Actual source code: ex5.c

petsc-3.6.1 2015-08-06
Report Typos and Errors
  1: static char help[] = "Nonlinear, time-dependent. Developed from radiative_surface_balance.c \n";
  2: /*
  3:   Contributed by Steve Froehlich, Illinois Institute of Technology

  5:    Usage:
  6:     mpiexec -n <np> ./ex5 [options]
  7:     ./ex5 -help  [view petsc options]
  8:     ./ex5 -ts_type sundials -ts_view
  9:     ./ex5 -da_grid_x 20 -da_grid_y 20 -log_summary
 10:     ./ex5 -da_grid_x 20 -da_grid_y 20 -ts_type rosw -ts_atol 1.e-6 -ts_rtol 1.e-6
 11:     ./ex5 -drawcontours -draw_pause 0.1 -draw_fields 0,1,2,3,4
 12: */

 14: /*
 15:    -----------------------------------------------------------------------

 17:    Governing equations:

 19:         R      = s*(Ea*Ta^4 - Es*Ts^4)
 20:         SH     = p*Cp*Ch*wind*(Ta - Ts)
 21:         LH     = p*L*Ch*wind*B(q(Ta) - q(Ts))
 22:         G      = k*(Tgnd - Ts)/dz

 24:         Fnet   = R + SH + LH + G

 26:         du/dt  = -u*(du/dx) - v*(du/dy) - 2*omeg*sin(lat)*v - (1/p)*(dP/dx)
 27:         dv/dt  = -u*(dv/dx) - v*(dv/dy) + 2*omeg*sin(lat)*u - (1/p)*(dP/dy)
 28:         dTs/dt = Fnet/(Cp*dz) - Div([u*Ts, v*Ts]) + D*Lap(Ts)
 29:                = Fnet/(Cs*dz) - u*(dTs/dx) - v*(dTs/dy) + D*(Ts_xx + Ts_yy)
 30:         dp/dt  = -Div([u*p,v*p])
 31:                = - u*dp/dx - v*dp/dy
 32:         dTa/dt = Fnet/Cp

 34:    Equation of State:

 36:         P = p*R*Ts

 38:    -----------------------------------------------------------------------

 40:    Program considers the evolution of a two dimensional atmosphere from
 41:    sunset to sunrise. There are two components:
 42:                 1. Surface energy balance model to compute diabatic dT (Fnet)
 43:                 2. Dynamical model using simplified primitive equations

 45:    Program is to be initiated at sunset and run to sunrise.

 47:    Inputs are:
 48:                 Surface temperature
 49:                 Dew point temperature
 50:                 Air temperature
 51:                 Temperature at cloud base (if clouds are present)
 52:                 Fraction of sky covered by clouds
 53:                 Wind speed
 54:                 Precipitable water in centimeters
 55:                 Wind direction

 57:    Inputs are are read in from the text file ex5_control.txt. To change an
 58:    input value use ex5_control.txt.

 60:    Solvers:
 61:             Backward Euler = default solver
 62:             Sundials = fastest and most accurate, requires Sundials libraries

 64:    This model is under development and should be used only as an example
 65:    and not as a predictive weather model.
 66: */

 68: #include <petscts.h>
 69: #include <petscdm.h>
 70: #include <petscdmda.h>

 72: /* stefan-boltzmann constant */
 73: #define SIG 0.000000056703
 74: /* absorption-emission constant for surface */
 75: #define EMMSFC   1
 76: /* amount of time (seconds) that passes before new flux is calculated */
 77: #define TIMESTEP 1

 79: /* variables of interest to be solved at each grid point */
 80: typedef struct {
 81:   PetscScalar Ts,Ta; /* surface and air temperature */
 82:   PetscScalar u,v;   /* wind speed */
 83:   PetscScalar p;     /* density */
 84: } Field;

 86: /* User defined variables. Used in solving for variables of interest */
 87: typedef struct {
 88:   DM          da;        /* grid */
 89:   PetscScalar csoil;     /* heat constant for layer */
 90:   PetscScalar dzlay;     /* thickness of top soil layer */
 91:   PetscScalar emma;      /* emission parameter */
 92:   PetscScalar wind;      /* wind speed */
 93:   PetscScalar dewtemp;   /* dew point temperature (moisture in air) */
 94:   PetscScalar pressure1; /* sea level pressure */
 95:   PetscScalar airtemp;   /* temperature of air near boundary layer inversion */
 96:   PetscScalar Ts;        /* temperature at the surface */
 97:   PetscScalar fract;     /* fraction of sky covered by clouds */
 98:   PetscScalar Tc;        /* temperature at base of lowest cloud layer */
 99:   PetscScalar lat;       /* Latitude in degrees */
100:   PetscScalar init;      /* initialization scenario */
101:   PetscScalar deep_grnd_temp; /* temperature of ground under top soil surface layer */
102: } AppCtx;

104: /* Struct for visualization */
105: typedef struct {
106:   PetscBool   drawcontours;   /* flag - 1 indicates drawing contours */
107:   PetscViewer drawviewer;
108:   PetscInt    interval;
109: } MonitorCtx;


112: /* Inputs read in from text file */
113: struct in {
114:   PetscScalar Ts;     /* surface temperature  */
115:   PetscScalar Td;     /* dewpoint temperature */
116:   PetscScalar Tc;     /* temperature of cloud base */
117:   PetscScalar fr;     /* fraction of sky covered by clouds */
118:   PetscScalar wnd;    /* wind speed */
119:   PetscScalar Ta;     /* air temperature */
120:   PetscScalar pwt;    /* precipitable water */
121:   PetscScalar wndDir; /* wind direction */
122:   PetscScalar lat;    /* latitude */
123:   PetscReal   time;   /* time in hours */
124:   PetscScalar init;
125: };

127: /* functions */
128: extern PetscScalar emission(PetscScalar);                           /* sets emission/absorption constant depending on water vapor content */
129: extern PetscScalar calc_q(PetscScalar);                             /* calculates specific humidity */
130: extern PetscScalar mph2mpers(PetscScalar);                          /* converts miles per hour to meters per second */
131: extern PetscScalar Lconst(PetscScalar);                             /* calculates latent heat constant taken from Satellite estimates of wind speed and latent heat flux over the global oceans., Bentamy et al. */
132: extern PetscScalar fahr_to_cel(PetscScalar);                        /* converts Fahrenheit to Celsius */
133: extern PetscScalar cel_to_fahr(PetscScalar);                        /* converts Celsius to Fahrenheit */
134: extern PetscScalar calcmixingr(PetscScalar, PetscScalar);           /* calculates mixing ratio */
135: extern PetscScalar cloud(PetscScalar);                              /* cloud radiative parameterization */
136: extern PetscErrorCode FormInitialSolution(DM,Vec,void*);            /* Specifies initial conditions for the system of equations (PETSc defined function) */
137: extern PetscErrorCode RhsFunc(TS,PetscReal,Vec,Vec,void*);          /* Specifies the user defined functions                     (PETSc defined function) */
138: extern PetscErrorCode Monitor(TS,PetscInt,PetscReal,Vec,void*);     /* Specifies output and visualization tools                 (PETSc defined function) */
139: extern void readinput(struct in *put);                              /* reads input from text file */
140: extern PetscErrorCode calcfluxs(PetscScalar, PetscScalar, PetscScalar, PetscScalar, PetscScalar, PetscScalar*); /* calculates upward IR from surface */
141: extern PetscErrorCode calcfluxa(PetscScalar, PetscScalar, PetscScalar, PetscScalar*);                           /* calculates downward IR from atmosphere */
142: extern PetscErrorCode sensibleflux(PetscScalar, PetscScalar, PetscScalar, PetscScalar*);                        /* calculates sensible heat flux */
143: extern PetscErrorCode potential_temperature(PetscScalar, PetscScalar, PetscScalar, PetscScalar, PetscScalar*);  /* calculates potential temperature */
144: extern PetscErrorCode latentflux(PetscScalar, PetscScalar, PetscScalar, PetscScalar, PetscScalar*);             /* calculates latent heat flux */
145: extern PetscErrorCode calc_gflux(PetscScalar, PetscScalar, PetscScalar*);                                       /* calculates flux between top soil layer and underlying earth */

149: int main(int argc,char **argv)
150: {
152:   int            time;           /* amount of loops */
153:   struct in      put;
154:   PetscScalar    rh;             /* relative humidity */
155:   PetscScalar    x;              /* memory varialbe for relative humidity calculation */
156:   PetscScalar    deep_grnd_temp; /* temperature of ground under top soil surface layer */
157:   PetscScalar    emma;           /* absorption-emission constant for air */
158:   PetscScalar    pressure1 = 101300; /* surface pressure */
159:   PetscScalar    mixratio;       /* mixing ratio */
160:   PetscScalar    airtemp;        /* temperature of air near boundary layer inversion */
161:   PetscScalar    dewtemp;        /* dew point temperature */
162:   PetscScalar    sfctemp;        /* temperature at surface */
163:   PetscScalar    pwat;           /* total column precipitable water */
164:   PetscScalar    cloudTemp;      /* temperature at base of cloud */
165:   AppCtx         user;           /*  user-defined work context */
166:   MonitorCtx     usermonitor;    /* user-defined monitor context */
167:   PetscMPIInt    rank,size;
168:   TS             ts;
169:   SNES           snes;
170:   DM             da;
171:   Vec            T,rhs;          /* solution vector */
172:   Mat            J;              /* Jacobian matrix */
173:   PetscReal      ftime,dt;
174:   PetscInt       steps,dof = 5;
175:   PetscBool      use_coloring  = PETSC_TRUE;
176:   MatFDColoring  matfdcoloring = 0;
177:   PetscBool      monitor_off = PETSC_FALSE;

179:   PetscInitialize(&argc,&argv,(char*)0,help);
180:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
181:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);

183:   /* Inputs */
184:   readinput(&put);

186:   sfctemp   = put.Ts;
187:   dewtemp   = put.Td;
188:   cloudTemp = put.Tc;
189:   airtemp   = put.Ta;
190:   pwat      = put.pwt;

192:   if (!rank) PetscPrintf(PETSC_COMM_SELF,"Initial Temperature = %g\n",(double)sfctemp); /* input surface temperature */

194:   deep_grnd_temp = sfctemp - 10;   /* set underlying ground layer temperature */
195:   emma           = emission(pwat); /* accounts for radiative effects of water vapor */

197:   /* Converts from Fahrenheit to Celsuis */
198:   sfctemp        = fahr_to_cel(sfctemp);
199:   airtemp        = fahr_to_cel(airtemp);
200:   dewtemp        = fahr_to_cel(dewtemp);
201:   cloudTemp      = fahr_to_cel(cloudTemp);
202:   deep_grnd_temp = fahr_to_cel(deep_grnd_temp);

204:   /* Converts from Celsius to Kelvin */
205:   sfctemp        += 273;
206:   airtemp        += 273;
207:   dewtemp        += 273;
208:   cloudTemp      += 273;
209:   deep_grnd_temp += 273;

211:   /* Calculates initial relative humidity */
212:   x        = calcmixingr(dewtemp,pressure1);
213:   mixratio = calcmixingr(sfctemp,pressure1);
214:   rh       = (x/mixratio)*100;

216:   if (!rank) printf("Initial RH = %.1f percent\n\n",(double)rh);   /* prints initial relative humidity */

218:   time = 3600*put.time;                         /* sets amount of timesteps to run model */

220:   /* Configure PETSc TS solver */
221:   /*------------------------------------------*/

223:   /* Create grid */
224:   DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_PERIODIC,DM_BOUNDARY_PERIODIC,DMDA_STENCIL_STAR,-20,-20,
225:                       PETSC_DECIDE,PETSC_DECIDE,dof,1,NULL,NULL,&da);
226:   DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0);

228:   /* Define output window for each variable of interest */
229:   DMDASetFieldName(da,0,"Ts");
230:   DMDASetFieldName(da,1,"Ta");
231:   DMDASetFieldName(da,2,"u");
232:   DMDASetFieldName(da,3,"v");
233:   DMDASetFieldName(da,4,"p");

235:   /* set values for appctx */
236:   user.da             = da;
237:   user.Ts             = sfctemp;
238:   user.fract          = put.fr;          /* fraction of sky covered by clouds */
239:   user.dewtemp        = dewtemp;         /* dew point temperature (mositure in air) */
240:   user.csoil          = 2000000;         /* heat constant for layer */
241:   user.dzlay          = 0.08;            /* thickness of top soil layer */
242:   user.emma           = emma;            /* emission parameter */
243:   user.wind           = put.wnd;         /* wind spped */
244:   user.pressure1      = pressure1;       /* sea level pressure */
245:   user.airtemp        = airtemp;         /* temperature of air near boundar layer inversion */
246:   user.Tc             = cloudTemp;       /* temperature at base of lowest cloud layer */
247:   user.init           = put.init;        /* user chosen initiation scenario */
248:   user.lat            = 70*0.0174532;    /* converts latitude degrees to latitude in radians */
249:   user.deep_grnd_temp = deep_grnd_temp;  /* temp in lowest ground layer */

251:   /* set values for MonitorCtx */
252:   usermonitor.drawcontours = PETSC_FALSE;
253:   PetscOptionsHasName(NULL,"-drawcontours",&usermonitor.drawcontours);
254:   if (usermonitor.drawcontours) {
255:     PetscReal bounds[] = {1000.0,-1000.,  -1000.,-1000.,  1000.,-1000.,  1000.,-1000.,  1000,-1000, 100700,100800};
256:     PetscViewerDrawOpen(PETSC_COMM_WORLD,0,0,0,0,300,300,&usermonitor.drawviewer);
257:     PetscViewerDrawSetBounds(usermonitor.drawviewer,dof,bounds);
258:   }
259:   usermonitor.interval = 1;
260:   PetscOptionsGetInt(NULL,"-monitor_interval",&usermonitor.interval,NULL);

262:   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
263:      Extract global vectors from DA;
264:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
265:   DMCreateGlobalVector(da,&T);
266:   VecDuplicate(T,&rhs); /* r: vector to put the computed right hand side */

268:   TSCreate(PETSC_COMM_WORLD,&ts);
269:   TSSetProblemType(ts,TS_NONLINEAR);
270:   TSSetType(ts,TSBEULER);
271:   TSSetRHSFunction(ts,rhs,RhsFunc,&user);

273:   /* Set Jacobian evaluation routine - use coloring to compute finite difference Jacobian efficiently */
274:   DMSetMatType(da,MATAIJ);
275:   DMCreateMatrix(da,&J);
276:   TSGetSNES(ts,&snes);
277:   if (use_coloring) {
278:     ISColoring iscoloring;
279:     DMCreateColoring(da,IS_COLORING_GLOBAL,&iscoloring);
280:     MatFDColoringCreate(J,iscoloring,&matfdcoloring);
281:     MatFDColoringSetFromOptions(matfdcoloring);
282:     MatFDColoringSetUp(J,iscoloring,matfdcoloring);
283:     ISColoringDestroy(&iscoloring);
284:     MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode (*)(void))SNESTSFormFunction,ts);
285:     SNESSetJacobian(snes,J,J,SNESComputeJacobianDefaultColor,matfdcoloring);
286:   } else {
287:     SNESSetJacobian(snes,J,J,SNESComputeJacobianDefault,NULL);
288:   }

290:   /* Define what to print for ts_monitor option */
291:   PetscOptionsHasName(NULL,"-monitor_off",&monitor_off);
292:   if (!monitor_off) {
293:     TSMonitorSet(ts,Monitor,&usermonitor,NULL);
294:   }
295:   FormInitialSolution(da,T,&user);
296:   dt    = TIMESTEP; /* initial time step */
297:   ftime = TIMESTEP*time;
298:   if (!rank) printf("time %d, ftime %g hour, TIMESTEP %g\n",time,(double)(ftime/3600),(double)dt);

300:   TSSetInitialTimeStep(ts,0.0,dt);
301:   TSSetDuration(ts,time,ftime);
302:   TSSetSolution(ts,T);
303:   TSSetDM(ts,da);

305:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
306:      Set runtime options
307:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
308:   TSSetFromOptions(ts);

310:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
311:      Solve nonlinear system
312:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
313:   TSSolve(ts,T);
314:   TSGetSolveTime(ts,&ftime);
315:   TSGetTimeStepNumber(ts,&steps);
316:   if (!rank) PetscPrintf(PETSC_COMM_WORLD,"Solution T after %g hours %d steps\n",(double)(ftime/3600),steps);


319:   if (matfdcoloring) {MatFDColoringDestroy(&matfdcoloring);}
320:   if (usermonitor.drawcontours) {
321:     PetscViewerDestroy(&usermonitor.drawviewer);
322:   }
323:   MatDestroy(&J);
324:   VecDestroy(&T);
325:   VecDestroy(&rhs);
326:   TSDestroy(&ts);
327:   DMDestroy(&da);

329:   PetscFinalize();
330:   return 0;
331: }
332: /*****************************end main program********************************/
333: /*****************************************************************************/
334: /*****************************************************************************/
335: /*****************************************************************************/
338: PetscErrorCode calcfluxs(PetscScalar sfctemp, PetscScalar airtemp, PetscScalar emma, PetscScalar fract, PetscScalar cloudTemp, PetscScalar *flux)
339: {
341:   *flux = SIG*((EMMSFC*emma*PetscPowScalarInt(airtemp,4)) + (EMMSFC*fract*(1 - emma)*PetscPowScalarInt(cloudTemp,4)) - (EMMSFC*PetscPowScalarInt(sfctemp,4)));   /* calculates flux using Stefan-Boltzmann relation */
342:   return(0);
343: }

347: PetscErrorCode calcfluxa(PetscScalar sfctemp, PetscScalar airtemp, PetscScalar emma, PetscScalar *flux)   /* this function is not currently called upon */
348: {
349:   PetscScalar emm = 0.001;

352:   *flux = SIG*(-emm*(PetscPowScalarInt(airtemp,4)));      /* calculates flux usinge Stefan-Boltzmann relation */
353:   return(0);
354: }
357: PetscErrorCode sensibleflux(PetscScalar sfctemp, PetscScalar airtemp, PetscScalar wind, PetscScalar *sheat)
358: {
359:   PetscScalar density = 1;    /* air density */
360:   PetscScalar Cp      = 1005; /* heat capicity for dry air */
361:   PetscScalar wndmix;         /* temperature change from wind mixing: wind*Ch */

364:   wndmix = 0.0025 + 0.0042*wind;                               /* regression equation valid for neutral and stable BL */
365:   *sheat = density*Cp*wndmix*(airtemp - sfctemp);              /* calculates sensible heat flux */
366:   return(0);
367: }

371: PetscErrorCode latentflux(PetscScalar sfctemp, PetscScalar dewtemp, PetscScalar wind, PetscScalar pressure1, PetscScalar *latentheat)
372: {
373:   PetscScalar density = 1;   /* density of dry air */
374:   PetscScalar q;             /* actual specific humitity */
375:   PetscScalar qs;            /* saturation specific humidity */
376:   PetscScalar wndmix;        /* temperature change from wind mixing: wind*Ch */
377:   PetscScalar beta = .4;     /* moisture availability */
378:   PetscScalar mr;            /* mixing ratio */
379:   PetscScalar lhcnst;        /* latent heat of vaporization constant = 2501000 J/kg at 0c */
380:                              /* latent heat of saturation const = 2834000 J/kg */
381:                              /* latent heat of fusion const = 333700 J/kg */

384:   wind   = mph2mpers(wind);                /* converts wind from mph to meters per second */
385:   wndmix = 0.0025 + 0.0042*wind;           /* regression equation valid for neutral BL */
386:   lhcnst = Lconst(sfctemp);                /* calculates latent heat of evaporation */
387:   mr     = calcmixingr(sfctemp,pressure1); /* calculates saturation mixing ratio */
388:   qs     = calc_q(mr);                     /* calculates saturation specific humidty */
389:   mr     = calcmixingr(dewtemp,pressure1); /* calculates mixing ratio */
390:   q      = calc_q(mr);                     /* calculates specific humidty */

392:   *latentheat = density*wndmix*beta*lhcnst*(q - qs); /* calculates latent heat flux */
393:   return(0);
394: }

398: PetscErrorCode potential_temperature(PetscScalar temp, PetscScalar pressure1, PetscScalar pressure2, PetscScalar sfctemp, PetscScalar *pottemp)
399: {
400:   PetscScalar kdry;     /* poisson constant for dry atmosphere */
401:   PetscScalar pavg;     /* average atmospheric pressure */
402:   /* PetscScalar mixratio; mixing ratio */
403:   /* PetscScalar kmoist;   poisson constant for moist atmosphere */

406:   /* mixratio = calcmixingr(sfctemp,pressure1); */

408: /* initialize poisson constant */
409:   kdry   = 0.2854;
410:   /* kmoist = 0.2854*(1 - 0.24*mixratio); */

412:   pavg     = ((0.7*pressure1)+pressure2)/2;     /* calculates simple average press */
413:   *pottemp = temp*(PetscPowScalar((pressure1/pavg),kdry)); /* calculates potential temperature */
414:   return(0);
415: }
416: extern PetscScalar calcmixingr(PetscScalar dtemp, PetscScalar pressure1)
417: {
418:   PetscScalar e;        /* vapor pressure */
419:   PetscScalar mixratio; /* mixing ratio */

421:   dtemp    = dtemp - 273;                                /* converts from Kelvin to Celsuis */
422:   e        = 6.11*(PetscPowScalar(10,((7.5*dtemp)/(237.7+dtemp)))); /* converts from dew point temp to vapor pressure */
423:   e        = e*100;                                      /* converts from hPa to Pa */
424:   mixratio = (0.622*e)/(pressure1 - e);                  /* computes mixing ratio */
425:   mixratio = mixratio*1;                                 /* convert to g/Kg */

427:   return mixratio;
428: }
429: extern PetscScalar calc_q(PetscScalar rv)
430: {
431:   PetscScalar specific_humidity;        /* define specific humidity variable */
432:   specific_humidity = rv/(1 + rv);      /* calculates specific humidity */
433:   return specific_humidity;
434: }

438: PetscErrorCode calc_gflux(PetscScalar sfctemp, PetscScalar deep_grnd_temp, PetscScalar* Gflux)
439: {
440:   PetscScalar k;                       /* thermal conductivity parameter */
441:   PetscScalar n                = 0.38; /* value of soil porosity */
442:   PetscScalar dz               = 1;    /* depth of layer between soil surface and deep soil layer */
443:   PetscScalar unit_soil_weight = 2700; /* unit soil weight in kg/m^3 */

446:   k      = ((0.135*(1-n)*unit_soil_weight) + 64.7)/(unit_soil_weight - (0.947*(1-n)*unit_soil_weight)); /* dry soil conductivity */
447:   *Gflux = (k*(deep_grnd_temp - sfctemp)/dz);   /* calculates flux from deep ground layer */
448:   return(0);
449: }
452: extern PetscScalar emission(PetscScalar pwat)
453: {
454:   PetscScalar emma;

456:   emma = 0.725 + 0.17*log10(pwat);

458:   return emma;
459: }
460: extern PetscScalar cloud(PetscScalar fract)
461: {
462:   PetscScalar emma = 0;

464:   /* modifies radiative balance depending on cloud cover */
465:   if (fract >= 0.9)                     emma = 1;
466:   else if (0.9 > fract && fract >= 0.8) emma = 0.9;
467:   else if (0.8 > fract && fract >= 0.7) emma = 0.85;
468:   else if (0.7 > fract && fract >= 0.6) emma = 0.75;
469:   else if (0.6 > fract && fract >= 0.5) emma = 0.65;
470:   else if (0.4 > fract && fract >= 0.3) emma = emma*1.086956;
471:   return emma;
472: }
473: extern PetscScalar Lconst(PetscScalar sfctemp)
474: {
475:   PetscScalar Lheat;
476:   sfctemp -=273;                               /* converts from kelvin to celsius */
477:   Lheat    = 4186.8*(597.31 - 0.5625*sfctemp); /* calculates latent heat constant */
478:   return Lheat;
479: }
480: extern PetscScalar mph2mpers(PetscScalar wind)
481: {
482:   wind = ((wind*1.6*1000)/3600);                 /* converts wind from mph to meters per second */
483:   return wind;
484: }
485: extern PetscScalar fahr_to_cel(PetscScalar temp)
486: {
487:   temp = (5*(temp-32))/9; /* converts from farhrenheit to celsuis */
488:   return temp;
489: }
490: extern PetscScalar cel_to_fahr(PetscScalar temp)
491: {
492:   temp = ((temp*9)/5) + 32; /* converts from celsuis to farhrenheit */
493:   return temp;
494: }
495: void readinput(struct in *put)
496: {
497:   int    i;
498:   char   x;
499:   FILE   *ifp;
500:   double tmp;

502:   ifp = fopen("ex5_control.txt", "r");

504:   for (i=0; i<110; i++) fscanf(ifp, "%c", &x);
505:   fscanf(ifp, "%lf", &tmp);
506:   put->Ts = tmp;

508:   for (i=0; i<43; i++) fscanf(ifp, "%c", &x);
509:   fscanf(ifp, "%lf", &tmp);
510:   put->Td = tmp;

512:   for (i=0; i<43; i++) fscanf(ifp, "%c", &x);
513:   fscanf(ifp, "%lf", &tmp);
514:   put->Ta = tmp;

516:   for (i=0; i<43; i++) fscanf(ifp, "%c", &x);
517:   fscanf(ifp, "%lf", &tmp);
518:   put->Tc = tmp;

520:   for (i=0; i<43; i++) fscanf(ifp, "%c", &x);
521:   fscanf(ifp, "%lf", &tmp);
522:   put->fr = tmp;

524:   for (i=0; i<43; i++) fscanf(ifp, "%c", &x);
525:   fscanf(ifp, "%lf", &tmp);
526:   put->wnd = tmp;

528:   for (i=0; i<43; i++) fscanf(ifp, "%c", &x);
529:   fscanf(ifp, "%lf", &tmp);
530:   put->pwt = tmp;

532:   for (i=0; i<43; i++) fscanf(ifp, "%c", &x);
533:   fscanf(ifp, "%lf", &tmp);
534:   put->wndDir = tmp;

536:   for (i=0; i<43; i++) fscanf(ifp, "%c", &x);
537:   fscanf(ifp, "%lf", &tmp);
538:   put->time = tmp;

540:   for (i=0; i<63; i++) fscanf(ifp, "%c", &x);
541:   fscanf(ifp, "%lf", &tmp);
542:   put->init = tmp;

544: }

546: /* ------------------------------------------------------------------- */
549: PetscErrorCode FormInitialSolution(DM da,Vec Xglobal,void *ctx)
550: {
552:   AppCtx         *user = (AppCtx*)ctx;       /* user-defined application context */
553:   PetscInt       i,j,xs,ys,xm,ym,Mx,My;
554:   Field          **X;

557:   DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
558:                      PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);

560:   /* Get pointers to vector data */
561:   DMDAVecGetArray(da,Xglobal,&X);

563:   /* Get local grid boundaries */
564:   DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);

566:   /* Compute function over the locally owned part of the grid */

568:   if (user->init == 1) {
569:     for (j=ys; j<ys+ym; j++) {
570:       for (i=xs; i<xs+xm; i++) {
571:         X[j][i].Ts = user->Ts - i*0.0001;
572:         X[j][i].Ta = X[j][i].Ts - 5;
573:         X[j][i].u  = 0;
574:         X[j][i].v  = 0;
575:         X[j][i].p  = 1.25;
576:         if ((j == 5 || j == 6) && (i == 4 || i == 5))   X[j][i].p += 0.00001;
577:         if ((j == 5 || j == 6) && (i == 12 || i == 13)) X[j][i].p += 0.00001;
578:       }
579:     }
580:   } else {
581:     for (j=ys; j<ys+ym; j++) {
582:       for (i=xs; i<xs+xm; i++) {
583:         X[j][i].Ts = user->Ts;
584:         X[j][i].Ta = X[j][i].Ts - 5;
585:         X[j][i].u  = 0;
586:         X[j][i].v  = 0;
587:         X[j][i].p  = 1.25;
588:       }
589:     }
590:   }

592:   /* Restore vectors */
593:   DMDAVecRestoreArray(da,Xglobal,&X);
594:   return(0);
595: }

599: /*
600:    RhsFunc - Evaluates nonlinear function F(u).

602:    Input Parameters:
603: .  ts - the TS context
604: .  t - current time
605: .  Xglobal - input vector
606: .  F - output vector
607: .  ptr - optional user-defined context, as set by SNESSetFunction()

609:    Output Parameter:
610: .  F - rhs function vector
611:  */
612: PetscErrorCode RhsFunc(TS ts,PetscReal t,Vec Xglobal,Vec F,void *ctx)
613: {
614:   AppCtx         *user = (AppCtx*)ctx;       /* user-defined application context */
615:   DM             da    = user->da;
617:   PetscInt       i,j,Mx,My,xs,ys,xm,ym;
618:   PetscReal      dhx,dhy;
619:   Vec            localT;
620:   Field          **X,**Frhs;                            /* structures that contain variables of interest and left hand side of governing equations respectively */
621:   PetscScalar    csoil          = user->csoil;          /* heat constant for layer */
622:   PetscScalar    dzlay          = user->dzlay;          /* thickness of top soil layer */
623:   PetscScalar    emma           = user->emma;           /* emission parameter */
624:   PetscScalar    wind           = user->wind;           /* wind speed */
625:   PetscScalar    dewtemp        = user->dewtemp;        /* dew point temperature (moisture in air) */
626:   PetscScalar    pressure1      = user->pressure1;      /* sea level pressure */
627:   PetscScalar    airtemp        = user->airtemp;        /* temperature of air near boundary layer inversion */
628:   PetscScalar    fract          = user->fract;          /* fraction of the sky covered by clouds */
629:   PetscScalar    Tc             = user->Tc;             /* temperature at base of lowest cloud layer */
630:   PetscScalar    lat            = user->lat;            /* latitude */
631:   PetscScalar    Cp             = 1005.7;               /* specific heat of air at constant pressure */
632:   PetscScalar    Rd             = 287.058;              /* gas constant for dry air */
633:   PetscScalar    diffconst      = 1000;                 /* diffusion coefficient */
634:   PetscScalar    f              = 2*0.0000727*PetscSinScalar(lat); /* coriolis force */
635:   PetscScalar    deep_grnd_temp = user->deep_grnd_temp; /* temp in lowest ground layer */
636:   PetscScalar    Ts,u,v,p;
637:   PetscScalar    u_abs,u_plus,u_minus,v_abs,v_plus,v_minus;

639:   PetscScalar sfctemp1,fsfc1,Ra;
640:   PetscScalar sheat;                   /* sensible heat flux */
641:   PetscScalar latentheat;              /* latent heat flux */
642:   PetscScalar groundflux;              /* flux from conduction of deep ground layer in contact with top soil */
643:   PetscInt    xend,yend;

646:   DMGetLocalVector(da,&localT);
647:   DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
648:                      PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);

650:   dhx = (PetscReal)(Mx-1)/(5000*(Mx-1));  /* dhx = 1/dx; assume 2D space domain: [0.0, 1.e5] x [0.0, 1.e5] */
651:   dhy = (PetscReal)(My-1)/(5000*(Mx-1));  /* dhy = 1/dy; */


654:   /*
655:      Scatter ghost points to local vector,using the 2-step process
656:         DAGlobalToLocalBegin(),DAGlobalToLocalEnd().
657:      By placing code between these two statements, computations can be
658:      done while messages are in transition.
659:   */
660:   DMGlobalToLocalBegin(da,Xglobal,INSERT_VALUES,localT);
661:   DMGlobalToLocalEnd(da,Xglobal,INSERT_VALUES,localT);

663:   /* Get pointers to vector data */
664:   DMDAVecGetArrayRead(da,localT,&X);
665:   DMDAVecGetArray(da,F,&Frhs);

667:   /* Get local grid boundaries */
668:   DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);

670:   /* Compute function over the locally owned part of the grid */
671:   /* the interior points */
672:   xend=xs+xm; yend=ys+ym;
673:   for (j=ys; j<yend; j++) {
674:     for (i=xs; i<xend; i++) {
675:       Ts = X[j][i].Ts; u = X[j][i].u; v = X[j][i].v; p = X[j][i].p; /*P = X[j][i].P; */

677:       sfctemp1 = (double)Ts;
678:       sfctemp1 = (double)X[j][i].Ts;
679:       calcfluxs(sfctemp1,airtemp,emma,fract,Tc,&fsfc1);        /* calculates surface net radiative flux */
680:       sensibleflux(sfctemp1,airtemp,wind,&sheat);              /* calculate sensible heat flux */
681:       latentflux(sfctemp1,dewtemp,wind,pressure1,&latentheat); /* calculates latent heat flux */
682:       calc_gflux(sfctemp1,deep_grnd_temp,&groundflux);         /* calculates flux from earth below surface soil layer by conduction */
683:       calcfluxa(sfctemp1,airtemp,emma,&Ra);                                  /* Calculates the change in downward radiative flux */
684:       fsfc1    = fsfc1 + latentheat + sheat + groundflux;                               /* adds radiative, sensible heat, latent heat, and ground heat flux yielding net flux */

686:       /* convective coefficients for upwinding */
687:       u_abs   = PetscAbsScalar(u);
688:       u_plus  = .5*(u + u_abs); /* u if u>0; 0 if u<0 */
689:       u_minus = .5*(u - u_abs); /* u if u <0; 0 if u>0 */

691:       v_abs   = PetscAbsScalar(v);
692:       v_plus  = .5*(v + v_abs); /* v if v>0; 0 if v<0 */
693:       v_minus = .5*(v - v_abs); /* v if v <0; 0 if v>0 */

695:       /* Solve governing equations */
696:       /* P = p*Rd*Ts; */

698:       /* du/dt -> time change of east-west component of the wind */
699:       Frhs[j][i].u = - u_plus*(u - X[j][i-1].u)*dhx - u_minus*(X[j][i+1].u - u)*dhx       /* - u(du/dx) */
700:                      - v_plus*(u - X[j-1][i].u)*dhy - v_minus*(X[j+1][i].u - u)*dhy       /* - v(du/dy) */
701:                      -(Rd/p)*(Ts*(X[j][i+1].p - X[j][i-1].p)*0.5*dhx  + p*0*(X[j][i+1].Ts - X[j][i-1].Ts)*0.5*dhx) /* -(R/p)[Ts(dp/dx)+ p(dTs/dx)] */
702: /*                     -(1/p)*(X[j][i+1].P - X[j][i-1].P)*dhx */
703:                      + f*v;

705:       /* dv/dt -> time change of north-south component of the wind */
706:       Frhs[j][i].v = - u_plus*(v - X[j][i-1].v)*dhx - u_minus*(X[j][i+1].v - v)*dhx       /* - u(dv/dx) */
707:                      - v_plus*(v - X[j-1][i].v)*dhy - v_minus*(X[j+1][i].v - v)*dhy       /* - v(dv/dy) */
708:                      -(Rd/p)*(Ts*(X[j+1][i].p - X[j-1][i].p)*0.5*dhy + p*0*(X[j+1][i].Ts - X[j-1][i].Ts)*0.5*dhy) /* -(R/p)[Ts(dp/dy)+ p(dTs/dy)] */
709: /*                     -(1/p)*(X[j+1][i].P - X[j-1][i].P)*dhy */
710:                      -f*u;

712:       /* dT/dt -> time change of temperature */
713:       Frhs[j][i].Ts = (fsfc1/(csoil*dzlay))                                            /* Fnet/(Cp*dz)  diabatic change in T */
714:                       -u_plus*(Ts - X[j][i-1].Ts)*dhx - u_minus*(X[j][i+1].Ts - Ts)*dhx  /* - u*(dTs/dx)  advection x */
715:                       -v_plus*(Ts - X[j-1][i].Ts)*dhy - v_minus*(X[j+1][i].Ts - Ts)*dhy  /* - v*(dTs/dy)  advection y */
716:                       + diffconst*((X[j][i+1].Ts - 2*Ts + X[j][i-1].Ts)*dhx*dhx               /* + D(Ts_xx + Ts_yy)  diffusion */
717:                                    + (X[j+1][i].Ts - 2*Ts + X[j-1][i].Ts)*dhy*dhy);

719:       /* dp/dt -> time change of */
720:       Frhs[j][i].p = -u_plus*(p - X[j][i-1].p)*dhx - u_minus*(X[j][i+1].p - p)*dhx     /* - u*(dp/dx) */
721:                      -v_plus*(p - X[j-1][i].p)*dhy - v_minus*(X[j+1][i].p - p)*dhy;    /* - v*(dp/dy) */

723:       Frhs[j][i].Ta = Ra/Cp;  /* dTa/dt time change of air temperature */
724:     }
725:   }

727:   /* Restore vectors */
728:   DMDAVecRestoreArrayRead(da,localT,&X);
729:   DMDAVecRestoreArray(da,F,&Frhs);
730:   DMRestoreLocalVector(da,&localT);
731:   return(0);
732: }

736: PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal time,Vec T,void *ctx)
737: {
738:   PetscErrorCode    ierr;
739:   const PetscScalar *array;
740:   MonitorCtx        *user  = (MonitorCtx*)ctx;
741:   PetscViewer       viewer = user->drawviewer;
742:   PetscMPIInt       rank;
743:   PetscReal         norm;

746:   MPI_Comm_rank(PetscObjectComm((PetscObject)ts),&rank);
747:   VecNorm(T,NORM_INFINITY,&norm);

749:   if (step%user->interval == 0) {
750:     VecGetArrayRead(T,&array);
751:     if (!rank) printf("step %4d, time %8.1f,  %6.4f, %6.4f, %6.4f, %6.4f, %6.4f, %6.4f\n",(int)step,(double)time,(double)(((array[0]-273)*9)/5 + 32),(double)(((array[1]-273)*9)/5 + 32),(double)array[2],(double)array[3],(double)array[4],(double)array[5]);
752:     VecRestoreArrayRead(T,&array);
753:   }

755:   if (user->drawcontours) {
756:     VecView(T,viewer);
757:   }
758:   return(0);
759: }