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

Last change on this file since 26482 was 26482, checked in by vverjans, 3 years ago

CHG: compliance with ISSM coding guidelines, random noise computation uniform accross CPUs in Smbautoregression, remove duplicate code in OceanExchangeDatax

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