source: issm/trunk/src/c/classes/Cfdragcoeffabsgrad.cpp

Last change on this file was 28013, checked in by Mathieu Morlighem, 17 months ago

merged trunk-jpl and trunk for revision 28011

File size: 5.8 KB
Line 
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 :*/
27Cfdragcoeffabsgrad::Cfdragcoeffabsgrad(){/*{{{*/
28
29 this->definitionenum = -1;
30 this->name = NULL;
31 this->J = 0.;
32 this->firsttimepassed = false;
33}
34/*}}}*/
35Cfdragcoeffabsgrad::Cfdragcoeffabsgrad(char* in_name, int in_definitionenum){/*{{{*/
36
37 this->definitionenum=in_definitionenum;
38
39 this->name = xNew<char>(strlen(in_name)+1);
40 xMemCpy<char>(this->name,in_name,strlen(in_name)+1);
41
42 this->J = 0.;
43 this->firsttimepassed = false;
44}
45/*}}}*/
46Cfdragcoeffabsgrad::Cfdragcoeffabsgrad(char* in_name, int in_definitionenum, IssmDouble in_J){/*{{{*/
47
48 this->definitionenum=in_definitionenum;
49
50 this->name = xNew<char>(strlen(in_name)+1);
51 xMemCpy<char>(this->name,in_name,strlen(in_name)+1);
52
53 this->J = in_J;
54 this->firsttimepassed = false;
55}
56/*}}}*/
57Cfdragcoeffabsgrad::Cfdragcoeffabsgrad(char* in_name, int in_definitionenum, IssmDouble in_J, bool in_firsttimepassed){/*{{{*/
58
59 this->definitionenum=in_definitionenum;
60
61 this->name = xNew<char>(strlen(in_name)+1);
62 xMemCpy<char>(this->name,in_name,strlen(in_name)+1);
63
64 this->J = in_J;
65 this->firsttimepassed = in_firsttimepassed;
66}
67/*}}}*/
68Cfdragcoeffabsgrad::~Cfdragcoeffabsgrad(){/*{{{*/
69 if(this->name)xDelete(this->name);
70}
71/*}}}*/
72/*Object virtual function resolutoin: */
73Object* Cfdragcoeffabsgrad::copy() {/*{{{*/
74 Cfdragcoeffabsgrad* mf = new Cfdragcoeffabsgrad(this->name,this->definitionenum, this->J, this->firsttimepassed);
75 return (Object*) mf;
76}
77/*}}}*/
78void Cfdragcoeffabsgrad::DeepEcho(void){/*{{{*/
79 this->Echo();
80}
81/*}}}*/
82void Cfdragcoeffabsgrad::Echo(void){/*{{{*/
83 _printf_(" Cfdragcoeffabsgrad: " << name << " " << this->definitionenum << "\n");
84}
85/*}}}*/
86int Cfdragcoeffabsgrad::Id(void){/*{{{*/
87 return -1;
88}
89/*}}}*/
90void Cfdragcoeffabsgrad::Marshall(MarshallHandle* marshallhandle){/*{{{*/
91
92 /*ok, marshall operations: */
93 int object_enum=CfdragcoeffabsgradEnum;
94 marshallhandle->call(object_enum);
95
96 marshallhandle->call(this->definitionenum);
97 marshallhandle->call(this->name);
98 marshallhandle->call(this->J);
99 marshallhandle->call(this->firsttimepassed);
100}
101/*}}}*/
102int Cfdragcoeffabsgrad::ObjectEnum(void){/*{{{*/
103 return CfdragcoeffabsgradEnum;
104}
105/*}}}*/
106/*Definition virtual function resolutoin: */
107int Cfdragcoeffabsgrad::DefinitionEnum(){/*{{{*/
108 return this->definitionenum;
109}
110/*}}}*/
111char* Cfdragcoeffabsgrad::Name(){/*{{{*/
112 char* name2=xNew<char>(strlen(this->name)+1);
113 xMemCpy(name2,this->name,strlen(this->name)+1);
114
115 return name2;
116}
117/*}}}*/
118IssmDouble Cfdragcoeffabsgrad::Response(FemModel* femmodel){/*{{{*/
119
120 /*recover parameters: */
121 IssmDouble J_part=0.;
122 IssmDouble J_sum=0.;
123
124 if (!this->firsttimepassed){
125 for(Object* & object : femmodel->elements->objects){
126 Element* element=xDynamicCast<Element*>(object);
127 J_part+=this->Cfdragcoeffabsgrad_Calculation(element);
128 }
129
130 ISSM_MPI_Allreduce ( (void*)&J_part,(void*)&J_sum,1,ISSM_MPI_DOUBLE,ISSM_MPI_SUM,IssmComm::GetComm());
131 ISSM_MPI_Bcast(&J_sum,1,ISSM_MPI_DOUBLE,0,IssmComm::GetComm());
132 this->J=J_sum;
133
134 this->firsttimepassed = true;
135 }
136 return this->J;
137}/*}}}*/
138IssmDouble Cfdragcoeffabsgrad::Cfdragcoeffabsgrad_Calculation(Element* element){/*{{{*/
139
140 int domaintype,numcomponents,frictionlaw;
141 IssmDouble Jelem=0.;
142 IssmDouble Jdet;
143 IssmDouble dp[2],weight;
144 IssmDouble* xyz_list = NULL;
145
146 /*Get basal element*/
147 if(!element->IsOnBase()) return 0.;
148
149 /*If on water, return 0: */
150 if(!element->IsIceInElement()) return 0.;
151
152 /*Get problem dimension*/
153 element->FindParam(&domaintype,DomainTypeEnum);
154 switch(domaintype){
155 case Domain2DverticalEnum: numcomponents = 1; break;
156 case Domain3DEnum: numcomponents = 2; break;
157 case Domain2DhorizontalEnum: numcomponents = 2; break;
158 default: _error_("not supported yet");
159 }
160
161 /*Spawn surface element*/
162 Element* basalelement = element->SpawnBasalElement();
163
164 /* Get node coordinates*/
165 basalelement->GetVerticesCoordinates(&xyz_list);
166
167 /*Get input if it already exists*/
168 DatasetInput *datasetinput = basalelement->GetDatasetInput(definitionenum); _assert_(datasetinput);
169 Input *drag_input = NULL;
170
171 /* get the friction law: if 2-Weertman, 11-Schoof, 14-RegularizedCoulomb 15-RegularizedCoulomb2, which has a different names of C */
172 element->FindParam(&frictionlaw, FrictionLawEnum);
173 switch(frictionlaw) {
174 case 2:
175 case 11:
176 case 14:
177 case 15:
178 drag_input = basalelement->GetInput(FrictionCEnum); _assert_(drag_input);
179 break;
180 default:
181 drag_input = basalelement->GetInput(FrictionCoefficientEnum); _assert_(drag_input);
182 }
183
184
185 /* Start looping on the number of gaussian points: */
186 Gauss* gauss=basalelement->NewGauss(2);
187 while(gauss->next()){
188
189 /* Get Jacobian determinant: */
190 basalelement->JacobianDeterminant(&Jdet,xyz_list,gauss);
191
192 /*Get all parameters at gaussian point*/
193 datasetinput->GetInputValue(&weight,gauss,WeightsSurfaceObservationEnum);
194 drag_input->GetInputDerivativeValue(&dp[0],xyz_list,gauss);
195
196 /*Add to cost function*/
197 Jelem+=weight*.5*dp[0]*dp[0]*Jdet*gauss->weight;
198 if(numcomponents==2) Jelem+=weight*.5*dp[1]*dp[1]*Jdet*gauss->weight;
199 }
200
201 /*clean up and Return: */
202 if(basalelement->IsSpawnedElement()){basalelement->DeleteMaterials(); delete basalelement;};
203 xDelete<IssmDouble>(xyz_list);
204 delete gauss;
205 return Jelem;
206}/*}}}*/
Note: See TracBrowser for help on using the repository browser.