source: issm/trunk-jpl/src/c/modules/SurfaceMassBalancex/SurfaceMassBalancex.cpp@ 26481

Last change on this file since 26481 was 26481, checked in by Mathieu Morlighem, 3 years ago

CHG: cosmetics

File size: 19.5 KB
RevLine 
[17085]1/*!\file SurfaceMassBalancex
[23366]2 * \brief: calculates SMB
[17085]3 */
4
[26479]5#include <config.h>
[17085]6#include "./SurfaceMassBalancex.h"
7#include "../../shared/shared.h"
8#include "../../toolkits/toolkits.h"
[23814]9#include "../modules.h"
[25379]10#include "../../classes/Inputs/TransientInput.h"
[26477]11#include "../../shared/Random/random.h"
[17085]12
[23814]13void SmbForcingx(FemModel* femmodel){/*{{{*/
14
15 // void SmbForcingx(smb,ni){
16 // INPUT parameters: ni: working size of arrays
17 // OUTPUT: mass-balance (m/yr ice): agd(NA)
18
19}/*}}}*/
[17085]20void SmbGradientsx(FemModel* femmodel){/*{{{*/
21
22 // void SurfaceMassBalancex(hd,agd,ni){
23 // INPUT parameters: ni: working size of arrays
24 // INPUT: surface elevation (m): hd(NA)
25 // OUTPUT: mass-balance (m/yr ice): agd(NA)
[18001]26 int v;
27 IssmDouble rho_water; // density of fresh water
28 IssmDouble rho_ice; // density of ice
[21469]29 IssmDouble yts; // conversion factor year to second
[17085]30
[18001]31 /*Loop over all the elements of this partition*/
[25539]32 for(Object* & object : femmodel->elements->objects){
33 Element* element=xDynamicCast<Element*>(object);
[18001]34
35 /*Allocate all arrays*/
36 int numvertices = element->GetNumberOfVertices();
37 IssmDouble* Href = xNew<IssmDouble>(numvertices); // reference elevation from which deviations are used to calculate the SMB adjustment
38 IssmDouble* Smbref = xNew<IssmDouble>(numvertices); // reference SMB to which deviations are added
39 IssmDouble* b_pos = xNew<IssmDouble>(numvertices); // Hs-SMB relation parameter
40 IssmDouble* b_neg = xNew<IssmDouble>(numvertices); // Hs-SMB relation paremeter
41 IssmDouble* s = xNew<IssmDouble>(numvertices); // surface elevation (m)
42 IssmDouble* smb = xNew<IssmDouble>(numvertices);
43
44 /*Recover SmbGradients*/
[19527]45 element->GetInputListOnVertices(Href,SmbHrefEnum);
46 element->GetInputListOnVertices(Smbref,SmbSmbrefEnum);
47 element->GetInputListOnVertices(b_pos,SmbBPosEnum);
48 element->GetInputListOnVertices(b_neg,SmbBNegEnum);
[18001]49
[18266]50 /*Recover surface elevation at vertices: */
[18001]51 element->GetInputListOnVertices(s,SurfaceEnum);
52
53 /*Get material parameters :*/
[23644]54 rho_ice=element->FindParam(MaterialsRhoIceEnum);
55 rho_water=element->FindParam(MaterialsRhoFreshwaterEnum);
[18001]56
[21469]57 /* Get constants */
58 femmodel->parameters->FindParam(&yts,ConstantsYtsEnum);
59
[18001]60 // loop over all vertices
61 for(v=0;v<numvertices;v++){
62 if(Smbref[v]>0){
63 smb[v]=Smbref[v]+b_pos[v]*(s[v]-Href[v]);
64 }
65 else{
66 smb[v]=Smbref[v]+b_neg[v]*(s[v]-Href[v]);
67 }
[21469]68
[18266]69 smb[v]=smb[v]/1000*rho_water/rho_ice; // SMB in m/y ice
[18001]70 } //end of the loop over the vertices
71
72 /*Add input to element and Free memory*/
[25379]73 element->AddInput(SmbMassBalanceEnum,smb,P1Enum);
[18001]74 xDelete<IssmDouble>(Href);
75 xDelete<IssmDouble>(Smbref);
76 xDelete<IssmDouble>(b_pos);
77 xDelete<IssmDouble>(b_neg);
78 xDelete<IssmDouble>(s);
79 xDelete<IssmDouble>(smb);
[17085]80 }
81
82}/*}}}*/
[21469]83void SmbGradientsElax(FemModel* femmodel){/*{{{*/
84
85 // void SurfaceMassBalancex(hd,agd,ni){
86 // INPUT parameters: ni: working size of arrays
87 // INPUT: surface elevation (m): hd(NA)
88 // OUTPUT: surface mass-balance (m/yr ice): agd(NA)
89 int v;
90
91 /*Loop over all the elements of this partition*/
[25539]92 for(Object* & object : femmodel->elements->objects){
93 Element* element=xDynamicCast<Element*>(object);
[21469]94
95 /*Allocate all arrays*/
96 int numvertices = element->GetNumberOfVertices();
97 IssmDouble* ela = xNew<IssmDouble>(numvertices); // Equilibrium Line Altitude (m a.s.l) to which deviations are used to calculate the SMB
98 IssmDouble* b_pos = xNew<IssmDouble>(numvertices); // SMB gradient above ELA (m ice eq. per m elevation change)
99 IssmDouble* b_neg = xNew<IssmDouble>(numvertices); // SMB gradient below ELA (m ice eq. per m elevation change)
100 IssmDouble* b_max = xNew<IssmDouble>(numvertices); // Upper cap on SMB rate (m/y ice eq.)
101 IssmDouble* b_min = xNew<IssmDouble>(numvertices); // Lower cap on SMB rate (m/y ice eq.)
102 IssmDouble* s = xNew<IssmDouble>(numvertices); // Surface elevation (m a.s.l.)
103 IssmDouble* smb = xNew<IssmDouble>(numvertices); // SMB (m/y ice eq.)
104
105 /*Recover ELA, SMB gradients, and caps*/
106 element->GetInputListOnVertices(ela,SmbElaEnum);
107 element->GetInputListOnVertices(b_pos,SmbBPosEnum);
108 element->GetInputListOnVertices(b_neg,SmbBNegEnum);
109 element->GetInputListOnVertices(b_max,SmbBMaxEnum);
110 element->GetInputListOnVertices(b_min,SmbBMinEnum);
111
112 /*Recover surface elevation at vertices: */
113 element->GetInputListOnVertices(s,SurfaceEnum);
114
115 /*Loop over all vertices, calculate SMB*/
116 for(v=0;v<numvertices;v++){
117 // if surface is above the ELA
[23366]118 if(s[v]>ela[v]){
[21469]119 smb[v]=b_pos[v]*(s[v]-ela[v]);
120 }
121 // if surface is below or equal to the ELA
122 else{
123 smb[v]=b_neg[v]*(s[v]-ela[v]);
124 }
125
126 // if SMB is larger than upper cap, set SMB to upper cap
127 if(smb[v]>b_max[v]){
128 smb[v]=b_max[v];
129 }
130 // if SMB is smaller than lower cap, set SMB to lower cap
131 if(smb[v]<b_min[v]){
132 smb[v]=b_min[v];
133 }
134 } //end of the loop over the vertices
135
136 /*Add input to element and Free memory*/
[25379]137 element->AddInput(SmbMassBalanceEnum,smb,P1Enum);
[21469]138 xDelete<IssmDouble>(ela);
139 xDelete<IssmDouble>(b_pos);
140 xDelete<IssmDouble>(b_neg);
141 xDelete<IssmDouble>(b_max);
142 xDelete<IssmDouble>(b_min);
143 xDelete<IssmDouble>(s);
144 xDelete<IssmDouble>(smb);
145
146 }
147
148}/*}}}*/
[26477]149void SmbautoregressionInitx(FemModel* femmodel){/*{{{*/
150 /*Initialization step of Smbautoregressionx*/
151 int M,N,Nphi,arorder,numbasins;
152 IssmDouble starttime,tstep_ar,tinit_ar;
153 IssmDouble* beta0 = xNew<IssmDouble>(numbasins);
154 IssmDouble* beta1 = xNew<IssmDouble>(numbasins);
155 IssmDouble* phi = xNew<IssmDouble>(numbasins*arorder);
156 IssmDouble* covmat = xNew<IssmDouble>(numbasins*numbasins);
157 femmodel->parameters->FindParam(&starttime,TimesteppingStartTimeEnum);
158 femmodel->parameters->FindParam(&tstep_ar,SmbAutoregressionTimestepEnum);
159 femmodel->parameters->FindParam(&tinit_ar,SmbAutoregressionInitialTimeEnum);
160 femmodel->parameters->FindParam(&numbasins,SmbNumBasinsEnum);
161 femmodel->parameters->FindParam(&arorder,SmbAutoregressiveOrderEnum);
162 femmodel->parameters->FindParam(&beta0,&M,SmbBeta0Enum); _assert_(M==numbasins);
163 femmodel->parameters->FindParam(&beta1,&M,SmbBeta1Enum); _assert_(M==numbasins);
164 femmodel->parameters->FindParam(&phi,&M,&Nphi,SmbPhiEnum); _assert_(M==numbasins); _assert_(Nphi==arorder);
165 femmodel->parameters->FindParam(&covmat,&M,&N,SmbCovmatEnum); _assert_(M==numbasins); _assert_(N==numbasins);
166 /*AR model spin-up*/
167 int nspin{2*arorder+5}; //number of spin-up steps to be executed
168 IssmDouble* noisespin = xNewZeroInit<IssmDouble>(numbasins*nspin); //sample of basin-specific noise at each spinup step
169 for(int ii{0};ii<nspin;ii++){
170 IssmDouble* temparray = xNew<IssmDouble>(numbasins);
171 multivariateNormal(&temparray,numbasins,0.0,covmat);
172 for(int jj{0};jj<numbasins;jj++){noisespin[ii*numbasins+jj]=temparray[jj];}
173 xDelete<IssmDouble>(temparray);
174 }
175 /*Initialize SmbValuesAutoregressionEnum*/
176 for(Object* &object:femmodel->elements->objects){
177 Element* element = xDynamicCast<Element*>(object); //generate element object
178 element->SmbautoregressionInit(numbasins,arorder,nspin,starttime,tstep_ar,tinit_ar,beta0,beta1,phi,noisespin);
179 }
180 /*Cleanup*/
181 xDelete<IssmDouble>(beta0);
182 xDelete<IssmDouble>(beta1);
183 xDelete<IssmDouble>(phi);
184 xDelete<IssmDouble>(noisespin);
185 xDelete<IssmDouble>(covmat);
186}/*}}}*/
187void Smbautoregressionx(FemModel* femmodel){/*{{{*/
188 /*Get time parameters*/
189 IssmDouble time,dt,starttime,tstep_ar;
190 femmodel->parameters->FindParam(&time,TimeEnum);
191 femmodel->parameters->FindParam(&dt,TimesteppingTimeStepEnum);
192 femmodel->parameters->FindParam(&starttime,TimesteppingStartTimeEnum);
193 femmodel->parameters->FindParam(&tstep_ar,SmbAutoregressionTimestepEnum);
194 /*Initialize module at first time step*/
195 if(time<=starttime+dt){SmbautoregressionInitx(femmodel);}
196 /*Determine if this is a time step for the AR model*/
197 bool isstepforar{false};
[26479]198
199 #ifndef _HAVE_AD_
200 if((fmod(time,tstep_ar)<fmod((time-dt),tstep_ar)) || (time<=starttime+dt) || tstep_ar==dt){isstepforar = true;}
201 #else
202 _error_("not implemented yet");
203 #endif
204
[26477]205 /*Load parameters*/
206 int M,N,Nphi,arorder,numbasins;
207 femmodel->parameters->FindParam(&numbasins,SmbNumBasinsEnum);
208 femmodel->parameters->FindParam(&arorder,SmbAutoregressiveOrderEnum);
209 IssmDouble tinit_ar;
[26481]210 IssmDouble* beta0 = xNew<IssmDouble>(numbasins);
211 IssmDouble* beta1 = xNew<IssmDouble>(numbasins);
212 IssmDouble* phi = xNew<IssmDouble>(numbasins*arorder);
213 IssmDouble* noiseterms = xNew<IssmDouble>(numbasins);
214 IssmDouble* covmat = xNew<IssmDouble>(numbasins*numbasins);
[26477]215 femmodel->parameters->FindParam(&tinit_ar,SmbAutoregressionInitialTimeEnum);
216 femmodel->parameters->FindParam(&beta0,&M,SmbBeta0Enum); _assert_(M==numbasins);
217 femmodel->parameters->FindParam(&beta1,&M,SmbBeta1Enum); _assert_(M==numbasins);
218 femmodel->parameters->FindParam(&phi,&M,&Nphi,SmbPhiEnum); _assert_(M==numbasins); _assert_(Nphi==arorder);
219 femmodel->parameters->FindParam(&covmat,&M,&N,SmbCovmatEnum); _assert_(M==numbasins); _assert_(N==numbasins);
[26481]220
221 /*time elapsed with respect to AR model initial time*/
222 IssmDouble telapsed_ar = time-tinit_ar;
223
[26477]224 /*Before looping through elements: compute noise term specific to each basin from covmat*/
225 multivariateNormal(&noiseterms,numbasins,0.0,covmat);
[26481]226
[26477]227 /*Loop over each element to compute SMB at vertices*/
228 for(Object* &object:femmodel->elements->objects){
[26481]229 Element* element = xDynamicCast<Element*>(object);
[26477]230 element->Smbautoregression(isstepforar,arorder,telapsed_ar,beta0,beta1,phi,noiseterms);
231 }
[26481]232
[26477]233 /*Cleanup*/
234 xDelete<IssmDouble>(beta0);
235 xDelete<IssmDouble>(beta1);
236 xDelete<IssmDouble>(phi);
237 xDelete<IssmDouble>(noiseterms);
238 xDelete<IssmDouble>(covmat);
239}/*}}}*/
[17085]240void Delta18oParameterizationx(FemModel* femmodel){/*{{{*/
241
[25539]242 for(Object* & object : femmodel->elements->objects){
243 Element* element=xDynamicCast<Element*>(object);
[17085]244 element->Delta18oParameterization();
245 }
246
247}/*}}}*/
[18968]248void MungsmtpParameterizationx(FemModel* femmodel){/*{{{*/
249
[25539]250 for(Object* & object : femmodel->elements->objects){
251 Element* element=xDynamicCast<Element*>(object);
[18968]252 element->MungsmtpParameterization();
253 }
254
255}/*}}}*/
[19172]256void Delta18opdParameterizationx(FemModel* femmodel){/*{{{*/
257
[25539]258 for(Object* & object : femmodel->elements->objects){
259 Element* element=xDynamicCast<Element*>(object);
[19172]260 element->Delta18opdParameterization();
261 }
262
263}/*}}}*/
[17085]264void PositiveDegreeDayx(FemModel* femmodel){/*{{{*/
265
266 // void PositiveDegreeDayx(hd,vTempsea,vPrec,agd,Tsurf,ni){
267 // note "v" prefix means 12 monthly means, ie time dimension
268 // INPUT parameters: ni: working size of arrays
269 // INPUT: surface elevation (m): hd(NA)
270 // monthly mean surface sealevel temperature (degrees C): vTempsea(NA
[23366]271 // ,NTIME)
[17085]272 // monthly mean precip rate (m/yr water equivalent): vPrec(NA,NTIME)
273 // OUTPUT: mass-balance (m/yr ice): agd(NA)
274 // mean annual surface temperature (degrees C): Tsurf(NA)
275
[25539]276 int it, jj, itm;
[17085]277 IssmDouble DT = 0.02, sigfac, snormfac;
[23366]278 IssmDouble signorm = 5.5; // signorm : sigma of the temperature distribution for a normal day
[17085]279 IssmDouble siglim; // sigma limit for the integration which is equal to 2.5 sigmanorm
280 IssmDouble signormc = signorm - 0.5; // sigma of the temperature distribution for cloudy day
281 IssmDouble siglimc, siglim0, siglim0c;
282 IssmDouble tstep, tsint, tint, tstepc;
283 int NPDMAX = 1504, NPDCMAX = 1454;
[23366]284 //IssmDouble pdds[NPDMAX]={0};
[17085]285 //IssmDouble pds[NPDCMAX]={0};
286 IssmDouble pddt, pd ; // pd : snow/precip fraction, precipitation falling as snow
287 IssmDouble PDup, PDCUT = 2.0; // PDcut: rain/snow cutoff temperature (C)
288 IssmDouble tstar; // monthly mean surface temp
289
[18968]290 bool ismungsm;
[22448]291 bool issetpddfac;
[18968]292
[17085]293 IssmDouble *pdds = NULL;
294 IssmDouble *pds = NULL;
295 Element *element = NULL;
296
[23366]297 pdds=xNew<IssmDouble>(NPDMAX+1);
298 pds=xNew<IssmDouble>(NPDCMAX+1);
[17085]299
[18968]300 // Get ismungsm parameter
[19527]301 femmodel->parameters->FindParam(&ismungsm,SmbIsmungsmEnum);
[18968]302
[22448]303 // Get issetpddfac parameter
304 femmodel->parameters->FindParam(&issetpddfac,SmbIssetpddfacEnum);
305
[17085]306 /* initialize PDD (creation of a lookup table)*/
307 tstep = 0.1;
308 tsint = tstep*0.5;
309 sigfac = -1.0/(2.0*pow(signorm,2));
310 snormfac = 1.0/(signorm*sqrt(2.0*acos(-1.0)));
311 siglim = 2.5*signorm;
312 siglimc = 2.5*signormc;
313 siglim0 = siglim/DT + 0.5;
314 siglim0c = siglimc/DT + 0.5;
315 PDup = siglimc+PDCUT;
316
317 itm = reCast<int,IssmDouble>((2*siglim/DT + 1.5));
318
319 if(itm >= NPDMAX) _error_("increase NPDMAX in massBalance.cpp");
[23366]320 for(it = 0; it < itm; it++){
[17085]321 // tstar = REAL(it)*DT-siglim;
322 tstar = it*DT-siglim;
323 tint = tsint;
324 pddt = 0.;
325 for ( jj = 0; jj < 600; jj++){
326 if (tint > (tstar+siglim)){break;}
327 pddt = pddt + tint*exp(sigfac*(pow((tint-tstar),2)))*tstep;
328 tint = tint+tstep;
329 }
330 pdds[it] = pddt*snormfac;
331 }
332 pdds[itm+1] = siglim + DT;
333
334 //*********compute PD(T) : snow/precip fraction. precipitation falling as snow
335 tstepc = 0.1;
336 tsint = PDCUT-tstepc*0.5;
337 signormc = signorm - 0.5;
338 sigfac = -1.0/(2.0*pow(signormc,2));
339 snormfac = 1.0/(signormc*sqrt(2.0*acos(-1.0)));
340 siglimc = 2.5*signormc ;
341 itm = reCast<int,IssmDouble>((PDCUT+2.*siglimc)/DT + 1.5);
342 if(itm >= NPDCMAX) _error_("increase NPDCMAX in p35com");
343 for(it = 0; it < itm; it++ ){
344 tstar = it*DT-siglimc;
345 // tstar = REAL(it)*DT-siglimc;
346 tint = tsint; // start against upper bound
347 pd = 0.;
348 for (jj = 0; jj < 600; jj++){
349 if (tint<(tstar-siglimc)) {break;}
350 pd = pd + exp(sigfac*(pow((tint-tstar),2)))*tstepc;
351 tint = tint-tstepc;
352 }
353 pds[it] = pd*snormfac; // gaussian integral lookup table for snow fraction
354 }
355 pds[itm+1] = 0.;
356 // *******END initialize PDD
357
[25539]358 for(Object* & object : femmodel->elements->objects){
359 element=xDynamicCast<Element*>(object);
[22448]360 element->PositiveDegreeDay(pdds,pds,signorm,ismungsm,issetpddfac);
[17085]361 }
362 /*free ressouces: */
363 xDelete<IssmDouble>(pdds);
364 xDelete<IssmDouble>(pds);
365}/*}}}*/
[23317]366void PositiveDegreeDaySicopolisx(FemModel* femmodel){/*{{{*/
[23366]367
[23317]368 bool isfirnwarming;
[23328]369 femmodel->parameters->FindParam(&isfirnwarming,SmbIsfirnwarmingEnum);
[23366]370
[25539]371 for(Object* & object : femmodel->elements->objects){
372 Element* element=xDynamicCast<Element*>(object);
[23317]373 element->PositiveDegreeDaySicopolis(isfirnwarming);
374 }
375
376}/*}}}*/
[17085]377void SmbHenningx(FemModel* femmodel){/*{{{*/
378
[17087]379 /*Intermediaries*/
[17403]380 IssmDouble z_critical = 1675.;
381 IssmDouble dz = 0;
382 IssmDouble a = -15.86;
383 IssmDouble b = 0.00969;
384 IssmDouble c = -0.235;
385 IssmDouble f = 1.;
386 IssmDouble g = -0.0011;
387 IssmDouble h = -1.54e-5;
388 IssmDouble smb,smbref,anomaly,yts,z;
[22249]389
390 /* Get constants */
391 femmodel->parameters->FindParam(&yts,ConstantsYtsEnum);
392 /*iomodel->FindConstant(&yts,"md.constants.yts");*/
393 /*this->parameters->FindParam(&yts,ConstantsYtsEnum);*/
394 /*Mathieu original*/
395 /*IssmDouble smb,smbref,z;*/
396
[17087]397 /*Loop over all the elements of this partition*/
[25539]398 for(Object* & object : femmodel->elements->objects){
399 Element* element=xDynamicCast<Element*>(object);
[17087]400
401 /*Get reference SMB (uncorrected) and allocate all arrays*/
402 int numvertices = element->GetNumberOfVertices();
403 IssmDouble* surfacelist = xNew<IssmDouble>(numvertices);
404 IssmDouble* smblistref = xNew<IssmDouble>(numvertices);
405 IssmDouble* smblist = xNew<IssmDouble>(numvertices);
406 element->GetInputListOnVertices(surfacelist,SurfaceEnum);
[19527]407 element->GetInputListOnVertices(smblistref,SmbSmbrefEnum);
[17087]408
409 /*Loop over all vertices of element and correct SMB as a function of altitude z*/
410 for(int v=0;v<numvertices;v++){
411
[17403]412 /*Get vertex elevation, anoma smb*/
[17087]413 z = surfacelist[v];
[17403]414 anomaly = smblistref[v];
[17087]415
[22249]416 /* Henning edited acc. to Riannes equations*/
417 /* Set SMB maximum elevation, if dz = 0 -> z_critical = 1675 */
418 z_critical = z_critical + dz;
419
420 /* Calculate smb acc. to the surface elevation z */
421 if(z<z_critical){
[17403]422 smb = a + b*z + c;
[17087]423 }
424 else{
[22249]425 smb = (a + b*z)*(f + g*(z-z_critical) + h*(z-z_critical)*(z-z_critical)) + c;
[17087]426 }
[22249]427
[18584]428 /* Compute smb including anomaly,
429 correct for number of seconds in a year [s/yr]*/
430 smb = smb/yts + anomaly;
431
[17087]432 /*Update array accordingly*/
433 smblist[v] = smb;
434
435 }
436
437 /*Add input to element and Free memory*/
[25379]438 element->AddInput(SmbMassBalanceEnum,smblist,P1Enum);
[17087]439 xDelete<IssmDouble>(surfacelist);
440 xDelete<IssmDouble>(smblistref);
441 xDelete<IssmDouble>(smblist);
[17085]442 }
443
444}/*}}}*/
[18001]445void SmbComponentsx(FemModel* femmodel){/*{{{*/
446
447 // void SmbComponentsx(acc,evap,runoff,ni){
448 // INPUT parameters: ni: working size of arrays
449 // INPUT: surface accumulation (m/yr water equivalent): acc
450 // surface evaporation (m/yr water equivalent): evap
451 // surface runoff (m/yr water equivalent): runoff
452 // OUTPUT: mass-balance (m/yr ice): agd(NA)
[23814]453
[18001]454 /*Loop over all the elements of this partition*/
[25539]455 for(Object* & object : femmodel->elements->objects){
456 Element* element=xDynamicCast<Element*>(object);
[18001]457
458 /*Allocate all arrays*/
459 int numvertices = element->GetNumberOfVertices();
[23366]460 IssmDouble* acc = xNew<IssmDouble>(numvertices);
[18001]461 IssmDouble* evap = xNew<IssmDouble>(numvertices);
[23366]462 IssmDouble* runoff = xNew<IssmDouble>(numvertices);
[18001]463 IssmDouble* smb = xNew<IssmDouble>(numvertices);
464
465 /*Recover Smb Components*/
[26208]466 element->GetInputListOnVertices(acc,SmbAccumulationEnum);
467 element->GetInputListOnVertices(evap,SmbEvaporationEnum);
468 element->GetInputListOnVertices(runoff,SmbRunoffEnum);
[18001]469
470 // loop over all vertices
[24335]471 for(int v=0;v<numvertices;v++) smb[v]=acc[v]-evap[v]-runoff[v];
[18001]472
473 /*Add input to element and Free memory*/
[25379]474 element->AddInput(SmbMassBalanceEnum,smb,P1Enum);
[18001]475 xDelete<IssmDouble>(acc);
476 xDelete<IssmDouble>(evap);
477 xDelete<IssmDouble>(runoff);
478 xDelete<IssmDouble>(smb);
479 }
480
481}/*}}}*/
482void SmbMeltComponentsx(FemModel* femmodel){/*{{{*/
483
484 // void SmbMeltComponentsx(acc,evap,melt,refreeze,ni){
485 // INPUT parameters: ni: working size of arrays
486 // INPUT: surface accumulation (m/yr water equivalent): acc
487 // surface evaporation (m/yr water equivalent): evap
488 // surface melt (m/yr water equivalent): melt
489 // refreeze of surface melt (m/yr water equivalent): refreeze
490 // OUTPUT: mass-balance (m/yr ice): agd(NA)
[23814]491
[18001]492 /*Loop over all the elements of this partition*/
[25539]493 for(Object* & object : femmodel->elements->objects){
494 Element* element=xDynamicCast<Element*>(object);
[18001]495
496 /*Allocate all arrays*/
497 int numvertices = element->GetNumberOfVertices();
498 IssmDouble* acc = xNew<IssmDouble>(numvertices);
[23366]499 IssmDouble* evap = xNew<IssmDouble>(numvertices);
[18001]500 IssmDouble* melt = xNew<IssmDouble>(numvertices);
501 IssmDouble* refreeze = xNew<IssmDouble>(numvertices);
502 IssmDouble* smb = xNew<IssmDouble>(numvertices);
503
504 /*Recover Smb Components*/
[26208]505 element->GetInputListOnVertices(acc,SmbAccumulationEnum);
506 element->GetInputListOnVertices(evap,SmbEvaporationEnum);
507 element->GetInputListOnVertices(melt,SmbMeltEnum);
508 element->GetInputListOnVertices(refreeze,SmbRefreezeEnum);
[18001]509
510 // loop over all vertices
[24335]511 for(int v=0;v<numvertices;v++) smb[v]=acc[v]-evap[v]-melt[v]+refreeze[v];
[18001]512
513 /*Add input to element and Free memory*/
[25379]514 element->AddInput(SmbMassBalanceEnum,smb,P1Enum);
[18001]515 xDelete<IssmDouble>(acc);
516 xDelete<IssmDouble>(evap);
517 xDelete<IssmDouble>(melt);
518 xDelete<IssmDouble>(refreeze);
519 xDelete<IssmDouble>(smb);
520 }
521
522}/*}}}*/
[23366]523void SmbGradientsComponentsx(FemModel* femmodel){/*{{{*/
524
[25539]525 for(Object* & object : femmodel->elements->objects){
526 Element* element=xDynamicCast<Element*>(object);
[23366]527 element->SmbGradCompParameterization();
528 }
529
530}/*}}}*/
[23540]531#ifdef _HAVE_SEMIC_
532void SmbSemicx(FemModel* femmodel){/*{{{*/
533
[25539]534 for(Object* & object : femmodel->elements->objects){
535 Element* element=xDynamicCast<Element*>(object);
[23540]536 element->SmbSemic();
537 }
538
539}/*}}}*/
540#else
541void SmbSemicx(FemModel* femmodel){_error_("SEMIC not installed");}
542#endif //_HAVE_SEMIC_
Note: See TracBrowser for help on using the repository browser.