1 | /*!\file Cfdragcoeffabsgrad.cpp
|
---|
2 | * \brief: Cfdragcoeffabsgrad Object
|
---|
3 | */
|
---|
4 |
|
---|
5 | /*Headers:*/
|
---|
6 | /*{{{*/
|
---|
7 | #ifdef HAVE_CONFIG_H
|
---|
8 | #include <config.h>
|
---|
9 | #else
|
---|
10 | #error "Cannot compile with HAVE_CONFIG_H symbol! run configure first!"
|
---|
11 | #endif
|
---|
12 |
|
---|
13 | #include "./classes.h"
|
---|
14 | #include "./ExternalResults/ExternalResult.h"
|
---|
15 | #include "./ExternalResults/Results.h"
|
---|
16 | #include "../datastructures/datastructures.h"
|
---|
17 | #include "./Elements/Element.h"
|
---|
18 | #include "./Elements/Elements.h"
|
---|
19 | #include "./FemModel.h"
|
---|
20 | #include "../modules/SurfaceAreax/SurfaceAreax.h"
|
---|
21 | #include "../classes/Params/Parameters.h"
|
---|
22 | #include "../classes/gauss/Gauss.h"
|
---|
23 | #include "./Inputs/DatasetInput.h"
|
---|
24 | /*}}}*/
|
---|
25 |
|
---|
26 | /*Cfdragcoeffabsgrad constructors, destructors :*/
|
---|
27 | Cfdragcoeffabsgrad::Cfdragcoeffabsgrad(){/*{{{*/
|
---|
28 |
|
---|
29 | this->definitionenum = -1;
|
---|
30 | this->name = NULL;
|
---|
31 | this->weights_enum = UNDEF;
|
---|
32 | this->misfit=0;
|
---|
33 | this->timepassedflag = false;
|
---|
34 | }
|
---|
35 | /*}}}*/
|
---|
36 | Cfdragcoeffabsgrad::Cfdragcoeffabsgrad(char* in_name, int in_definitionenum, int in_weights_enum, bool in_timepassedflag){/*{{{*/
|
---|
37 |
|
---|
38 | this->definitionenum=in_definitionenum;
|
---|
39 |
|
---|
40 | this->name = xNew<char>(strlen(in_name)+1);
|
---|
41 | xMemCpy<char>(this->name,in_name,strlen(in_name)+1);
|
---|
42 |
|
---|
43 | this->weights_enum=in_weights_enum;
|
---|
44 | this->timepassedflag=in_timepassedflag;
|
---|
45 |
|
---|
46 | this->misfit=0;
|
---|
47 | }
|
---|
48 | /*}}}*/
|
---|
49 | Cfdragcoeffabsgrad::~Cfdragcoeffabsgrad(){/*{{{*/
|
---|
50 | if(this->name)xDelete(this->name);
|
---|
51 | this->misfit=0;
|
---|
52 | }
|
---|
53 | /*}}}*/
|
---|
54 | /*Object virtual function resolutoin: */
|
---|
55 | Object* Cfdragcoeffabsgrad::copy() {/*{{{*/
|
---|
56 | Cfdragcoeffabsgrad* mf = new Cfdragcoeffabsgrad(this->name,this->definitionenum, this->weights_enum,this->timepassedflag);
|
---|
57 | mf->misfit=this->misfit;
|
---|
58 | return (Object*) mf;
|
---|
59 | }
|
---|
60 | /*}}}*/
|
---|
61 | void Cfdragcoeffabsgrad::DeepEcho(void){/*{{{*/
|
---|
62 | this->Echo();
|
---|
63 | }
|
---|
64 | /*}}}*/
|
---|
65 | void Cfdragcoeffabsgrad::Echo(void){/*{{{*/
|
---|
66 | _printf_(" Cfdragcoeffabsgrad: " << name << " " << this->definitionenum << "\n");
|
---|
67 | _printf_(" weights_enum: " << weights_enum << " " << EnumToStringx(weights_enum) << "\n");
|
---|
68 | _printf_(" timepassedflag: "<<timepassedflag<<"\n");
|
---|
69 | }
|
---|
70 | /*}}}*/
|
---|
71 | int Cfdragcoeffabsgrad::Id(void){/*{{{*/
|
---|
72 | return -1;
|
---|
73 | }
|
---|
74 | /*}}}*/
|
---|
75 | void Cfdragcoeffabsgrad::Marshall(MarshallHandle* marshallhandle){/*{{{*/
|
---|
76 | _error_("not implemented yet!");
|
---|
77 | }
|
---|
78 | /*}}}*/
|
---|
79 | int Cfdragcoeffabsgrad::ObjectEnum(void){/*{{{*/
|
---|
80 | return CfdragcoeffabsgradEnum;
|
---|
81 | }
|
---|
82 | /*}}}*/
|
---|
83 | /*Definition virtual function resolutoin: */
|
---|
84 | int Cfdragcoeffabsgrad::DefinitionEnum(){/*{{{*/
|
---|
85 | return this->definitionenum;
|
---|
86 | }
|
---|
87 | /*}}}*/
|
---|
88 | char* Cfdragcoeffabsgrad::Name(){/*{{{*/
|
---|
89 | char* name2=xNew<char>(strlen(this->name)+1);
|
---|
90 | xMemCpy(name2,this->name,strlen(this->name)+1);
|
---|
91 |
|
---|
92 | return name2;
|
---|
93 | }
|
---|
94 | /*}}}*/
|
---|
95 | IssmDouble Cfdragcoeffabsgrad::Response(FemModel* femmodel){/*{{{*/
|
---|
96 |
|
---|
97 | /*recover parameters: */
|
---|
98 | IssmDouble J=0.;
|
---|
99 | IssmDouble J_sum=0.;
|
---|
100 |
|
---|
101 | for(Object* & object : femmodel->elements->objects){
|
---|
102 | Element* element=xDynamicCast<Element*>(object);
|
---|
103 | J+=this->Cfdragcoeffabsgrad_Calculation(element,weights_enum);
|
---|
104 | }
|
---|
105 |
|
---|
106 | ISSM_MPI_Allreduce ( (void*)&J,(void*)&J_sum,1,ISSM_MPI_DOUBLE,ISSM_MPI_SUM,IssmComm::GetComm());
|
---|
107 | ISSM_MPI_Bcast(&J_sum,1,ISSM_MPI_DOUBLE,0,IssmComm::GetComm());
|
---|
108 | J=J_sum;
|
---|
109 |
|
---|
110 | timepassedflag = true;
|
---|
111 | return J;
|
---|
112 | }/*}}}*/
|
---|
113 | IssmDouble Cfdragcoeffabsgrad::Cfdragcoeffabsgrad_Calculation(Element* element, int weights_enum){/*{{{*/
|
---|
114 |
|
---|
115 | int domaintype,numcomponents;
|
---|
116 | IssmDouble Jelem=0.;
|
---|
117 | IssmDouble misfit,Jdet;
|
---|
118 | IssmDouble dp[2],weight;
|
---|
119 | IssmDouble* xyz_list = NULL;
|
---|
120 |
|
---|
121 | /*Get basal element*/
|
---|
122 | if(!element->IsOnBase()) return 0.;
|
---|
123 |
|
---|
124 | /*If on water, return 0: */
|
---|
125 | if(!element->IsIceInElement()) return 0.;
|
---|
126 |
|
---|
127 | /*Get problem dimension*/
|
---|
128 | element->FindParam(&domaintype,DomainTypeEnum);
|
---|
129 | switch(domaintype){
|
---|
130 | case Domain2DverticalEnum: numcomponents = 1; break;
|
---|
131 | case Domain3DEnum: numcomponents = 2; break;
|
---|
132 | case Domain2DhorizontalEnum: numcomponents = 2; break;
|
---|
133 | default: _error_("not supported yet");
|
---|
134 | }
|
---|
135 |
|
---|
136 | /*Spawn surface element*/
|
---|
137 | Element* basalelement = element->SpawnBasalElement();
|
---|
138 |
|
---|
139 | /* Get node coordinates*/
|
---|
140 | basalelement->GetVerticesCoordinates(&xyz_list);
|
---|
141 |
|
---|
142 | /*Get input if it already exists*/
|
---|
143 | DatasetInput *datasetinput = basalelement->GetDatasetInput(definitionenum); _assert_(datasetinput);
|
---|
144 | Input *drag_input = basalelement->GetInput(FrictionCoefficientEnum); _assert_(drag_input);
|
---|
145 |
|
---|
146 | /* Start looping on the number of gaussian points: */
|
---|
147 | Gauss* gauss=basalelement->NewGauss(2);
|
---|
148 | while(gauss->next()){
|
---|
149 |
|
---|
150 | /* Get Jacobian determinant: */
|
---|
151 | basalelement->JacobianDeterminant(&Jdet,xyz_list,gauss);
|
---|
152 |
|
---|
153 | /*Get all parameters at gaussian point*/
|
---|
154 | datasetinput->GetInputValue(&weight,gauss,WeightsSurfaceObservationEnum);
|
---|
155 | drag_input->GetInputDerivativeValue(&dp[0],xyz_list,gauss);
|
---|
156 |
|
---|
157 | /*Add to cost function*/
|
---|
158 | Jelem+=weight*.5*dp[0]*dp[0]*Jdet*gauss->weight;
|
---|
159 | if(numcomponents==2) Jelem+=weight*.5*dp[1]*dp[1]*Jdet*gauss->weight;
|
---|
160 | }
|
---|
161 |
|
---|
162 | /*clean up and Return: */
|
---|
163 | if(domaintype!=Domain2DhorizontalEnum){basalelement->DeleteMaterials(); delete basalelement;};
|
---|
164 | xDelete<IssmDouble>(xyz_list);
|
---|
165 | delete gauss;
|
---|
166 | return Jelem;
|
---|
167 | }/*}}}*/
|
---|