source: issm/trunk-jpl/src/c/analyses/EnthalpyAnalysis.cpp@ 21542

Last change on this file since 21542 was 21542, checked in by Mathieu Morlighem, 8 years ago

NEW: added Higher order finite elements for enthalpy

File size: 63.2 KB
Line 
1#include "./EnthalpyAnalysis.h"
2#include "../toolkits/toolkits.h"
3#include "../classes/classes.h"
4#include "../shared/shared.h"
5#include "../modules/modules.h"
6#include "../solutionsequences/solutionsequences.h"
7
8/*Model processing*/
9void EnthalpyAnalysis::CreateConstraints(Constraints* constraints,IoModel* iomodel){/*{{{*/
10
11 /*Intermediary*/
12 int count;
13 int M,N;
14 bool spcpresent = false;
15 IssmDouble heatcapacity;
16 IssmDouble referencetemperature;
17
18 /*Output*/
19 IssmDouble *spcvector = NULL;
20 IssmDouble* times=NULL;
21 IssmDouble* values=NULL;
22
23 /*Fetch parameters: */
24 iomodel->FindConstant(&heatcapacity,"md.materials.heatcapacity");
25 iomodel->FindConstant(&referencetemperature,"md.constants.referencetemperature");
26
27 /*return if 2d mesh*/
28 if(iomodel->domaintype==Domain2DhorizontalEnum) return;
29
30 /*Fetch data: */
31 iomodel->FetchData(&spcvector,&M,&N,"md.thermal.spctemperature");
32
33 //FIX ME: SHOULD USE IOMODELCREATECONSTRAINTS
34 /*Transient or static?:*/
35 if(M==iomodel->numberofvertices){
36 /*static: just create Constraints objects*/
37 count=0;
38
39 for(int i=0;i<iomodel->numberofvertices;i++){
40 /*keep only this partition's nodes:*/
41 if((iomodel->my_vertices[i])){
42
43 if (!xIsNan<IssmDouble>(spcvector[i])){
44
45 constraints->AddObject(new SpcStatic(iomodel->constraintcounter+count+1,iomodel->nodecounter+i+1,0,heatcapacity*(spcvector[i]-referencetemperature),EnthalpyAnalysisEnum));
46 count++;
47
48 }
49 }
50 }
51 }
52 else if (M==(iomodel->numberofvertices+1)){
53 /*transient: create transient SpcTransient objects. Same logic, except we need to retrieve
54 * various times and values to initialize an SpcTransient object: */
55 count=0;
56
57 /*figure out times: */
58 times=xNew<IssmDouble>(N);
59 for(int j=0;j<N;j++){
60 times[j]=spcvector[(M-1)*N+j];
61 }
62
63 /*Create constraints from x,y,z: */
64 for(int i=0;i<iomodel->numberofvertices;i++){
65
66 /*keep only this partition's nodes:*/
67 if((iomodel->my_vertices[i])){
68
69 /*figure out times and values: */
70 values=xNew<IssmDouble>(N);
71 spcpresent=false;
72 for(int j=0;j<N;j++){
73 values[j]=heatcapacity*(spcvector[i*N+j]-referencetemperature);
74 if(!xIsNan<IssmDouble>(values[j]))spcpresent=true; //NaN means no spc by default
75 }
76
77 if(spcpresent){
78 constraints->AddObject(new SpcTransient(iomodel->constraintcounter+count+1,iomodel->nodecounter+i+1,0,N,times,values,EnthalpyAnalysisEnum));
79 count++;
80 }
81 xDelete<IssmDouble>(values);
82 }
83 }
84 }
85 else{
86 _error_("Size of field " << EnumToStringx(ThermalSpctemperatureEnum) << " not supported");
87 }
88
89 /*Free ressources:*/
90 iomodel->DeleteData(spcvector,"md.thermal.spctemperature");
91 xDelete<IssmDouble>(times);
92 xDelete<IssmDouble>(values);
93}/*}}}*/
94void EnthalpyAnalysis::CreateLoads(Loads* loads, IoModel* iomodel){/*{{{*/
95
96 /*No loads */
97}/*}}}*/
98void EnthalpyAnalysis::CreateNodes(Nodes* nodes,IoModel* iomodel){/*{{{*/
99
100 int finiteelement;
101 iomodel->FindConstant(&finiteelement,"md.thermal.fe");
102
103 if(iomodel->domaintype==Domain3DEnum) iomodel->FetchData(2,"md.mesh.vertexonbase","md.mesh.vertexonsurface");
104 ::CreateNodes(nodes,iomodel,EnthalpyAnalysisEnum,finiteelement);
105 iomodel->DeleteData(2,"md.mesh.vertexonbase","md.mesh.vertexonsurface");
106}/*}}}*/
107int EnthalpyAnalysis::DofsPerNode(int** doflist,int domaintype,int approximation){/*{{{*/
108 return 1;
109}/*}}}*/
110void EnthalpyAnalysis::UpdateElements(Elements* elements,IoModel* iomodel,int analysis_counter,int analysis_type){/*{{{*/
111
112 bool dakota_analysis,ismovingfront,isenthalpy;
113 int frictionlaw,basalforcing_model,materialstype;
114 int FrictionCoupling;
115
116 /*Now, is the model 3d? otherwise, do nothing: */
117 if(iomodel->domaintype==Domain2DhorizontalEnum)return;
118
119 /*Is enthalpy requested?*/
120 iomodel->FindConstant(&isenthalpy,"md.thermal.isenthalpy");
121 if(!isenthalpy) return;
122
123 /*Fetch data needed: */
124 iomodel->FetchData(3,"md.initialization.temperature","md.initialization.waterfraction","md.initialization.pressure");
125
126 /*Finite element type*/
127 int finiteelement;
128 iomodel->FindConstant(&finiteelement,"md.thermal.fe");
129
130 /*Update elements: */
131 int counter=0;
132 for(int i=0;i<iomodel->numberofelements;i++){
133 if(iomodel->my_elements[i]){
134 Element* element=(Element*)elements->GetObjectByOffset(counter);
135 element->Update(i,iomodel,analysis_counter,analysis_type,finiteelement);
136 counter++;
137 }
138 }
139
140 iomodel->FindConstant(&dakota_analysis,"md.qmu.isdakota");
141 iomodel->FindConstant(&ismovingfront,"md.transient.ismovingfront");
142 iomodel->FindConstant(&frictionlaw,"md.friction.law");
143 iomodel->FindConstant(&materialstype,"md.materials.type");
144
145 iomodel->FetchDataToInput(elements,"md.geometry.thickness",ThicknessEnum);
146 iomodel->FetchDataToInput(elements,"md.geometry.surface",SurfaceEnum);
147 iomodel->FetchDataToInput(elements,"md.slr.sealevel",SealevelEnum,0);
148 iomodel->FetchDataToInput(elements,"md.geometry.base",BaseEnum);
149 iomodel->FetchDataToInput(elements,"md.mask.ice_levelset",MaskIceLevelsetEnum);
150 iomodel->FetchDataToInput(elements,"md.mask.groundedice_levelset",MaskGroundediceLevelsetEnum);
151 if(iomodel->domaintype!=Domain2DhorizontalEnum){
152 iomodel->FetchDataToInput(elements,"md.mesh.vertexonbase",MeshVertexonbaseEnum);
153 iomodel->FetchDataToInput(elements,"md.mesh.vertexonsurface",MeshVertexonsurfaceEnum);
154 }
155 iomodel->FetchDataToInput(elements,"md.initialization.pressure",PressureEnum);
156 iomodel->FetchDataToInput(elements,"md.initialization.temperature",TemperatureEnum);
157 iomodel->FetchDataToInput(elements,"md.initialization.waterfraction",WaterfractionEnum);
158 iomodel->FetchDataToInput(elements,"md.initialization.enthalpy",EnthalpyEnum);
159 iomodel->FetchDataToInput(elements,"md.initialization.watercolumn",WatercolumnEnum);
160 iomodel->FetchDataToInput(elements,"md.basalforcings.groundedice_melting_rate",BasalforcingsGroundediceMeltingRateEnum);
161 iomodel->FetchDataToInput(elements,"md.initialization.vx",VxEnum);
162 iomodel->FetchDataToInput(elements,"md.initialization.vy",VyEnum);
163 iomodel->FetchDataToInput(elements,"md.initialization.vz",VzEnum);
164 InputUpdateFromConstantx(elements,0.,VxMeshEnum);
165 InputUpdateFromConstantx(elements,0.,VyMeshEnum);
166 InputUpdateFromConstantx(elements,0.,VzMeshEnum);
167 if(ismovingfront){
168 iomodel->FetchDataToInput(elements,"md.mesh.vertexonbase",MeshVertexonbaseEnum); // required for updating active nodes
169 }
170
171 /*Basal forcings variables*/
172 iomodel->FindConstant(&basalforcing_model,"md.basalforcings.model");
173 switch(basalforcing_model){
174 case MantlePlumeGeothermalFluxEnum:
175 break;
176 default:
177 iomodel->FetchDataToInput(elements,"md.basalforcings.geothermalflux",BasalforcingsGeothermalfluxEnum);
178 break;
179 }
180
181 /*Rheology type*/
182 iomodel->FetchDataToInput(elements,"md.materials.rheology_B",MaterialsRheologyBEnum);
183 switch(materialstype){
184 case MatenhancediceEnum:
185 iomodel->FetchDataToInput(elements,"md.materials.rheology_n",MaterialsRheologyNEnum);
186 iomodel->FetchDataToInput(elements,"md.materials.rheology_E",MaterialsRheologyEEnum);
187 break;
188 case MatdamageiceEnum:
189 iomodel->FetchDataToInput(elements,"md.materials.rheology_n",MaterialsRheologyNEnum);
190 break;
191 case MatestarEnum:
192 iomodel->FetchDataToInput(elements,"md.materials.rheology_Ec",MaterialsRheologyEcEnum);
193 iomodel->FetchDataToInput(elements,"md.materials.rheology_Es",MaterialsRheologyEsEnum);
194 break;
195 case MaticeEnum:
196 iomodel->FetchDataToInput(elements,"md.materials.rheology_n",MaterialsRheologyNEnum);
197 break;
198 default:
199 _error_("not supported");
200 }
201
202 /*Friction law variables*/
203 switch(frictionlaw){
204 case 1:
205 iomodel->FetchDataToInput(elements,"md.friction.coefficient",FrictionCoefficientEnum);
206 iomodel->FetchDataToInput(elements,"md.friction.p",FrictionPEnum);
207 iomodel->FetchDataToInput(elements,"md.friction.q",FrictionQEnum);
208 break;
209 case 2:
210 iomodel->FetchDataToInput(elements,"md.friction.C",FrictionCEnum);
211 iomodel->FetchDataToInput(elements,"md.friction.m",FrictionMEnum);
212 break;
213 case 3:
214 iomodel->FindConstant(&FrictionCoupling,"md.friction.coupling");
215 iomodel->FetchDataToInput(elements,"md.friction.C",FrictionCEnum);
216 iomodel->FetchDataToInput(elements,"md.friction.As",FrictionAsEnum);
217 iomodel->FetchDataToInput(elements,"md.friction.q",FrictionQEnum);
218 if (FrictionCoupling==0){
219 iomodel->FetchDataToInput(elements,"md.friction.effective_pressure",FrictionEffectivePressureEnum);
220 }
221 break;
222 case 4:
223 iomodel->FetchDataToInput(elements,"md.friction.coefficient",FrictionCoefficientEnum);
224 iomodel->FetchDataToInput(elements,"md.friction.p",FrictionPEnum);
225 iomodel->FetchDataToInput(elements,"md.friction.q",FrictionQEnum);
226 iomodel->FetchDataToInput(elements,"md.initialization.pressure",PressureEnum);
227 iomodel->FetchDataToInput(elements,"md.initialization.temperature",TemperatureEnum);
228 break;
229 case 5:
230 iomodel->FetchDataToInput(elements,"md.friction.coefficient",FrictionCoefficientEnum);
231 iomodel->FetchDataToInput(elements,"md.friction.p",FrictionPEnum);
232 iomodel->FetchDataToInput(elements,"md.friction.q",FrictionQEnum);
233 iomodel->FetchDataToInput(elements,"md.friction.water_layer",FrictionWaterLayerEnum);
234 break;
235 case 6:
236 iomodel->FetchDataToInput(elements,"md.friction.C",FrictionCEnum);
237 iomodel->FetchDataToInput(elements,"md.friction.m",FrictionMEnum);
238 iomodel->FetchDataToInput(elements,"md.initialization.pressure",PressureEnum);
239 iomodel->FetchDataToInput(elements,"md.initialization.temperature",TemperatureEnum);
240 break;
241 default:
242 _error_("not supported");
243 }
244
245 /*Free data: */
246 iomodel->DeleteData(3,"md.initialization.temperature","md.initialization.waterfraction","md.initialization.pressure");
247
248}/*}}}*/
249void EnthalpyAnalysis::UpdateParameters(Parameters* parameters,IoModel* iomodel,int solution_enum,int analysis_enum){/*{{{*/
250
251 int numoutputs;
252 char** requestedoutputs = NULL;
253
254 parameters->AddObject(iomodel->CopyConstantObject("md.thermal.stabilization",ThermalStabilizationEnum));
255 parameters->AddObject(iomodel->CopyConstantObject("md.thermal.maxiter",ThermalMaxiterEnum));
256 parameters->AddObject(iomodel->CopyConstantObject("md.thermal.reltol",ThermalReltolEnum));
257 parameters->AddObject(iomodel->CopyConstantObject("md.thermal.isenthalpy",ThermalIsenthalpyEnum));
258 parameters->AddObject(iomodel->CopyConstantObject("md.thermal.isdynamicbasalspc",ThermalIsdynamicbasalspcEnum));
259 parameters->AddObject(iomodel->CopyConstantObject("md.friction.law",FrictionLawEnum));
260
261 iomodel->FindConstant(&requestedoutputs,&numoutputs,"md.thermal.requested_outputs");
262 parameters->AddObject(new IntParam(ThermalNumRequestedOutputsEnum,numoutputs));
263 if(numoutputs)parameters->AddObject(new StringArrayParam(ThermalRequestedOutputsEnum,requestedoutputs,numoutputs));
264 iomodel->DeleteData(&requestedoutputs,numoutputs,"md.thermal.requested_outputs");
265
266 /*Deal with friction parameters*/
267 int frictionlaw;
268 iomodel->FindConstant(&frictionlaw,"md.friction.law");
269 if(frictionlaw==4 || frictionlaw==6) parameters->AddObject(iomodel->CopyConstantObject("md.friction.gamma",FrictionGammaEnum));
270 if(frictionlaw==3) parameters->AddObject(iomodel->CopyConstantObject("md.friction.coupling",FrictionCouplingEnum));
271}/*}}}*/
272
273/*Finite Element Analysis*/
274void EnthalpyAnalysis::ApplyBasalConstraints(IssmDouble* serial_spc,Element* element){/*{{{*/
275
276 /* Do not check if ice in element, this may lead to inconsistencies between cpu partitions */
277 /* Only update constraints at the base. */
278 if(!(element->IsOnBase())) return;
279
280 /*Intermediary*/
281 bool isdynamicbasalspc;
282 int numindices;
283 int *indices = NULL;
284 Node* node = NULL;
285 IssmDouble pressure;
286
287 /*Check wether dynamic basal boundary conditions are activated */
288 element->FindParam(&isdynamicbasalspc,ThermalIsdynamicbasalspcEnum);
289 if(!isdynamicbasalspc) return;
290
291 /*Get parameters and inputs: */
292 Input* pressure_input = element->GetInput(PressureEnum); _assert_(pressure_input);
293
294 /*Fetch indices of basal & surface nodes for this finite element*/
295 Penta *penta = (Penta *) element; // TODO: add Basal-/SurfaceNodeIndices to element.h, and change this to Element*
296 penta->BasalNodeIndices(&numindices,&indices,element->GetElementType());
297
298 GaussPenta* gauss=new GaussPenta();
299 for(int i=0;i<numindices;i++){
300 gauss->GaussNode(element->GetElementType(),indices[i]);
301
302 pressure_input->GetInputValue(&pressure,gauss);
303
304 /*apply or release spc*/
305 node=element->GetNode(indices[i]);
306 if(!node->IsActive()) continue;
307 if(serial_spc[node->Sid()]==1.){
308 pressure_input->GetInputValue(&pressure, gauss);
309 node->ApplyConstraint(0,PureIceEnthalpy(element,pressure));
310 }
311 else {
312 node->DofInFSet(0);
313 }
314 }
315
316 /*Free ressources:*/
317 xDelete<int>(indices);
318 delete gauss;
319}/*}}}*/
320void EnthalpyAnalysis::ComputeBasalMeltingrate(FemModel* femmodel){/*{{{*/
321 /*Compute basal melting rates: */
322 for(int i=0;i<femmodel->elements->Size();i++){
323 Element* element=xDynamicCast<Element*>(femmodel->elements->GetObjectByOffset(i));
324 ComputeBasalMeltingrate(element);
325 }
326}/*}}}*/
327void EnthalpyAnalysis::ComputeBasalMeltingrate(Element* element){/*{{{*/
328 /*Calculate the basal melt rates of the enthalpy model after Aschwanden 2012*/
329 /* melting rate is positive when melting, negative when refreezing*/
330
331 /* Check if ice in element */
332 if(!element->IsIceInElement()) return;
333
334 /* Only compute melt rates at the base of grounded ice*/
335 if(!element->IsOnBase() || element->IsFloating()) return;
336
337 /* Intermediaries */
338 bool converged;
339 const int dim=3;
340 int i,is,state;
341 int vertexdown,vertexup,numvertices,numsegments;
342 int enthalpy_enum;
343 IssmDouble vec_heatflux[dim],normal_base[dim],d1enthalpy[dim],d1pressure[dim];
344 IssmDouble basalfriction,alpha2,geothermalflux,heatflux;
345 IssmDouble dt,yts;
346 IssmDouble melting_overshoot,lambda;
347 IssmDouble vx,vy,vz;
348 IssmDouble *xyz_list = NULL;
349 IssmDouble *xyz_list_base = NULL;
350 int *pairindices = NULL;
351
352 /*Fetch parameters*/
353 element->GetVerticesCoordinates(&xyz_list);
354 element->GetVerticesCoordinatesBase(&xyz_list_base);
355 element->GetInputValue(&converged,ConvergedEnum);
356 element->FindParam(&dt,TimesteppingTimeStepEnum);
357 element->FindParam(&yts, ConstantsYtsEnum);
358
359 if(dt==0. && !converged) enthalpy_enum=EnthalpyPicardEnum;
360 else enthalpy_enum=EnthalpyEnum;
361
362 IssmDouble latentheat = element->GetMaterialParameter(MaterialsLatentheatEnum);
363 IssmDouble rho_ice = element->GetMaterialParameter(MaterialsRhoIceEnum);
364 IssmDouble rho_water = element->GetMaterialParameter(MaterialsRhoFreshwaterEnum);
365 IssmDouble beta = element->GetMaterialParameter(MaterialsBetaEnum);
366 IssmDouble kappa = EnthalpyDiffusionParameterVolume(element,enthalpy_enum); _assert_(kappa>=0.);
367 IssmDouble kappa_mix;
368
369 /*retrieve inputs*/
370 Input* enthalpy_input = element->GetInput(enthalpy_enum); _assert_(enthalpy_input);
371 Input* pressure_input = element->GetInput(PressureEnum); _assert_(pressure_input);
372 Input* geothermalflux_input = element->GetInput(BasalforcingsGeothermalfluxEnum); _assert_(geothermalflux_input);
373 Input* vx_input = element->GetInput(VxEnum); _assert_(vx_input);
374 Input* vy_input = element->GetInput(VyEnum); _assert_(vy_input);
375 Input* vz_input = element->GetInput(VzEnum); _assert_(vz_input);
376
377 /*Build friction element, needed later: */
378 Friction* friction=new Friction(element,dim);
379
380 /******** MELTING RATES ************************************//*{{{*/
381 element->NormalBase(&normal_base[0],xyz_list_base);
382 element->VerticalSegmentIndices(&pairindices,&numsegments);
383 IssmDouble* meltingrate_enthalpy = xNew<IssmDouble>(numsegments);
384 IssmDouble* heating = xNew<IssmDouble>(numsegments);
385
386 numvertices=element->GetNumberOfVertices();
387 IssmDouble* enthalpies = xNew<IssmDouble>(numvertices);
388 IssmDouble* pressures = xNew<IssmDouble>(numvertices);
389 IssmDouble* watercolumns = xNew<IssmDouble>(numvertices);
390 IssmDouble* basalmeltingrates = xNew<IssmDouble>(numvertices);
391 element->GetInputListOnVertices(enthalpies,enthalpy_enum);
392 element->GetInputListOnVertices(pressures,PressureEnum);
393 element->GetInputListOnVertices(watercolumns,WatercolumnEnum);
394 element->GetInputListOnVertices(basalmeltingrates,BasalforcingsGroundediceMeltingRateEnum);
395
396 Gauss* gauss=element->NewGauss();
397 for(is=0;is<numsegments;is++){
398 vertexdown = pairindices[is*2+0];
399 vertexup = pairindices[is*2+1];
400 gauss->GaussVertex(vertexdown);
401
402 state=GetThermalBasalCondition(element, enthalpies[vertexdown], enthalpies[vertexup], pressures[vertexdown], pressures[vertexup], watercolumns[vertexdown], basalmeltingrates[vertexdown]);
403 switch (state) {
404 case 0:
405 // cold, dry base: apply basal surface forcing
406 for(i=0;i<3;i++) vec_heatflux[i]=0.;
407 break;
408 case 1: case 2: case 3:
409 // case 1 : cold, wet base: keep at pressure melting point
410 // case 2: temperate, thin refreezing base: release spc
411 // case 3: temperate, thin melting base: set spc
412 enthalpy_input->GetInputDerivativeValue(&d1enthalpy[0],xyz_list,gauss);
413 for(i=0;i<3;i++) vec_heatflux[i]=-kappa*d1enthalpy[i];
414 break;
415 case 4:
416 // temperate, thick melting base: set grad H*n=0
417 kappa_mix=GetWetIceConductivity(element, enthalpies[vertexdown], pressures[vertexdown]);
418 pressure_input->GetInputDerivativeValue(&d1pressure[0],xyz_list,gauss);
419 for(i=0;i<3;i++) vec_heatflux[i]=kappa_mix*beta*d1pressure[i];
420 break;
421 default:
422 _printf0_(" unknown thermal basal state found!");
423 }
424 if(state==0) meltingrate_enthalpy[is]=0.;
425 else{
426 /*heat flux along normal*/
427 heatflux=0.;
428 for(i=0;i<3;i++) heatflux+=(vec_heatflux[i])*normal_base[i];
429
430 /*basal friction*/
431 friction->GetAlpha2(&alpha2,gauss);
432 vx_input->GetInputValue(&vx,gauss); vy_input->GetInputValue(&vy,gauss); vz_input->GetInputValue(&vz,gauss);
433 basalfriction=alpha2*(vx*vx + vy*vy + vz*vz);
434 geothermalflux_input->GetInputValue(&geothermalflux,gauss);
435 /* -Mb= Fb-(q-q_geo)/((1-w)*L*rho), and (1-w)*rho=rho_ice, cf Aschwanden 2012, eqs.1, 2, 66*/
436 heating[is]=(heatflux+basalfriction+geothermalflux);
437 meltingrate_enthalpy[is]=heating[is]/(latentheat*rho_ice); // m/s water equivalent
438 }
439 }/*}}}*/
440
441 /******** UPDATE MELTINGRATES AND WATERCOLUMN **************//*{{{*/
442 for(is=0;is<numsegments;is++){
443 vertexdown = pairindices[is*2+0];
444 vertexup = pairindices[is*2+1];
445 if(dt!=0.){
446 if(watercolumns[vertexdown]+meltingrate_enthalpy[is]*dt<0.){ // prevent too much freeze on
447 lambda = -watercolumns[vertexdown]/(dt*meltingrate_enthalpy[is]); _assert_(lambda>=0.); _assert_(lambda<1.);
448 watercolumns[vertexdown]=0.;
449 basalmeltingrates[vertexdown]=lambda*meltingrate_enthalpy[is]; // restrict freeze on only to size of watercolumn
450 enthalpies[vertexdown]+=(1.-lambda)*dt/yts*meltingrate_enthalpy[is]*latentheat*rho_ice; // use rest of energy to cool down base: dE=L*m, m=(1-lambda)*meltingrate*rho_ice
451 }
452 else{
453 basalmeltingrates[vertexdown]=meltingrate_enthalpy[is];
454 watercolumns[vertexdown]+=dt*meltingrate_enthalpy[is];
455 }
456 }
457 else{
458 basalmeltingrates[vertexdown]=meltingrate_enthalpy[is];
459 if(watercolumns[vertexdown]+meltingrate_enthalpy[is]<0.)
460 watercolumns[vertexdown]=0.;
461 else
462 watercolumns[vertexdown]+=meltingrate_enthalpy[is];
463 }
464 basalmeltingrates[vertexdown]*=rho_water/rho_ice; // convert meltingrate from water to ice equivalent
465 _assert_(watercolumns[vertexdown]>=0.);
466 }/*}}}*/
467
468 /*feed updated variables back into model*/
469 if(dt!=0.){
470 element->AddInput(enthalpy_enum,enthalpies,P1Enum); //TODO: distinguis for steadystate and transient run
471 element->AddInput(WatercolumnEnum,watercolumns,P1Enum);
472 }
473 element->AddInput(BasalforcingsGroundediceMeltingRateEnum,basalmeltingrates,P1Enum);
474
475 /*Clean up and return*/
476 delete gauss;
477 delete friction;
478 xDelete<int>(pairindices);
479 xDelete<IssmDouble>(enthalpies);
480 xDelete<IssmDouble>(pressures);
481 xDelete<IssmDouble>(watercolumns);
482 xDelete<IssmDouble>(basalmeltingrates);
483 xDelete<IssmDouble>(meltingrate_enthalpy);
484 xDelete<IssmDouble>(heating);
485 xDelete<IssmDouble>(xyz_list);
486 xDelete<IssmDouble>(xyz_list_base);
487}/*}}}*/
488void EnthalpyAnalysis::Core(FemModel* femmodel){/*{{{*/
489
490 IssmDouble dt;
491 bool isdynamicbasalspc;
492
493 femmodel->parameters->FindParam(&dt,TimesteppingTimeStepEnum);
494 femmodel->parameters->FindParam(&isdynamicbasalspc,ThermalIsdynamicbasalspcEnum);
495
496 if(VerboseSolution()) _printf0_(" computing enthalpy\n");
497 femmodel->SetCurrentConfiguration(EnthalpyAnalysisEnum);
498 if((dt>0.) && isdynamicbasalspc) UpdateBasalConstraints(femmodel);
499 solutionsequence_thermal_nonlinear(femmodel);
500
501 /*transfer enthalpy to enthalpy picard for the next step: */
502 InputDuplicatex(femmodel,EnthalpyEnum,EnthalpyPicardEnum);
503
504 PostProcessing(femmodel);
505
506}/*}}}*/
507ElementVector* EnthalpyAnalysis::CreateDVector(Element* element){/*{{{*/
508 /*Default, return NULL*/
509 return NULL;
510}/*}}}*/
511ElementMatrix* EnthalpyAnalysis::CreateJacobianMatrix(Element* element){/*{{{*/
512_error_("Not implemented");
513}/*}}}*/
514ElementMatrix* EnthalpyAnalysis::CreateKMatrix(Element* element){/*{{{*/
515
516 /* Check if ice in element */
517 if(!element->IsIceInElement()) return NULL;
518
519 /*compute all stiffness matrices for this element*/
520 ElementMatrix* Ke1=CreateKMatrixVolume(element);
521 ElementMatrix* Ke2=CreateKMatrixShelf(element);
522 ElementMatrix* Ke =new ElementMatrix(Ke1,Ke2);
523
524 /*clean-up and return*/
525 delete Ke1;
526 delete Ke2;
527 return Ke;
528}/*}}}*/
529ElementMatrix* EnthalpyAnalysis::CreateKMatrixVolume(Element* element){/*{{{*/
530
531 /* Check if ice in element */
532 if(!element->IsIceInElement()) return NULL;
533
534 /*Intermediaries */
535 int stabilization;
536 IssmDouble Jdet,dt,u,v,w,um,vm,wm,vel;
537 IssmDouble h,hx,hy,hz,vx,vy,vz;
538 IssmDouble tau_parameter,diameter;
539 IssmDouble D_scalar;
540 IssmDouble* xyz_list = NULL;
541
542 /*Fetch number of nodes and dof for this finite element*/
543 int numnodes = element->GetNumberOfNodes();
544
545 /*Initialize Element vector and other vectors*/
546 ElementMatrix* Ke = element->NewElementMatrix();
547 IssmDouble* basis = xNew<IssmDouble>(numnodes);
548 IssmDouble* dbasis = xNew<IssmDouble>(3*numnodes);
549 IssmDouble* B = xNew<IssmDouble>(3*numnodes);
550 IssmDouble* Bprime = xNew<IssmDouble>(3*numnodes);
551 IssmDouble D[3][3] = {0.};
552 IssmDouble K[3][3];
553
554 /*Retrieve all inputs and parameters*/
555 element->GetVerticesCoordinates(&xyz_list);
556 element->FindParam(&dt,TimesteppingTimeStepEnum);
557 element->FindParam(&stabilization,ThermalStabilizationEnum);
558 IssmDouble rho_water = element->GetMaterialParameter(MaterialsRhoSeawaterEnum);
559 IssmDouble rho_ice = element->GetMaterialParameter(MaterialsRhoIceEnum);
560 IssmDouble gravity = element->GetMaterialParameter(ConstantsGEnum);
561 IssmDouble heatcapacity = element->GetMaterialParameter(MaterialsHeatcapacityEnum);
562 IssmDouble thermalconductivity = element->GetMaterialParameter(MaterialsThermalconductivityEnum);
563 Input* vx_input = element->GetInput(VxEnum); _assert_(vx_input);
564 Input* vy_input = element->GetInput(VyEnum); _assert_(vy_input);
565 Input* vz_input = element->GetInput(VzEnum); _assert_(vz_input);
566 Input* vxm_input = element->GetInput(VxMeshEnum); _assert_(vxm_input);
567 Input* vym_input = element->GetInput(VyMeshEnum); _assert_(vym_input);
568 Input* vzm_input = element->GetInput(VzMeshEnum); _assert_(vzm_input);
569 if(stabilization==2) diameter=element->MinEdgeLength(xyz_list);
570
571 /*Enthalpy diffusion parameter*/
572 IssmDouble kappa=this->EnthalpyDiffusionParameterVolume(element,EnthalpyPicardEnum); _assert_(kappa>=0.);
573
574 /* Start looping on the number of gaussian points: */
575 Gauss* gauss=element->NewGauss(4);
576 for(int ig=gauss->begin();ig<gauss->end();ig++){
577 gauss->GaussPoint(ig);
578
579 element->JacobianDeterminant(&Jdet,xyz_list,gauss);
580 D_scalar=gauss->weight*Jdet;
581 if(dt!=0.) D_scalar=D_scalar*dt;
582
583 /*Conduction: */
584 GetBConduct(B,element,xyz_list,gauss);
585 D[0][0]=D_scalar*kappa/rho_ice;
586 D[1][1]=D_scalar*kappa/rho_ice;
587 D[2][2]=D_scalar*kappa/rho_ice;
588 TripleMultiply(B,3,numnodes,1,
589 &D[0][0],3,3,0,
590 B,3,numnodes,0,
591 &Ke->values[0],1);
592
593 /*Advection: */
594 GetBAdvec(B,element,xyz_list,gauss);
595 GetBAdvecprime(Bprime,element,xyz_list,gauss);
596 vx_input->GetInputValue(&u,gauss); vxm_input->GetInputValue(&um,gauss); vx=u-um;
597 vy_input->GetInputValue(&v,gauss); vym_input->GetInputValue(&vm,gauss); vy=v-vm;
598 vz_input->GetInputValue(&w,gauss); vzm_input->GetInputValue(&wm,gauss); vz=w-wm;
599 D[0][0]=D_scalar*vx;
600 D[1][1]=D_scalar*vy;
601 D[2][2]=D_scalar*vz;
602 TripleMultiply(B,3,numnodes,1,
603 &D[0][0],3,3,0,
604 Bprime,3,numnodes,0,
605 &Ke->values[0],1);
606
607 /*Transient: */
608 if(dt!=0.){
609 D_scalar=gauss->weight*Jdet;
610 element->NodalFunctions(basis,gauss);
611 TripleMultiply(basis,numnodes,1,0,
612 &D_scalar,1,1,0,
613 basis,1,numnodes,0,
614 &Ke->values[0],1);
615 D_scalar=D_scalar*dt;
616 }
617
618 /*Artificial diffusivity*/
619 if(stabilization==1){
620 element->ElementSizes(&hx,&hy,&hz);
621 vel=sqrt(vx*vx + vy*vy + vz*vz)+1.e-14;
622 h=sqrt( pow(hx*vx/vel,2) + pow(hy*vy/vel,2) + pow(hz*vz/vel,2));
623 K[0][0]=h/(2.*vel)*fabs(vx*vx); K[0][1]=h/(2.*vel)*fabs(vx*vy); K[0][2]=h/(2.*vel)*fabs(vx*vz);
624 K[1][0]=h/(2.*vel)*fabs(vy*vx); K[1][1]=h/(2.*vel)*fabs(vy*vy); K[1][2]=h/(2.*vel)*fabs(vy*vz);
625 K[2][0]=h/(2.*vel)*fabs(vz*vx); K[2][1]=h/(2.*vel)*fabs(vz*vy); K[2][2]=h/(2.*vel)*fabs(vz*vz);
626 for(int i=0;i<3;i++) for(int j=0;j<3;j++) K[i][j] = D_scalar*K[i][j];
627
628 GetBAdvecprime(Bprime,element,xyz_list,gauss);
629 TripleMultiply(Bprime,3,numnodes,1,
630 &K[0][0],3,3,0,
631 Bprime,3,numnodes,0,
632 &Ke->values[0],1);
633 }
634 else if(stabilization==2){
635 element->NodalFunctionsDerivatives(dbasis,xyz_list,gauss);
636 tau_parameter=element->StabilizationParameter(u-um,v-vm,w-wm,diameter,kappa/rho_ice);
637 for(int i=0;i<numnodes;i++){
638 for(int j=0;j<numnodes;j++){
639 Ke->values[i*numnodes+j]+=tau_parameter*D_scalar*
640 ((u-um)*dbasis[0*numnodes+i]+(v-vm)*dbasis[1*numnodes+i]+(w-wm)*dbasis[2*numnodes+i])*((u-um)*dbasis[0*numnodes+j]+(v-vm)*dbasis[1*numnodes+j]+(w-wm)*dbasis[2*numnodes+j]);
641 }
642 }
643 if(dt!=0.){
644 D_scalar=gauss->weight*Jdet;
645 for(int i=0;i<numnodes;i++){
646 for(int j=0;j<numnodes;j++){
647 Ke->values[i*numnodes+j]+=tau_parameter*D_scalar*basis[j]*((u-um)*dbasis[0*numnodes+i]+(v-vm)*dbasis[1*numnodes+i]+(w-wm)*dbasis[2*numnodes+i]);
648 }
649 }
650 }
651 }
652 }
653
654 /*Clean up and return*/
655 xDelete<IssmDouble>(xyz_list);
656 xDelete<IssmDouble>(basis);
657 xDelete<IssmDouble>(dbasis);
658 xDelete<IssmDouble>(B);
659 xDelete<IssmDouble>(Bprime);
660 delete gauss;
661 return Ke;
662}/*}}}*/
663ElementMatrix* EnthalpyAnalysis::CreateKMatrixShelf(Element* element){/*{{{*/
664
665 /* Check if ice in element */
666 if(!element->IsIceInElement()) return NULL;
667
668 /*Initialize Element matrix and return if necessary*/
669 if(!element->IsOnBase() || !element->IsFloating()) return NULL;
670
671 /*Intermediaries*/
672 IssmDouble dt,Jdet,D;
673 IssmDouble *xyz_list_base = NULL;
674
675 /*Fetch number of nodes for this finite element*/
676 int numnodes = element->GetNumberOfNodes();
677
678 /*Initialize vectors*/
679 ElementMatrix* Ke = element->NewElementMatrix();
680 IssmDouble* basis = xNew<IssmDouble>(numnodes);
681
682 /*Retrieve all inputs and parameters*/
683 element->GetVerticesCoordinatesBase(&xyz_list_base);
684 element->FindParam(&dt,TimesteppingTimeStepEnum);
685 IssmDouble gravity = element->GetMaterialParameter(ConstantsGEnum);
686 IssmDouble rho_water = element->GetMaterialParameter(MaterialsRhoSeawaterEnum);
687 IssmDouble rho_ice = element->GetMaterialParameter(MaterialsRhoIceEnum);
688 IssmDouble heatcapacity = element->GetMaterialParameter(MaterialsHeatcapacityEnum);
689 IssmDouble mixed_layer_capacity= element->GetMaterialParameter(MaterialsMixedLayerCapacityEnum);
690 IssmDouble thermal_exchange_vel= element->GetMaterialParameter(MaterialsThermalExchangeVelocityEnum);
691
692 /* Start looping on the number of gaussian points: */
693 Gauss* gauss=element->NewGaussBase(4);
694 for(int ig=gauss->begin();ig<gauss->end();ig++){
695 gauss->GaussPoint(ig);
696
697 element->JacobianDeterminantBase(&Jdet,xyz_list_base,gauss);
698 element->NodalFunctions(basis,gauss);
699
700 D=gauss->weight*Jdet*rho_water*mixed_layer_capacity*thermal_exchange_vel/(heatcapacity*rho_ice);
701 if(reCast<bool,IssmDouble>(dt)) D=dt*D;
702 TripleMultiply(basis,numnodes,1,0,
703 &D,1,1,0,
704 basis,1,numnodes,0,
705 &Ke->values[0],1);
706
707 }
708
709 /*Clean up and return*/
710 delete gauss;
711 xDelete<IssmDouble>(basis);
712 xDelete<IssmDouble>(xyz_list_base);
713 return Ke;
714}/*}}}*/
715ElementVector* EnthalpyAnalysis::CreatePVector(Element* element){/*{{{*/
716
717 /* Check if ice in element */
718 if(!element->IsIceInElement()) return NULL;
719
720 /*compute all load vectors for this element*/
721 ElementVector* pe1=CreatePVectorVolume(element);
722 ElementVector* pe2=CreatePVectorSheet(element);
723 ElementVector* pe3=CreatePVectorShelf(element);
724 ElementVector* pe =new ElementVector(pe1,pe2,pe3);
725
726 /*clean-up and return*/
727 delete pe1;
728 delete pe2;
729 delete pe3;
730 return pe;
731}/*}}}*/
732ElementVector* EnthalpyAnalysis::CreatePVectorVolume(Element* element){/*{{{*/
733
734 /* Check if ice in element */
735 if(!element->IsIceInElement()) return NULL;
736
737 /*Intermediaries*/
738 int i, stabilization;
739 IssmDouble Jdet,phi,dt;
740 IssmDouble enthalpy, Hpmp;
741 IssmDouble enthalpypicard, d1enthalpypicard[3];
742 IssmDouble pressure, d1pressure[3], d2pressure;
743 IssmDouble waterfractionpicard;
744 IssmDouble kappa,tau_parameter,diameter,kappa_w;
745 IssmDouble u,v,w;
746 IssmDouble scalar_def, scalar_sens ,scalar_transient;
747 IssmDouble* xyz_list = NULL;
748 IssmDouble d1H_d1P, d1P2;
749
750 /*Fetch number of nodes and dof for this finite element*/
751 int numnodes = element->GetNumberOfNodes();
752 int numvertices = element->GetNumberOfVertices();
753
754 /*Initialize Element vector*/
755 ElementVector* pe = element->NewElementVector();
756 IssmDouble* basis = xNew<IssmDouble>(numnodes);
757 IssmDouble* dbasis = xNew<IssmDouble>(3*numnodes);
758
759 /*Retrieve all inputs and parameters*/
760 element->GetVerticesCoordinates(&xyz_list);
761 IssmDouble rho_ice = element->GetMaterialParameter(MaterialsRhoIceEnum);
762 IssmDouble heatcapacity = element->GetMaterialParameter(MaterialsHeatcapacityEnum);
763 IssmDouble thermalconductivity = element->GetMaterialParameter(MaterialsThermalconductivityEnum);
764 IssmDouble temperateiceconductivity = element->GetMaterialParameter(MaterialsTemperateiceconductivityEnum);
765 IssmDouble beta = element->GetMaterialParameter(MaterialsBetaEnum);
766 IssmDouble latentheat = element->GetMaterialParameter(MaterialsLatentheatEnum);
767 element->FindParam(&dt,TimesteppingTimeStepEnum);
768 element->FindParam(&stabilization,ThermalStabilizationEnum);
769 Input* vx_input=element->GetInput(VxEnum); _assert_(vx_input);
770 Input* vy_input=element->GetInput(VyEnum); _assert_(vy_input);
771 Input* vz_input=element->GetInput(VzEnum); _assert_(vz_input);
772 Input* enthalpypicard_input=element->GetInput(EnthalpyPicardEnum); _assert_(enthalpypicard_input);
773 Input* pressure_input=element->GetInput(PressureEnum); _assert_(pressure_input);
774 Input* enthalpy_input=NULL;
775 if(reCast<bool,IssmDouble>(dt)){enthalpy_input = element->GetInput(EnthalpyEnum); _assert_(enthalpy_input);}
776 if(stabilization==2){
777 diameter=element->MinEdgeLength(xyz_list);
778 kappa=this->EnthalpyDiffusionParameterVolume(element,EnthalpyPicardEnum); _assert_(kappa>=0.);
779 }
780
781 /* Start looping on the number of gaussian points: */
782 Gauss* gauss=element->NewGauss(4);
783 for(int ig=gauss->begin();ig<gauss->end();ig++){
784 gauss->GaussPoint(ig);
785
786 element->JacobianDeterminant(&Jdet,xyz_list,gauss);
787 element->NodalFunctions(basis,gauss);
788
789 /*viscous dissipation*/
790 element->ViscousHeating(&phi,xyz_list,gauss,vx_input,vy_input,vz_input);
791
792 scalar_def=phi/rho_ice*Jdet*gauss->weight;
793 if(dt!=0.) scalar_def=scalar_def*dt;
794
795 for(i=0;i<numnodes;i++) pe->values[i]+=scalar_def*basis[i];
796
797 /*sensible heat flux in temperate ice*/
798 enthalpypicard_input->GetInputValue(&enthalpypicard,gauss);
799 pressure_input->GetInputValue(&pressure,gauss);
800 Hpmp=this->PureIceEnthalpy(element, pressure);
801
802 if(enthalpypicard>=Hpmp){
803 enthalpypicard_input->GetInputDerivativeValue(&d1enthalpypicard[0],xyz_list,gauss);
804 pressure_input->GetInputDerivativeValue(&d1pressure[0],xyz_list,gauss);
805 d2pressure=0.; // for linear elements, 2nd derivative is zero
806
807 d1H_d1P=0.;
808 for(i=0;i<3;i++) d1H_d1P+=d1enthalpypicard[i]*d1pressure[i];
809 d1P2=0.;
810 for(i=0;i<3;i++) d1P2+=pow(d1pressure[i],2.);
811
812 scalar_sens=-beta*((temperateiceconductivity - thermalconductivity)/latentheat*(d1H_d1P + beta*heatcapacity*d1P2))/rho_ice;
813 if(dt!=0.) scalar_sens=scalar_sens*dt;
814 for(i=0;i<numnodes;i++) pe->values[i]+=scalar_sens*basis[i];
815 }
816
817 /* Build transient now */
818 if(reCast<bool,IssmDouble>(dt)){
819 enthalpy_input->GetInputValue(&enthalpy, gauss);
820 scalar_transient=enthalpy*Jdet*gauss->weight;
821 for(i=0;i<numnodes;i++) pe->values[i]+=scalar_transient*basis[i];
822 }
823
824 if(stabilization==2){
825 element->NodalFunctionsDerivatives(dbasis,xyz_list,gauss);
826
827 vx_input->GetInputValue(&u,gauss);
828 vy_input->GetInputValue(&v,gauss);
829 vz_input->GetInputValue(&w,gauss);
830 tau_parameter=element->StabilizationParameter(u,v,w,diameter,kappa/rho_ice);
831
832 for(i=0;i<numnodes;i++) pe->values[i]+=tau_parameter*scalar_def*(u*dbasis[0*numnodes+i]+v*dbasis[1*numnodes+i]+w*dbasis[2*numnodes+i]);
833
834 if(dt!=0.){
835 for(i=0;i<numnodes;i++) pe->values[i]+=tau_parameter*scalar_transient*(u*dbasis[0*numnodes+i]+v*dbasis[1*numnodes+i]+w*dbasis[2*numnodes+i]);
836 }
837 }
838 }
839
840 /*Clean up and return*/
841 xDelete<IssmDouble>(basis);
842 xDelete<IssmDouble>(dbasis);
843 xDelete<IssmDouble>(xyz_list);
844 delete gauss;
845 return pe;
846
847}/*}}}*/
848ElementVector* EnthalpyAnalysis::CreatePVectorSheet(Element* element){/*{{{*/
849
850 /* Check if ice in element */
851 if(!element->IsIceInElement()) return NULL;
852
853 /* implementation of the basal condition decision chart of Aschwanden 2012, Fig.5 */
854 if(!element->IsOnBase() || element->IsFloating()) return NULL;
855
856 bool converged, isdynamicbasalspc;
857 int i, state;
858 int enthalpy_enum;
859 IssmDouble dt,Jdet,scalar;
860 IssmDouble enthalpy, enthalpyup, pressure, pressureup, watercolumn, meltingrate;
861 IssmDouble vx,vy,vz;
862 IssmDouble alpha2,basalfriction,geothermalflux,heatflux;
863 IssmDouble *xyz_list_base = NULL;
864
865 /*Fetch number of nodes for this finite element*/
866 int numnodes = element->GetNumberOfNodes();
867
868 /*Initialize vectors*/
869 ElementVector* pe = element->NewElementVector();
870 IssmDouble* basis = xNew<IssmDouble>(numnodes);
871
872 /*Retrieve all inputs and parameters*/
873 element->GetVerticesCoordinatesBase(&xyz_list_base);
874 element->FindParam(&dt,TimesteppingTimeStepEnum);
875 element->FindParam(&isdynamicbasalspc,ThermalIsdynamicbasalspcEnum);
876 element->GetInputValue(&converged,ConvergedEnum);
877 if(dt==0. && !converged) enthalpy_enum=EnthalpyPicardEnum; // use enthalpy from last iteration
878 else enthalpy_enum=EnthalpyEnum; // use enthalpy from last time step
879 Input* vx_input = element->GetInput(VxEnum); _assert_(vx_input);
880 Input* vy_input = element->GetInput(VyEnum); _assert_(vy_input);
881 Input* vz_input = element->GetInput(VzEnum); _assert_(vz_input);
882 Input* enthalpy_input = element->GetInput(enthalpy_enum); _assert_(enthalpy_input);
883 Input* pressure_input = element->GetInput(PressureEnum); _assert_(pressure_input);
884 Input* watercolumn_input = element->GetInput(WatercolumnEnum); _assert_(watercolumn_input);
885 Input* meltingrate_input = element->GetInput(BasalforcingsGroundediceMeltingRateEnum); _assert_(meltingrate_input);
886 Input* geothermalflux_input = element->GetInput(BasalforcingsGeothermalfluxEnum); _assert_(geothermalflux_input);
887 IssmDouble rho_ice = element->GetMaterialParameter(MaterialsRhoIceEnum);
888
889 /*Build friction element, needed later: */
890 Friction* friction=new Friction(element,3);
891
892 /* Start looping on the number of gaussian points: */
893 Gauss* gauss=element->NewGaussBase(4);
894 Gauss* gaussup=element->NewGaussTop(4);
895 for(int ig=gauss->begin();ig<gauss->end();ig++){
896 gauss->GaussPoint(ig);
897 gaussup->GaussPoint(ig);
898
899 element->JacobianDeterminantBase(&Jdet,xyz_list_base,gauss);
900 element->NodalFunctions(basis,gauss);
901
902 if(isdynamicbasalspc){
903 enthalpy_input->GetInputValue(&enthalpy,gauss);
904 enthalpy_input->GetInputValue(&enthalpyup,gaussup);
905 pressure_input->GetInputValue(&pressure,gauss);
906 pressure_input->GetInputValue(&pressureup,gaussup);
907 watercolumn_input->GetInputValue(&watercolumn,gauss);
908 meltingrate_input->GetInputValue(&meltingrate,gauss);
909 state=GetThermalBasalCondition(element, enthalpy, enthalpyup, pressure, pressureup, watercolumn, meltingrate);
910 }
911 else
912 state=0;
913
914 switch (state) {
915 case 0: case 1: case 2: case 3:
916 // cold, dry base; cold, wet base; refreezing temperate base; thin temperate base:
917 // Apply basal surface forcing.
918 // Interpolated values of enthalpy on gauss nodes may indicate cold base,
919 // although one node might have become temperate. So keep heat flux switched on.
920 geothermalflux_input->GetInputValue(&geothermalflux,gauss);
921 friction->GetAlpha2(&alpha2,gauss);
922 vx_input->GetInputValue(&vx,gauss);
923 vy_input->GetInputValue(&vy,gauss);
924 vz_input->GetInputValue(&vz,gauss);
925 basalfriction=alpha2*(vx*vx+vy*vy+vz*vz);
926 heatflux=(basalfriction+geothermalflux)/(rho_ice);
927 scalar=gauss->weight*Jdet*heatflux;
928 if(dt!=0.) scalar=dt*scalar;
929 for(i=0;i<numnodes;i++)
930 pe->values[i]+=scalar*basis[i];
931 break;
932 case 4:
933 // temperate, thick melting base: set grad H*n=0
934 for(i=0;i<numnodes;i++)
935 pe->values[i]+=0.;
936 break;
937 default:
938 _printf0_(" unknown thermal basal state found!");
939 }
940 }
941
942 /*Clean up and return*/
943 delete gauss;
944 delete gaussup;
945 delete friction;
946 xDelete<IssmDouble>(basis);
947 xDelete<IssmDouble>(xyz_list_base);
948 return pe;
949
950}/*}}}*/
951ElementVector* EnthalpyAnalysis::CreatePVectorShelf(Element* element){/*{{{*/
952
953 /* Check if ice in element */
954 if(!element->IsIceInElement()) return NULL;
955
956 /*Get basal element*/
957 if(!element->IsOnBase() || !element->IsFloating()) return NULL;
958
959 IssmDouble Hpmp,dt,Jdet,scalar_ocean,pressure;
960 IssmDouble *xyz_list_base = NULL;
961
962 /*Fetch number of nodes for this finite element*/
963 int numnodes = element->GetNumberOfNodes();
964
965 /*Initialize vectors*/
966 ElementVector* pe = element->NewElementVector();
967 IssmDouble* basis = xNew<IssmDouble>(numnodes);
968
969 /*Retrieve all inputs and parameters*/
970 element->GetVerticesCoordinatesBase(&xyz_list_base);
971 element->FindParam(&dt,TimesteppingTimeStepEnum);
972 Input* pressure_input=element->GetInput(PressureEnum); _assert_(pressure_input);
973 IssmDouble gravity = element->GetMaterialParameter(ConstantsGEnum);
974 IssmDouble rho_water = element->GetMaterialParameter(MaterialsRhoSeawaterEnum);
975 IssmDouble rho_ice = element->GetMaterialParameter(MaterialsRhoIceEnum);
976 IssmDouble heatcapacity = element->GetMaterialParameter(MaterialsHeatcapacityEnum);
977 IssmDouble mixed_layer_capacity= element->GetMaterialParameter(MaterialsMixedLayerCapacityEnum);
978 IssmDouble thermal_exchange_vel= element->GetMaterialParameter(MaterialsThermalExchangeVelocityEnum);
979
980 /* Start looping on the number of gaussian points: */
981 Gauss* gauss=element->NewGaussBase(4);
982 for(int ig=gauss->begin();ig<gauss->end();ig++){
983 gauss->GaussPoint(ig);
984
985 element->JacobianDeterminantBase(&Jdet,xyz_list_base,gauss);
986 element->NodalFunctions(basis,gauss);
987
988 pressure_input->GetInputValue(&pressure,gauss);
989 Hpmp=element->PureIceEnthalpy(pressure);
990
991 scalar_ocean=gauss->weight*Jdet*rho_water*mixed_layer_capacity*thermal_exchange_vel*Hpmp/(heatcapacity*rho_ice);
992 if(reCast<bool,IssmDouble>(dt)) scalar_ocean=dt*scalar_ocean;
993
994 for(int i=0;i<numnodes;i++) pe->values[i]+=scalar_ocean*basis[i];
995 }
996
997 /*Clean up and return*/
998 delete gauss;
999 xDelete<IssmDouble>(basis);
1000 xDelete<IssmDouble>(xyz_list_base);
1001 return pe;
1002}/*}}}*/
1003void EnthalpyAnalysis::DrainWaterfraction(FemModel* femmodel){/*{{{*/
1004 /*Drain excess water fraction in ice column: */
1005 for(int i=0;i<femmodel->elements->Size();i++){
1006 Element* element=xDynamicCast<Element*>(femmodel->elements->GetObjectByOffset(i));
1007 DrainWaterfractionIcecolumn(element);
1008 }
1009}/*}}}*/
1010void EnthalpyAnalysis::DrainWaterfraction(Element* element, IssmDouble* pdrainrate_element){/*{{{*/
1011
1012 /* Check if ice in element */
1013 if(!element->IsIceInElement()) return;
1014
1015 /*Intermediaries*/
1016 int iv,is,vertexdown,vertexup,numsegments;
1017 IssmDouble dt, height_element;
1018 IssmDouble rho_water, rho_ice;
1019 int numvertices = element->GetNumberOfVertices();
1020
1021 IssmDouble* xyz_list = NULL;
1022 IssmDouble* enthalpies = xNew<IssmDouble>(numvertices);
1023 IssmDouble* pressures = xNew<IssmDouble>(numvertices);
1024 IssmDouble* temperatures = xNew<IssmDouble>(numvertices);
1025 IssmDouble* waterfractions = xNew<IssmDouble>(numvertices);
1026 IssmDouble* deltawaterfractions = xNew<IssmDouble>(numvertices);
1027 int *pairindices = NULL;
1028
1029 rho_ice=element->GetMaterialParameter(MaterialsRhoIceEnum);
1030 rho_water=element->GetMaterialParameter(MaterialsRhoSeawaterEnum);
1031
1032 element->GetVerticesCoordinates(&xyz_list);
1033 element->GetInputListOnVertices(enthalpies,EnthalpyEnum);
1034 element->GetInputListOnVertices(pressures,PressureEnum);
1035
1036 element->FindParam(&dt,TimesteppingTimeStepEnum);
1037 for(iv=0;iv<numvertices;iv++){
1038 element->EnthalpyToThermal(&temperatures[iv],&waterfractions[iv], enthalpies[iv],pressures[iv]);
1039 deltawaterfractions[iv]=DrainageFunctionWaterfraction(waterfractions[iv], dt);
1040 }
1041
1042 /*drain waterfraction, feed updated variables back into model*/
1043 for(iv=0;iv<numvertices;iv++){
1044 if(reCast<bool,IssmDouble>(dt))
1045 waterfractions[iv]-=deltawaterfractions[iv]*dt;
1046 else
1047 waterfractions[iv]-=deltawaterfractions[iv];
1048 element->ThermalToEnthalpy(&enthalpies[iv], temperatures[iv], waterfractions[iv], pressures[iv]);
1049 }
1050 element->AddInput(EnthalpyEnum,enthalpies,P1Enum);
1051 element->AddInput(WaterfractionEnum,waterfractions,P1Enum);
1052
1053 /*return meltwater column equivalent to drained water*/
1054 element->VerticalSegmentIndices(&pairindices,&numsegments);
1055 for(is=0;is<numsegments;is++){
1056 vertexdown = pairindices[is*2+0];
1057 vertexup = pairindices[is*2+1];
1058 height_element=fabs(xyz_list[vertexup*3+2]-xyz_list[vertexdown*3+2]);
1059 pdrainrate_element[is]=(deltawaterfractions[vertexdown]+deltawaterfractions[vertexup])/2.*height_element; // return water equivalent of drainage
1060 _assert_(pdrainrate_element[is]>=0.);
1061 }
1062
1063 /*Clean up and return*/
1064 xDelete<int>(pairindices);
1065 xDelete<IssmDouble>(xyz_list);
1066 xDelete<IssmDouble>(enthalpies);
1067 xDelete<IssmDouble>(pressures);
1068 xDelete<IssmDouble>(temperatures);
1069 xDelete<IssmDouble>(waterfractions);
1070 xDelete<IssmDouble>(deltawaterfractions);
1071}/*}}}*/
1072void EnthalpyAnalysis::DrainWaterfractionIcecolumn(Element* element){/*{{{*/
1073
1074 /* Check if ice in element */
1075 if(!element->IsIceInElement()) return;
1076
1077 /* Only drain waterfraction of ice column from element at base*/
1078 if(!element->IsOnBase()) return; //FIXME: allow freeze on for floating elements
1079
1080 /* Intermediaries*/
1081 int is, numvertices, numsegments;
1082 int *pairindices = NULL;
1083
1084 numvertices=element->GetNumberOfVertices();
1085 element->VerticalSegmentIndices(&pairindices,&numsegments);
1086
1087 IssmDouble* watercolumn = xNew<IssmDouble>(numvertices);
1088 IssmDouble* drainrate_column = xNew<IssmDouble>(numsegments);
1089 IssmDouble* drainrate_element = xNew<IssmDouble>(numsegments);
1090
1091 element->GetInputListOnVertices(watercolumn,WatercolumnEnum);
1092
1093 for(is=0;is<numsegments;is++) drainrate_column[is]=0.;
1094 Element* elementi = element;
1095 for(;;){
1096 for(is=0;is<numsegments;is++) drainrate_element[is]=0.;
1097 DrainWaterfraction(elementi,drainrate_element); // TODO: make sure every vertex is only drained once
1098 for(is=0;is<numsegments;is++) drainrate_column[is]+=drainrate_element[is];
1099
1100 if(elementi->IsOnSurface()) break;
1101 elementi=elementi->GetUpperElement();
1102 }
1103 /* add drained water to water column*/
1104 for(is=0;is<numsegments;is++) watercolumn[is]+=drainrate_column[is];
1105 /* Feed updated water column back into model */
1106 element->AddInput(WatercolumnEnum,watercolumn,P1Enum);
1107
1108 xDelete<int>(pairindices);
1109 xDelete<IssmDouble>(drainrate_column);
1110 xDelete<IssmDouble>(drainrate_element);
1111 xDelete<IssmDouble>(watercolumn);
1112}/*}}}*/
1113IssmDouble EnthalpyAnalysis::EnthalpyDiffusionParameter(Element* element,IssmDouble enthalpy,IssmDouble pressure){/*{{{*/
1114
1115 IssmDouble heatcapacity = element->GetMaterialParameter(MaterialsHeatcapacityEnum);
1116 IssmDouble temperateiceconductivity = element->GetMaterialParameter(MaterialsTemperateiceconductivityEnum);
1117 IssmDouble thermalconductivity = element->GetMaterialParameter(MaterialsThermalconductivityEnum);
1118
1119 if(enthalpy < PureIceEnthalpy(element,pressure)){
1120 return thermalconductivity/heatcapacity;
1121 }
1122 else{
1123 return temperateiceconductivity/heatcapacity;
1124 }
1125}/*}}}*/
1126IssmDouble EnthalpyAnalysis::EnthalpyDiffusionParameterVolume(Element* element,int enthalpy_enum){/*{{{*/
1127
1128 int iv;
1129 IssmDouble lambda; /* fraction of cold ice */
1130 IssmDouble kappa,kappa_c,kappa_t; /* enthalpy conductivities */
1131 IssmDouble Hc,Ht;
1132
1133 /*Get pressures and enthalpies on vertices*/
1134 int numvertices = element->GetNumberOfVertices();
1135 IssmDouble* pressures = xNew<IssmDouble>(numvertices);
1136 IssmDouble* enthalpies = xNew<IssmDouble>(numvertices);
1137 IssmDouble* PIE = xNew<IssmDouble>(numvertices);
1138 IssmDouble* dHpmp = xNew<IssmDouble>(numvertices);
1139 element->GetInputListOnVertices(pressures,PressureEnum);
1140 element->GetInputListOnVertices(enthalpies,enthalpy_enum);
1141 for(iv=0;iv<numvertices;iv++){
1142 PIE[iv] = PureIceEnthalpy(element,pressures[iv]);
1143 dHpmp[iv] = enthalpies[iv]-PIE[iv];
1144 }
1145
1146 bool allequalsign = true;
1147 if(dHpmp[0]<0.){
1148 for(iv=1; iv<numvertices;iv++) allequalsign=(allequalsign && (dHpmp[iv]<0.));
1149 }
1150 else{
1151 for(iv=1; iv<numvertices;iv++) allequalsign=(allequalsign && (dHpmp[iv]>=0.));
1152 }
1153
1154 if(allequalsign){
1155 kappa = EnthalpyDiffusionParameter(element,enthalpies[0],pressures[0]);
1156 }
1157 else{
1158 /* return harmonic mean of thermal conductivities, weighted by fraction of cold/temperate ice,
1159 cf Patankar 1980, pp44 */
1160 kappa_c = EnthalpyDiffusionParameter(element,PureIceEnthalpy(element,0.)-1.,0.);
1161 kappa_t = EnthalpyDiffusionParameter(element,PureIceEnthalpy(element,0.)+1.,0.);
1162 Hc=0.; Ht=0.;
1163 for(iv=0; iv<numvertices;iv++){
1164 if(enthalpies[iv]<PIE[iv])
1165 Hc+=(PIE[iv]-enthalpies[iv]);
1166 else
1167 Ht+=(enthalpies[iv]-PIE[iv]);
1168 }
1169 _assert_((Hc+Ht)>0.);
1170 lambda = Hc/(Hc+Ht);
1171 kappa = kappa_c*kappa_t/(lambda*kappa_t+(1.-lambda)*kappa_c); // ==(lambda/kappa_c + (1.-lambda)/kappa_t)^-1
1172 }
1173
1174 /*Clean up and return*/
1175 xDelete<IssmDouble>(PIE);
1176 xDelete<IssmDouble>(dHpmp);
1177 xDelete<IssmDouble>(pressures);
1178 xDelete<IssmDouble>(enthalpies);
1179 return kappa;
1180}/*}}}*/
1181void EnthalpyAnalysis::GetBAdvec(IssmDouble* B,Element* element,IssmDouble* xyz_list,Gauss* gauss){/*{{{*/
1182 /*Compute B matrix. B=[B1 B2 B3 B4 B5 B6] where Bi is of size 5*NDOF1.
1183 * For node i, Bi' can be expressed in the actual coordinate system
1184 * by:
1185 * Bi_advec =[ h ]
1186 * [ h ]
1187 * [ h ]
1188 * where h is the interpolation function for node i.
1189 *
1190 * We assume B has been allocated already, of size: 3x(NDOF1*NUMNODESP1)
1191 */
1192
1193 /*Fetch number of nodes for this finite element*/
1194 int numnodes = element->GetNumberOfNodes();
1195
1196 /*Get nodal functions*/
1197 IssmDouble* basis=xNew<IssmDouble>(numnodes);
1198 element->NodalFunctions(basis,gauss);
1199
1200 /*Build B: */
1201 for(int i=0;i<numnodes;i++){
1202 B[numnodes*0+i] = basis[i];
1203 B[numnodes*1+i] = basis[i];
1204 B[numnodes*2+i] = basis[i];
1205 }
1206
1207 /*Clean-up*/
1208 xDelete<IssmDouble>(basis);
1209}/*}}}*/
1210void EnthalpyAnalysis::GetBAdvecprime(IssmDouble* B,Element* element,IssmDouble* xyz_list,Gauss* gauss){/*{{{*/
1211 /*Compute B matrix. B=[B1 B2 B3 B4 B5 B6] where Bi is of size 5*NDOF1.
1212 * For node i, Bi' can be expressed in the actual coordinate system
1213 * by:
1214 * Biprime_advec=[ dh/dx ]
1215 * [ dh/dy ]
1216 * [ dh/dz ]
1217 * where h is the interpolation function for node i.
1218 *
1219 * We assume B has been allocated already, of size: 3x(NDOF1*numnodes)
1220 */
1221
1222 /*Fetch number of nodes for this finite element*/
1223 int numnodes = element->GetNumberOfNodes();
1224
1225 /*Get nodal functions derivatives*/
1226 IssmDouble* dbasis=xNew<IssmDouble>(3*numnodes);
1227 element->NodalFunctionsDerivatives(dbasis,xyz_list,gauss);
1228
1229 /*Build B: */
1230 for(int i=0;i<numnodes;i++){
1231 B[numnodes*0+i] = dbasis[0*numnodes+i];
1232 B[numnodes*1+i] = dbasis[1*numnodes+i];
1233 B[numnodes*2+i] = dbasis[2*numnodes+i];
1234 }
1235
1236 /*Clean-up*/
1237 xDelete<IssmDouble>(dbasis);
1238}/*}}}*/
1239void EnthalpyAnalysis::GetBasalConstraints(Vector<IssmDouble>* vec_spc,Element* element){/*{{{*/
1240
1241 /*Intermediary*/
1242 bool isdynamicbasalspc;
1243 IssmDouble dt;
1244
1245 /*Check wether dynamic basal boundary conditions are activated */
1246 element->FindParam(&isdynamicbasalspc,ThermalIsdynamicbasalspcEnum);
1247 if(!isdynamicbasalspc) return;
1248
1249 element->FindParam(&dt,TimesteppingTimeStepEnum);
1250 if(dt==0.){
1251 GetBasalConstraintsSteadystate(vec_spc,element);
1252 }
1253 else{
1254 GetBasalConstraintsTransient(vec_spc,element);
1255 }
1256}/*}}}*/
1257void EnthalpyAnalysis::GetBasalConstraintsSteadystate(Vector<IssmDouble>* vec_spc,Element* element){/*{{{*/
1258
1259 /* Check if ice in element */
1260 if(!element->IsIceInElement()) return;
1261
1262 /* Only update constraints at the base.
1263 * Floating ice is not affected by basal BC decision chart. */
1264 if(!(element->IsOnBase()) || element->IsFloating()) return;
1265
1266 /*Intermediary*/
1267 int numindices, numindicesup, state;
1268 int *indices = NULL, *indicesup = NULL;
1269 IssmDouble enthalpy, enthalpyup, pressure, pressureup, watercolumn, meltingrate;
1270
1271 /*Get parameters and inputs: */
1272 Input* enthalpy_input = element->GetInput(EnthalpyPicardEnum); _assert_(enthalpy_input);
1273 Input* pressure_input = element->GetInput(PressureEnum); _assert_(pressure_input);
1274 Input* watercolumn_input = element->GetInput(WatercolumnEnum); _assert_(watercolumn_input);
1275 Input* meltingrate_input = element->GetInput(BasalforcingsGroundediceMeltingRateEnum); _assert_(meltingrate_input);
1276
1277 /*Fetch indices of basal & surface nodes for this finite element*/
1278 Penta *penta = (Penta *) element; // TODO: add Basal-/SurfaceNodeIndices to element.h, and change this to Element*
1279 penta->BasalNodeIndices(&numindices,&indices,element->GetElementType());
1280 penta->SurfaceNodeIndices(&numindicesup,&indicesup,element->GetElementType()); _assert_(numindices==numindicesup);
1281
1282 GaussPenta* gauss=new GaussPenta();
1283 GaussPenta* gaussup=new GaussPenta();
1284 for(int i=0;i<numindices;i++){
1285 gauss->GaussNode(element->GetElementType(),indices[i]);
1286 gaussup->GaussNode(element->GetElementType(),indicesup[i]);
1287
1288 enthalpy_input->GetInputValue(&enthalpy,gauss);
1289 enthalpy_input->GetInputValue(&enthalpyup,gaussup);
1290 pressure_input->GetInputValue(&pressure,gauss);
1291 pressure_input->GetInputValue(&pressureup,gaussup);
1292 watercolumn_input->GetInputValue(&watercolumn,gauss);
1293 meltingrate_input->GetInputValue(&meltingrate,gauss);
1294
1295 state=GetThermalBasalCondition(element, enthalpy, enthalpyup, pressure, pressureup, watercolumn, meltingrate);
1296 switch (state) {
1297 case 0:
1298 // cold, dry base: apply basal surface forcing
1299 vec_spc->SetValue(element->nodes[i]->Sid(),0.,INS_VAL);
1300 break;
1301 case 1:
1302 // cold, wet base: keep at pressure melting point
1303 vec_spc->SetValue(element->nodes[i]->Sid(),1.,INS_VAL);
1304 break;
1305 case 2:
1306 // temperate, thin refreezing base:
1307 vec_spc->SetValue(element->nodes[i]->Sid(),1.,INS_VAL);
1308 break;
1309 case 3:
1310 // temperate, thin melting base: set spc
1311 vec_spc->SetValue(element->nodes[i]->Sid(),1.,INS_VAL);
1312 break;
1313 case 4:
1314 // temperate, thick melting base:
1315 vec_spc->SetValue(element->nodes[i]->Sid(),1.,INS_VAL);
1316 break;
1317 default:
1318 _printf0_(" unknown thermal basal state found!");
1319 }
1320 }
1321
1322 /*Free ressources:*/
1323 xDelete<int>(indices);
1324 xDelete<int>(indicesup);
1325 delete gauss;
1326 delete gaussup;
1327}/*}}}*/
1328void EnthalpyAnalysis::GetBasalConstraintsTransient(Vector<IssmDouble>* vec_spc,Element* element){/*{{{*/
1329
1330 /* Check if ice in element */
1331 if(!element->IsIceInElement()) return;
1332
1333 /* Only update constraints at the base.
1334 * Floating ice is not affected by basal BC decision chart.*/
1335 if(!(element->IsOnBase()) || element->IsFloating()) return;
1336
1337 /*Intermediary*/
1338 int numindices, numindicesup, state;
1339 int *indices = NULL, *indicesup = NULL;
1340 IssmDouble enthalpy, enthalpyup, pressure, pressureup, watercolumn, meltingrate;
1341
1342 /*Get parameters and inputs: */
1343 Input* enthalpy_input = element->GetInput(EnthalpyEnum); _assert_(enthalpy_input); //TODO: check EnthalpyPicard?
1344 Input* pressure_input = element->GetInput(PressureEnum); _assert_(pressure_input);
1345 Input* watercolumn_input = element->GetInput(WatercolumnEnum); _assert_(watercolumn_input);
1346 Input* meltingrate_input = element->GetInput(BasalforcingsGroundediceMeltingRateEnum); _assert_(meltingrate_input);
1347
1348 /*Fetch indices of basal & surface nodes for this finite element*/
1349 Penta *penta = (Penta *) element; // TODO: add Basal-/SurfaceNodeIndices to element.h, and change this to Element*
1350 penta->BasalNodeIndices(&numindices,&indices,element->GetElementType());
1351 penta->SurfaceNodeIndices(&numindicesup,&indicesup,element->GetElementType()); _assert_(numindices==numindicesup);
1352
1353 GaussPenta* gauss=new GaussPenta();
1354 GaussPenta* gaussup=new GaussPenta();
1355
1356 for(int i=0;i<numindices;i++){
1357 gauss->GaussNode(element->GetElementType(),indices[i]);
1358 gaussup->GaussNode(element->GetElementType(),indicesup[i]);
1359
1360 enthalpy_input->GetInputValue(&enthalpy,gauss);
1361 enthalpy_input->GetInputValue(&enthalpyup,gaussup);
1362 pressure_input->GetInputValue(&pressure,gauss);
1363 pressure_input->GetInputValue(&pressureup,gaussup);
1364 watercolumn_input->GetInputValue(&watercolumn,gauss);
1365 meltingrate_input->GetInputValue(&meltingrate,gauss);
1366
1367 state=GetThermalBasalCondition(element, enthalpy, enthalpyup, pressure, pressureup, watercolumn, meltingrate);
1368
1369 switch (state) {
1370 case 0:
1371 // cold, dry base: apply basal surface forcing
1372 vec_spc->SetValue(element->nodes[i]->Sid(),0.,INS_VAL);
1373 break;
1374 case 1:
1375 // cold, wet base: keep at pressure melting point
1376 vec_spc->SetValue(element->nodes[i]->Sid(),1.,INS_VAL);
1377 break;
1378 case 2:
1379 // temperate, thin refreezing base: release spc
1380 vec_spc->SetValue(element->nodes[i]->Sid(),0.,INS_VAL);
1381 break;
1382 case 3:
1383 // temperate, thin melting base: set spc
1384 vec_spc->SetValue(element->nodes[i]->Sid(),1.,INS_VAL);
1385 break;
1386 case 4:
1387 // temperate, thick melting base: set grad H*n=0
1388 vec_spc->SetValue(element->nodes[i]->Sid(),0.,INS_VAL);
1389 break;
1390 default:
1391 _printf0_(" unknown thermal basal state found!");
1392 }
1393
1394 }
1395
1396 /*Free ressources:*/
1397 xDelete<int>(indices);
1398 xDelete<int>(indicesup);
1399 delete gauss;
1400 delete gaussup;
1401}/*}}}*/
1402void EnthalpyAnalysis::GetBConduct(IssmDouble* B,Element* element,IssmDouble* xyz_list,Gauss* gauss){/*{{{*/
1403 /*Compute B matrix. B=[B1 B2 B3 B4 B5 B6] where Bi is of size 5*NDOF1.
1404 * For node i, Bi' can be expressed in the actual coordinate system
1405 * by:
1406 * Bi_conduct=[ dh/dx ]
1407 * [ dh/dy ]
1408 * [ dh/dz ]
1409 * where h is the interpolation function for node i.
1410 *
1411 * We assume B has been allocated already, of size: 3x(NDOF1*numnodes)
1412 */
1413
1414 /*Fetch number of nodes for this finite element*/
1415 int numnodes = element->GetNumberOfNodes();
1416
1417 /*Get nodal functions derivatives*/
1418 IssmDouble* dbasis=xNew<IssmDouble>(3*numnodes);
1419 element->NodalFunctionsDerivatives(dbasis,xyz_list,gauss);
1420
1421 /*Build B: */
1422 for(int i=0;i<numnodes;i++){
1423 B[numnodes*0+i] = dbasis[0*numnodes+i];
1424 B[numnodes*1+i] = dbasis[1*numnodes+i];
1425 B[numnodes*2+i] = dbasis[2*numnodes+i];
1426 }
1427
1428 /*Clean-up*/
1429 xDelete<IssmDouble>(dbasis);
1430}/*}}}*/
1431void EnthalpyAnalysis::GetSolutionFromInputs(Vector<IssmDouble>* solution,Element* element){/*{{{*/
1432 element->GetSolutionFromInputsOneDof(solution,EnthalpyEnum);
1433}/*}}}*/
1434int EnthalpyAnalysis::GetThermalBasalCondition(Element* element, IssmDouble enthalpy, IssmDouble enthalpyup, IssmDouble pressure, IssmDouble pressureup, IssmDouble watercolumn, IssmDouble meltingrate){/*{{{*/
1435
1436 /* Check if ice in element */
1437 if(!element->IsIceInElement()) return -1;
1438
1439 /* Only update Constraints at the base of grounded ice*/
1440 if(!(element->IsOnBase())) return -1;
1441
1442 /*Intermediary*/
1443 int state=-1;
1444 IssmDouble dt;
1445
1446 /*Get parameters and inputs: */
1447 element->FindParam(&dt,TimesteppingTimeStepEnum);
1448
1449 if(enthalpy<PureIceEnthalpy(element,pressure)){
1450 if(watercolumn<=0.) state=0; // cold, dry base
1451 else state=1; // cold, wet base (refreezing)
1452 }
1453 else{
1454 if(enthalpyup<PureIceEnthalpy(element,pressureup)){
1455 if((dt==0.) && (meltingrate<0.)) state=2; // refreezing temperate base (non-physical, only for steadystate solver)
1456 else state=3; // temperate base, but no temperate layer
1457 }
1458 else state=4; // temperate layer with positive thickness
1459 }
1460
1461 _assert_(state>=0);
1462 return state;
1463}/*}}}*/
1464IssmDouble EnthalpyAnalysis::GetWetIceConductivity(Element* element, IssmDouble enthalpy, IssmDouble pressure){/*{{{*/
1465
1466 IssmDouble temperature, waterfraction;
1467 IssmDouble kappa_w = 0.6; // thermal conductivity of water (in W/m/K)
1468 IssmDouble kappa_i = element->GetMaterialParameter(MaterialsThermalconductivityEnum);
1469 element->EnthalpyToThermal(&temperature, &waterfraction, enthalpy, pressure);
1470
1471 return (1.-waterfraction)*kappa_i + waterfraction*kappa_w;
1472}/*}}}*/
1473void EnthalpyAnalysis::GradientJ(Vector<IssmDouble>* gradient,Element* element,int control_type,int control_index){/*{{{*/
1474 _error_("Not implemented yet");
1475}/*}}}*/
1476void EnthalpyAnalysis::InputUpdateFromSolution(IssmDouble* solution,Element* element){/*{{{*/
1477
1478 bool converged;
1479 int i,rheology_law;
1480 IssmDouble B_average,s_average,T_average=0.,P_average=0.;
1481 IssmDouble n=3.0;
1482 int *doflist = NULL;
1483 IssmDouble *xyz_list = NULL;
1484
1485 /*Fetch number of nodes and dof for this finite element*/
1486 int numnodes = element->GetNumberOfNodes();
1487
1488 /*Fetch dof list and allocate solution vector*/
1489 element->GetDofList(&doflist,NoneApproximationEnum,GsetEnum);
1490 IssmDouble* values = xNew<IssmDouble>(numnodes);
1491 IssmDouble* pressure = xNew<IssmDouble>(numnodes);
1492 IssmDouble* surface = xNew<IssmDouble>(numnodes);
1493 IssmDouble* B = xNew<IssmDouble>(numnodes);
1494 IssmDouble* temperature = xNew<IssmDouble>(numnodes);
1495 IssmDouble* waterfraction = xNew<IssmDouble>(numnodes);
1496
1497 /*Use the dof list to index into the solution vector: */
1498 for(i=0;i<numnodes;i++){
1499 values[i]=solution[doflist[i]];
1500
1501 /*Check solution*/
1502 if(xIsNan<IssmDouble>(values[i])) _error_("NaN found in solution vector");
1503 if(xIsInf<IssmDouble>(values[i])) _error_("Inf found in solution vector");
1504 }
1505
1506 /*Get all inputs and parameters*/
1507 if(element->material->ObjectEnum()!=MatestarEnum) n=element->GetMaterialParameter(MaterialsRheologyNEnum);
1508 element->GetInputValue(&converged,ConvergedEnum);
1509 element->GetInputListOnNodes(&pressure[0],PressureEnum);
1510 if(converged){
1511 for(i=0;i<numnodes;i++){
1512 element->EnthalpyToThermal(&temperature[i],&waterfraction[i],values[i],pressure[i]);
1513 if(waterfraction[i]<0.) _error_("Negative water fraction found in solution vector");
1514 //if(waterfraction[i]>1.) _error_("Water fraction >1 found in solution vector");
1515 }
1516 element->AddInput(EnthalpyEnum,values,element->GetElementType());
1517 element->AddInput(WaterfractionEnum,waterfraction,element->GetElementType());
1518 element->AddInput(TemperatureEnum,temperature,element->GetElementType());
1519
1520 /*Update Rheology only if converged (we must make sure that the temperature is below melting point
1521 * otherwise the rheology could be negative*/
1522 element->FindParam(&rheology_law,MaterialsRheologyLawEnum);
1523 element->GetInputListOnNodes(&surface[0],SurfaceEnum);
1524 switch(rheology_law){
1525 case NoneEnum:
1526 /*Do nothing: B is not temperature dependent*/
1527 break;
1528 case BuddJackaEnum:
1529 for(i=0;i<numnodes;i++) B[i]=BuddJacka(temperature[i]);
1530 element->AddInput(MaterialsRheologyBEnum,&B[0],element->GetElementType());
1531 break;
1532 case CuffeyEnum:
1533 for(i=0;i<numnodes;i++) B[i]=Cuffey(temperature[i]);
1534 element->AddInput(MaterialsRheologyBEnum,&B[0],element->GetElementType());
1535 break;
1536 case CuffeyTemperateEnum:
1537 for(i=0;i<numnodes;i++) B[i]=CuffeyTemperate(temperature[i], waterfraction[i],n);
1538 element->AddInput(MaterialsRheologyBEnum,&B[0],element->GetElementType());
1539 break;
1540 case PatersonEnum:
1541 for(i=0;i<numnodes;i++) B[i]=Paterson(temperature[i]);
1542 element->AddInput(MaterialsRheologyBEnum,&B[0],element->GetElementType());
1543 break;
1544 case ArrheniusEnum:
1545 element->GetVerticesCoordinates(&xyz_list);
1546 for(i=0;i<numnodes;i++) B[i]=Arrhenius(temperature[i],surface[i]-xyz_list[i*3+2],n);
1547 element->AddInput(MaterialsRheologyBEnum,&B[0],element->GetElementType());
1548 break;
1549 case LliboutryDuvalEnum:
1550 for(i=0;i<numnodes;i++) B[i]=LliboutryDuval(values[i],pressure[i],n,element->GetMaterialParameter(MaterialsBetaEnum),element->GetMaterialParameter(ConstantsReferencetemperatureEnum),element->GetMaterialParameter(MaterialsHeatcapacityEnum),element->GetMaterialParameter(MaterialsLatentheatEnum));
1551 element->AddInput(MaterialsRheologyBEnum,&B[0],element->GetElementType());
1552 break;
1553 default: _error_("Rheology law " << EnumToStringx(rheology_law) << " not supported yet");
1554 }
1555 }
1556 else{
1557 element->AddInput(EnthalpyPicardEnum,values,element->GetElementType());
1558 }
1559
1560 /*Free ressources:*/
1561 xDelete<IssmDouble>(values);
1562 xDelete<IssmDouble>(pressure);
1563 xDelete<IssmDouble>(surface);
1564 xDelete<IssmDouble>(B);
1565 xDelete<IssmDouble>(temperature);
1566 xDelete<IssmDouble>(waterfraction);
1567 xDelete<IssmDouble>(xyz_list);
1568 xDelete<int>(doflist);
1569}/*}}}*/
1570void EnthalpyAnalysis::PostProcessing(FemModel* femmodel){/*{{{*/
1571
1572 /*Intermediaries*/
1573 bool computebasalmeltingrates=true;
1574 bool drainicecolumn=true;
1575 IssmDouble dt;
1576
1577 femmodel->parameters->FindParam(&dt,TimesteppingTimeStepEnum);
1578
1579 if(drainicecolumn && (dt>0.)) DrainWaterfraction(femmodel);
1580 if(computebasalmeltingrates) ComputeBasalMeltingrate(femmodel);
1581
1582}/*}}}*/
1583IssmDouble EnthalpyAnalysis::PureIceEnthalpy(Element* element,IssmDouble pressure){/*{{{*/
1584
1585 IssmDouble heatcapacity = element->GetMaterialParameter(MaterialsHeatcapacityEnum);
1586 IssmDouble referencetemperature = element->GetMaterialParameter(ConstantsReferencetemperatureEnum);
1587
1588 return heatcapacity*(TMeltingPoint(element,pressure)-referencetemperature);
1589}/*}}}*/
1590IssmDouble EnthalpyAnalysis::TMeltingPoint(Element* element,IssmDouble pressure){/*{{{*/
1591
1592 IssmDouble meltingpoint = element->GetMaterialParameter(MaterialsMeltingpointEnum);
1593 IssmDouble beta = element->GetMaterialParameter(MaterialsBetaEnum);
1594
1595 return meltingpoint-beta*pressure;
1596}/*}}}*/
1597void EnthalpyAnalysis::UpdateBasalConstraints(FemModel* femmodel){/*{{{*/
1598
1599 /*Update basal dirichlet BCs for enthalpy: */
1600 Vector<IssmDouble>* spc = NULL;
1601 IssmDouble* serial_spc = NULL;
1602
1603 spc=new Vector<IssmDouble>(femmodel->nodes->NumberOfNodes(EnthalpyAnalysisEnum));
1604 /*First create a vector to figure out what elements should be constrained*/
1605 for(int i=0;i<femmodel->elements->Size();i++){
1606 Element* element=xDynamicCast<Element*>(femmodel->elements->GetObjectByOffset(i));
1607 GetBasalConstraints(spc,element);
1608 }
1609
1610 /*Assemble and serialize*/
1611 spc->Assemble();
1612 serial_spc=spc->ToMPISerial();
1613 delete spc;
1614
1615 /*Then update basal constraints nodes accordingly*/
1616 for(int i=0;i<femmodel->elements->Size();i++){
1617 Element* element=xDynamicCast<Element*>(femmodel->elements->GetObjectByOffset(i));
1618 ApplyBasalConstraints(serial_spc,element);
1619 }
1620
1621 femmodel->UpdateConstraintsx();
1622
1623 /*Delete*/
1624 xDelete<IssmDouble>(serial_spc);
1625}/*}}}*/
1626void EnthalpyAnalysis::UpdateConstraints(FemModel* femmodel){/*{{{*/
1627 SetActiveNodesLSMx(femmodel);
1628}/*}}}*/
Note: See TracBrowser for help on using the repository browser.