Actual source code: ex5.c
petsc-dev 2014-02-02
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 <petscdmda.h>
71: /* stefan-boltzmann constant */
72: #define SIG 0.000000056703
73: /* absorption-emission constant for surface */
74: #define EMMSFC 1
75: /* amount of time (seconds) that passes before new flux is calculated */
76: #define TIMESTEP 1
78: /* variables of interest to be solved at each grid point */
79: typedef struct {
80: PetscScalar Ts,Ta; /* surface and air temperature */
81: PetscScalar u,v; /* wind speed */
82: PetscScalar p; /* density */
83: } Field;
85: /* User defined variables. Used in solving for variables of interest */
86: typedef struct {
87: DM da; /* grid */
88: PetscScalar csoil; /* heat constant for layer */
89: PetscScalar dzlay; /* thickness of top soil layer */
90: PetscScalar emma; /* emission parameter */
91: PetscScalar wind; /* wind speed */
92: PetscScalar dewtemp; /* dew point temperature (moisture in air) */
93: PetscScalar pressure1; /* sea level pressure */
94: PetscScalar airtemp; /* temperature of air near boundary layer inversion */
95: PetscScalar Ts; /* temperature at the surface */
96: PetscScalar fract; /* fraction of sky covered by clouds */
97: PetscScalar Tc; /* temperature at base of lowest cloud layer */
98: PetscScalar lat; /* Latitude in degrees */
99: PetscScalar init; /* initialization scenario */
100: PetscScalar deep_grnd_temp; /* temperature of ground under top soil surface layer */
101: } AppCtx;
103: /* Struct for visualization */
104: typedef struct {
105: PetscBool drawcontours; /* flag - 1 indicates drawing contours */
106: PetscViewer drawviewer;
107: PetscInt interval;
108: } MonitorCtx;
111: /* Inputs read in from text file */
112: struct in {
113: PetscScalar Ts; /* surface temperature */
114: PetscScalar Td; /* dewpoint temperature */
115: PetscScalar Tc; /* temperature of cloud base */
116: PetscScalar fr; /* fraction of sky covered by clouds */
117: PetscScalar wnd; /* wind speed */
118: PetscScalar Ta; /* air temperature */
119: PetscScalar pwt; /* precipitable water */
120: PetscScalar wndDir; /* wind direction */
121: PetscScalar lat; /* latitude */
122: PetscReal time; /* time in hours */
123: PetscScalar init;
124: };
126: /* functions */
127: extern PetscScalar emission(PetscScalar); /* sets emission/absorption constant depending on water vapor content */
128: extern PetscScalar calc_q(PetscScalar); /* calculates specific humidity */
129: extern PetscScalar mph2mpers(PetscScalar); /* converts miles per hour to meters per second */
130: 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. */
131: extern PetscScalar fahr_to_cel(PetscScalar); /* converts Fahrenheit to Celsius */
132: extern PetscScalar cel_to_fahr(PetscScalar); /* converts Celsius to Fahrenheit */
133: extern PetscScalar calcmixingr(PetscScalar, PetscScalar); /* calculates mixing ratio */
134: extern PetscScalar cloud(PetscScalar); /* cloud radiative parameterization */
135: extern PetscErrorCode FormInitialSolution(DM,Vec,void*); /* Specifies initial conditions for the system of equations (PETSc defined function) */
136: extern PetscErrorCode RhsFunc(TS,PetscReal,Vec,Vec,void*); /* Specifies the user defined functions (PETSc defined function) */
137: extern PetscErrorCode Monitor(TS,PetscInt,PetscReal,Vec,void*); /* Specifies output and visualization tools (PETSc defined function) */
138: extern void readinput(struct in *put); /* reads input from text file */
139: extern PetscErrorCode calcfluxs(PetscScalar, PetscScalar, PetscScalar, PetscScalar, PetscScalar, PetscScalar*); /* calculates upward IR from surface */
140: extern PetscErrorCode calcfluxa(PetscScalar, PetscScalar, PetscScalar, PetscScalar*); /* calculates downward IR from atmosphere */
141: extern PetscErrorCode sensibleflux(PetscScalar, PetscScalar, PetscScalar, PetscScalar*); /* calculates sensible heat flux */
142: extern PetscErrorCode potential_temperature(PetscScalar, PetscScalar, PetscScalar, PetscScalar, PetscScalar*); /* calculates potential temperature */
143: extern PetscErrorCode latentflux(PetscScalar, PetscScalar, PetscScalar, PetscScalar, PetscScalar*); /* calculates latent heat flux */
144: extern PetscErrorCode calc_gflux(PetscScalar, PetscScalar, PetscScalar*); /* calculates flux between top soil layer and underlying earth */
148: int main(int argc,char **argv)
149: {
151: int time; /* amount of loops */
152: struct in put;
153: PetscScalar rh; /* relative humidity */
154: PetscScalar x; /* memory varialbe for relative humidity calculation */
155: PetscScalar deep_grnd_temp; /* temperature of ground under top soil surface layer */
156: PetscScalar emma; /* absorption-emission constant for air */
157: PetscScalar pressure1 = 101300; /* surface pressure */
158: PetscScalar mixratio; /* mixing ratio */
159: PetscScalar airtemp; /* temperature of air near boundary layer inversion */
160: PetscScalar dewtemp; /* dew point temperature */
161: PetscScalar sfctemp; /* temperature at surface */
162: PetscScalar pwat; /* total column precipitable water */
163: PetscScalar cloudTemp; /* temperature at base of cloud */
164: AppCtx user; /* user-defined work context */
165: MonitorCtx usermonitor; /* user-defined monitor context */
166: PetscMPIInt rank,size;
167: TS ts;
168: SNES snes;
169: DM da;
170: Vec T,rhs; /* solution vector */
171: Mat J; /* Jacobian matrix */
172: PetscReal ftime,dt;
173: PetscInt steps,dof = 5;
174: PetscBool use_coloring = PETSC_TRUE;
175: MatFDColoring matfdcoloring = 0;
176: PetscBool monitor_off = PETSC_FALSE;
178: PetscInitialize(&argc,&argv,(char*)0,help);
179: MPI_Comm_size(PETSC_COMM_WORLD,&size);
180: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
182: /* Inputs */
183: readinput(&put);
185: sfctemp = put.Ts;
186: dewtemp = put.Td;
187: cloudTemp = put.Tc;
188: airtemp = put.Ta;
189: pwat = put.pwt;
191: if (!rank) PetscPrintf(PETSC_COMM_SELF,"Initial Temperature = %g\n",sfctemp); /* input surface temperature */
193: deep_grnd_temp = sfctemp - 10; /* set underlying ground layer temperature */
194: emma = emission(pwat); /* accounts for radiative effects of water vapor */
196: /* Converts from Fahrenheit to Celsuis */
197: sfctemp = fahr_to_cel(sfctemp);
198: airtemp = fahr_to_cel(airtemp);
199: dewtemp = fahr_to_cel(dewtemp);
200: cloudTemp = fahr_to_cel(cloudTemp);
201: deep_grnd_temp = fahr_to_cel(deep_grnd_temp);
203: /* Converts from Celsius to Kelvin */
204: sfctemp += 273;
205: airtemp += 273;
206: dewtemp += 273;
207: cloudTemp += 273;
208: deep_grnd_temp += 273;
210: /* Calculates initial relative humidity */
211: x = calcmixingr(dewtemp,pressure1);
212: mixratio = calcmixingr(sfctemp,pressure1);
213: rh = (x/mixratio)*100;
215: if (!rank) printf("Initial RH = %.1f percent\n\n",(double)rh); /* prints initial relative humidity */
217: time = 3600*put.time; /* sets amount of timesteps to run model */
219: /* Configure PETSc TS solver */
220: /*------------------------------------------*/
222: /* Create grid */
223: DMDACreate2d(PETSC_COMM_WORLD,DMDA_BOUNDARY_PERIODIC,DMDA_BOUNDARY_PERIODIC,DMDA_STENCIL_STAR,-20,-20,
224: PETSC_DECIDE,PETSC_DECIDE,dof,1,NULL,NULL,&da);
225: DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0);
227: /* Define output window for each variable of interest */
228: DMDASetFieldName(da,0,"Ts");
229: DMDASetFieldName(da,1,"Ta");
230: DMDASetFieldName(da,2,"u");
231: DMDASetFieldName(da,3,"v");
232: DMDASetFieldName(da,4,"p");
234: /* set values for appctx */
235: user.da = da;
236: user.Ts = sfctemp;
237: user.fract = put.fr; /* fraction of sky covered by clouds */
238: user.dewtemp = dewtemp; /* dew point temperature (mositure in air) */
239: user.csoil = 2000000; /* heat constant for layer */
240: user.dzlay = 0.08; /* thickness of top soil layer */
241: user.emma = emma; /* emission parameter */
242: user.wind = put.wnd; /* wind spped */
243: user.pressure1 = pressure1; /* sea level pressure */
244: user.airtemp = airtemp; /* temperature of air near boundar layer inversion */
245: user.Tc = cloudTemp; /* temperature at base of lowest cloud layer */
246: user.init = put.init; /* user chosen initiation scenario */
247: user.lat = 70*0.0174532; /* converts latitude degrees to latitude in radians */
248: user.deep_grnd_temp = deep_grnd_temp; /* temp in lowest ground layer */
250: /* set values for MonitorCtx */
251: usermonitor.drawcontours = PETSC_FALSE;
252: PetscOptionsHasName(NULL,"-drawcontours",&usermonitor.drawcontours);
253: if (usermonitor.drawcontours) {
254: PetscReal bounds[] = {1000.0,-1000., -1000.,-1000., 1000.,-1000., 1000.,-1000., 1000,-1000, 100700,100800};
255: PetscViewerDrawOpen(PETSC_COMM_WORLD,0,0,0,0,300,300,&usermonitor.drawviewer);
256: PetscViewerDrawSetBounds(usermonitor.drawviewer,dof,bounds);
257: }
258: usermonitor.interval = 1;
259: PetscOptionsGetInt(NULL,"-monitor_interval",&usermonitor.interval,NULL);
261: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
262: Extract global vectors from DA;
263: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
264: DMCreateGlobalVector(da,&T);
265: VecDuplicate(T,&rhs); /* r: vector to put the computed right hand side */
267: TSCreate(PETSC_COMM_WORLD,&ts);
268: TSSetProblemType(ts,TS_NONLINEAR);
269: TSSetType(ts,TSBEULER);
270: TSSetRHSFunction(ts,rhs,RhsFunc,&user);
272: /* Set Jacobian evaluation routine - use coloring to compute finite difference Jacobian efficiently */
273: DMSetMatType(da,MATAIJ);
274: DMCreateMatrix(da,&J);
275: TSGetSNES(ts,&snes);
276: if (use_coloring) {
277: ISColoring iscoloring;
278: DMCreateColoring(da,IS_COLORING_GLOBAL,&iscoloring);
279: MatFDColoringCreate(J,iscoloring,&matfdcoloring);
280: MatFDColoringSetFromOptions(matfdcoloring);
281: MatFDColoringSetUp(J,iscoloring,matfdcoloring);
282: ISColoringDestroy(&iscoloring);
283: MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode (*)(void))SNESTSFormFunction,ts);
284: SNESSetJacobian(snes,J,J,SNESComputeJacobianDefaultColor,matfdcoloring);
285: } else {
286: SNESSetJacobian(snes,J,J,SNESComputeJacobianDefault,NULL);
287: }
289: /* Define what to print for ts_monitor option */
290: PetscOptionsHasName(NULL,"-monitor_off",&monitor_off);
291: if (!monitor_off) {
292: TSMonitorSet(ts,Monitor,&usermonitor,NULL);
293: }
294: FormInitialSolution(da,T,&user);
295: dt = TIMESTEP; /* initial time step */
296: ftime = TIMESTEP*time;
297: if (!rank) printf("time %d, ftime %g hour, TIMESTEP %g\n",time,(double)(ftime/3600),(double)dt);
299: TSSetInitialTimeStep(ts,0.0,dt);
300: TSSetDuration(ts,time,ftime);
301: TSSetSolution(ts,T);
302: TSSetDM(ts,da);
304: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
305: Set runtime options
306: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
307: TSSetFromOptions(ts);
309: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
310: Solve nonlinear system
311: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
312: TSSolve(ts,T);
313: TSGetSolveTime(ts,&ftime);
314: TSGetTimeStepNumber(ts,&steps);
315: if (!rank) PetscPrintf(PETSC_COMM_WORLD,"Solution T after %g hours %d steps\n",(double)(ftime/3600),steps);
318: if (matfdcoloring) {MatFDColoringDestroy(&matfdcoloring);}
319: if (usermonitor.drawcontours) {
320: PetscViewerDestroy(&usermonitor.drawviewer);
321: }
322: MatDestroy(&J);
323: VecDestroy(&T);
324: VecDestroy(&rhs);
325: TSDestroy(&ts);
326: DMDestroy(&da);
328: PetscFinalize();
329: return 0;
330: }
331: /*****************************end main program********************************/
332: /*****************************************************************************/
333: /*****************************************************************************/
334: /*****************************************************************************/
337: PetscErrorCode calcfluxs(PetscScalar sfctemp, PetscScalar airtemp, PetscScalar emma, PetscScalar fract, PetscScalar cloudTemp, PetscScalar *flux)
338: {
340: *flux = SIG*((EMMSFC*emma*PetscPowScalarInt(airtemp,4)) + (EMMSFC*fract*(1 - emma)*PetscPowScalarInt(cloudTemp,4)) - (EMMSFC*PetscPowScalarInt(sfctemp,4))); /* calculates flux using Stefan-Boltzmann relation */
341: return(0);
342: }
346: PetscErrorCode calcfluxa(PetscScalar sfctemp, PetscScalar airtemp, PetscScalar emma, PetscScalar *flux) /* this function is not currently called upon */
347: {
348: PetscScalar emm = 0.001;
351: *flux = SIG*(-emm*(PetscPowScalarInt(airtemp,4))); /* calculates flux usinge Stefan-Boltzmann relation */
352: return(0);
353: }
356: PetscErrorCode sensibleflux(PetscScalar sfctemp, PetscScalar airtemp, PetscScalar wind, PetscScalar *sheat)
357: {
358: PetscScalar density = 1; /* air density */
359: PetscScalar Cp = 1005; /* heat capicity for dry air */
360: PetscScalar wndmix; /* temperature change from wind mixing: wind*Ch */
363: wndmix = 0.0025 + 0.0042*wind; /* regression equation valid for neutral and stable BL */
364: *sheat = density*Cp*wndmix*(airtemp - sfctemp); /* calculates sensible heat flux */
365: return(0);
366: }
370: PetscErrorCode latentflux(PetscScalar sfctemp, PetscScalar dewtemp, PetscScalar wind, PetscScalar pressure1, PetscScalar *latentheat)
371: {
372: PetscScalar density = 1; /* density of dry air */
373: PetscScalar q; /* actual specific humitity */
374: PetscScalar qs; /* saturation specific humidity */
375: PetscScalar wndmix; /* temperature change from wind mixing: wind*Ch */
376: PetscScalar beta = .4; /* moisture availability */
377: PetscScalar mr; /* mixing ratio */
378: PetscScalar lhcnst; /* latent heat of vaporization constant = 2501000 J/kg at 0c */
379: /* latent heat of saturation const = 2834000 J/kg */
380: /* latent heat of fusion const = 333700 J/kg */
383: wind = mph2mpers(wind); /* converts wind from mph to meters per second */
384: wndmix = 0.0025 + 0.0042*wind; /* regression equation valid for neutral BL */
385: lhcnst = Lconst(sfctemp); /* calculates latent heat of evaporation */
386: mr = calcmixingr(sfctemp,pressure1); /* calculates saturation mixing ratio */
387: qs = calc_q(mr); /* calculates saturation specific humidty */
388: mr = calcmixingr(dewtemp,pressure1); /* calculates mixing ratio */
389: q = calc_q(mr); /* calculates specific humidty */
391: *latentheat = density*wndmix*beta*lhcnst*(q - qs); /* calculates latent heat flux */
392: return(0);
393: }
397: PetscErrorCode potential_temperature(PetscScalar temp, PetscScalar pressure1, PetscScalar pressure2, PetscScalar sfctemp, PetscScalar *pottemp)
398: {
399: PetscScalar kdry; /* poisson constant for dry atmosphere */
400: PetscScalar pavg; /* average atmospheric pressure */
401: /* PetscScalar mixratio; mixing ratio */
402: /* PetscScalar kmoist; poisson constant for moist atmosphere */
405: /* mixratio = calcmixingr(sfctemp,pressure1); */
407: /* initialize poisson constant */
408: kdry = 0.2854;
409: /* kmoist = 0.2854*(1 - 0.24*mixratio); */
411: pavg = ((0.7*pressure1)+pressure2)/2; /* calculates simple average press */
412: *pottemp = temp*(PetscPowScalar((pressure1/pavg),kdry)); /* calculates potential temperature */
413: return(0);
414: }
415: extern PetscScalar calcmixingr(PetscScalar dtemp, PetscScalar pressure1)
416: {
417: PetscScalar e; /* vapor pressure */
418: PetscScalar mixratio; /* mixing ratio */
420: dtemp = dtemp - 273; /* converts from Kelvin to Celsuis */
421: e = 6.11*(PetscPowScalar(10,((7.5*dtemp)/(237.7+dtemp)))); /* converts from dew point temp to vapor pressure */
422: e = e*100; /* converts from hPa to Pa */
423: mixratio = (0.622*e)/(pressure1 - e); /* computes mixing ratio */
424: mixratio = mixratio*1; /* convert to g/Kg */
426: return mixratio;
427: }
428: extern PetscScalar calc_q(PetscScalar rv)
429: {
430: PetscScalar specific_humidity; /* define specific humidity variable */
431: specific_humidity = rv/(1 + rv); /* calculates specific humidity */
432: return specific_humidity;
433: }
437: PetscErrorCode calc_gflux(PetscScalar sfctemp, PetscScalar deep_grnd_temp, PetscScalar* Gflux)
438: {
439: PetscScalar k; /* thermal conductivity parameter */
440: PetscScalar n = 0.38; /* value of soil porosity */
441: PetscScalar dz = 1; /* depth of layer between soil surface and deep soil layer */
442: PetscScalar unit_soil_weight = 2700; /* unit soil weight in kg/m^3 */
445: k = ((0.135*(1-n)*unit_soil_weight) + 64.7)/(unit_soil_weight - (0.947*(1-n)*unit_soil_weight)); /* dry soil conductivity */
446: *Gflux = (k*(deep_grnd_temp - sfctemp)/dz); /* calculates flux from deep ground layer */
447: return(0);
448: }
451: extern PetscScalar emission(PetscScalar pwat)
452: {
453: PetscScalar emma;
455: emma = 0.725 + 0.17*log10(pwat);
457: return emma;
458: }
459: extern PetscScalar cloud(PetscScalar fract)
460: {
461: PetscScalar emma = 0;
463: /* modifies radiative balance depending on cloud cover */
464: if (fract >= 0.9) emma = 1;
465: else if (0.9 > fract && fract >= 0.8) emma = 0.9;
466: else if (0.8 > fract && fract >= 0.7) emma = 0.85;
467: else if (0.7 > fract && fract >= 0.6) emma = 0.75;
468: else if (0.6 > fract && fract >= 0.5) emma = 0.65;
469: else if (0.4 > fract && fract >= 0.3) emma = emma*1.086956;
470: return emma;
471: }
472: extern PetscScalar Lconst(PetscScalar sfctemp)
473: {
474: PetscScalar Lheat;
475: sfctemp -=273; /* converts from kelvin to celsius */
476: Lheat = 4186.8*(597.31 - 0.5625*sfctemp); /* calculates latent heat constant */
477: return Lheat;
478: }
479: extern PetscScalar mph2mpers(PetscScalar wind)
480: {
481: wind = ((wind*1.6*1000)/3600); /* converts wind from mph to meters per second */
482: return wind;
483: }
484: extern PetscScalar fahr_to_cel(PetscScalar temp)
485: {
486: temp = (5*(temp-32))/9; /* converts from farhrenheit to celsuis */
487: return temp;
488: }
489: extern PetscScalar cel_to_fahr(PetscScalar temp)
490: {
491: temp = ((temp*9)/5) + 32; /* converts from celsuis to farhrenheit */
492: return temp;
493: }
494: void readinput(struct in *put)
495: {
496: int i;
497: char x;
498: FILE *ifp;
500: ifp = fopen("ex5_control.txt", "r");
502: for (i=0; i<110; i++) fscanf(ifp, "%c", &x);
503: fscanf(ifp, "%lf", &put->Ts);
505: for (i=0; i<43; i++) fscanf(ifp, "%c", &x);
506: fscanf(ifp, "%lf", &put->Td);
508: for (i=0; i<43; i++) fscanf(ifp, "%c", &x);
509: fscanf(ifp, "%lf", &put->Ta);
511: for (i=0; i<43; i++) fscanf(ifp, "%c", &x);
512: fscanf(ifp, "%lf", &put->Tc);
514: for (i=0; i<43; i++) fscanf(ifp, "%c", &x);
515: fscanf(ifp, "%lf", &put->fr);
517: for (i=0; i<43; i++) fscanf(ifp, "%c", &x);
518: fscanf(ifp, "%lf", &put->wnd);
520: for (i=0; i<43; i++) fscanf(ifp, "%c", &x);
521: fscanf(ifp, "%lf", &put->pwt);
523: for (i=0; i<43; i++) fscanf(ifp, "%c", &x);
524: fscanf(ifp, "%lf", &put->wndDir);
526: for (i=0; i<43; i++) fscanf(ifp, "%c", &x);
527: fscanf(ifp, "%lf", &put->time);
529: for (i=0; i<63; i++) fscanf(ifp, "%c", &x);
530: fscanf(ifp, "%lf", &put->init);
531: }
533: /* ------------------------------------------------------------------- */
536: PetscErrorCode FormInitialSolution(DM da,Vec Xglobal,void *ctx)
537: {
539: AppCtx *user = (AppCtx*)ctx; /* user-defined application context */
540: PetscInt i,j,xs,ys,xm,ym,Mx,My;
541: Field **X;
544: DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
545: PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
547: /* Get pointers to vector data */
548: DMDAVecGetArray(da,Xglobal,&X);
550: /* Get local grid boundaries */
551: DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);
553: /* Compute function over the locally owned part of the grid */
555: if (user->init == 1) {
556: for (j=ys; j<ys+ym; j++) {
557: for (i=xs; i<xs+xm; i++) {
558: X[j][i].Ts = user->Ts - i*0.0001;
559: X[j][i].Ta = X[j][i].Ts - 5;
560: X[j][i].u = 0;
561: X[j][i].v = 0;
562: X[j][i].p = 1.25;
563: if ((j == 5 || j == 6) && (i == 4 || i == 5)) X[j][i].p += 0.00001;
564: if ((j == 5 || j == 6) && (i == 12 || i == 13)) X[j][i].p += 0.00001;
565: }
566: }
567: } else {
568: for (j=ys; j<ys+ym; j++) {
569: for (i=xs; i<xs+xm; i++) {
570: X[j][i].Ts = user->Ts;
571: X[j][i].Ta = X[j][i].Ts - 5;
572: X[j][i].u = 0;
573: X[j][i].v = 0;
574: X[j][i].p = 1.25;
575: }
576: }
577: }
579: /* Restore vectors */
580: DMDAVecRestoreArray(da,Xglobal,&X);
581: return(0);
582: }
586: /*
587: RhsFunc - Evaluates nonlinear function F(u).
589: Input Parameters:
590: . ts - the TS context
591: . t - current time
592: . Xglobal - input vector
593: . F - output vector
594: . ptr - optional user-defined context, as set by SNESSetFunction()
596: Output Parameter:
597: . F - rhs function vector
598: */
599: PetscErrorCode RhsFunc(TS ts,PetscReal t,Vec Xglobal,Vec F,void *ctx)
600: {
601: AppCtx *user = (AppCtx*)ctx; /* user-defined application context */
602: DM da = user->da;
604: PetscInt i,j,Mx,My,xs,ys,xm,ym;
605: PetscReal dhx,dhy;
606: Vec localT;
607: Field **X,**Frhs; /* structures that contain variables of interest and left hand side of governing equations respectively */
608: PetscScalar csoil = user->csoil; /* heat constant for layer */
609: PetscScalar dzlay = user->dzlay; /* thickness of top soil layer */
610: PetscScalar emma = user->emma; /* emission parameter */
611: PetscScalar wind = user->wind; /* wind speed */
612: PetscScalar dewtemp = user->dewtemp; /* dew point temperature (moisture in air) */
613: PetscScalar pressure1 = user->pressure1; /* sea level pressure */
614: PetscScalar airtemp = user->airtemp; /* temperature of air near boundary layer inversion */
615: PetscScalar fract = user->fract; /* fraction of the sky covered by clouds */
616: PetscScalar Tc = user->Tc; /* temperature at base of lowest cloud layer */
617: PetscScalar lat = user->lat; /* latitude */
618: PetscScalar Cp = 1005.7; /* specific heat of air at constant pressure */
619: PetscScalar Rd = 287.058; /* gas constant for dry air */
620: PetscScalar diffconst = 1000; /* diffusion coefficient */
621: PetscScalar f = 2*0.0000727*PetscSinScalar(lat); /* coriolis force */
622: PetscScalar deep_grnd_temp = user->deep_grnd_temp; /* temp in lowest ground layer */
623: PetscScalar Ts,u,v,p;
624: PetscScalar u_abs,u_plus,u_minus,v_abs,v_plus,v_minus;
626: PetscScalar sfctemp1,fsfc1,Ra;
627: PetscScalar sheat; /* sensible heat flux */
628: PetscScalar latentheat; /* latent heat flux */
629: PetscScalar groundflux; /* flux from conduction of deep ground layer in contact with top soil */
630: PetscInt xend,yend;
633: DMGetLocalVector(da,&localT);
634: DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
635: PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
637: dhx = (PetscReal)(Mx-1)/(5000*(Mx-1)); /* dhx = 1/dx; assume 2D space domain: [0.0, 1.e5] x [0.0, 1.e5] */
638: dhy = (PetscReal)(My-1)/(5000*(Mx-1)); /* dhy = 1/dy; */
641: /*
642: Scatter ghost points to local vector,using the 2-step process
643: DAGlobalToLocalBegin(),DAGlobalToLocalEnd().
644: By placing code between these two statements, computations can be
645: done while messages are in transition.
646: */
647: DMGlobalToLocalBegin(da,Xglobal,INSERT_VALUES,localT);
648: DMGlobalToLocalEnd(da,Xglobal,INSERT_VALUES,localT);
650: /* Get pointers to vector data */
651: DMDAVecGetArray(da,localT,&X);
652: DMDAVecGetArray(da,F,&Frhs);
654: /* Get local grid boundaries */
655: DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);
657: /* Compute function over the locally owned part of the grid */
658: /* the interior points */
659: xend=xs+xm; yend=ys+ym;
660: for (j=ys; j<yend; j++) {
661: for (i=xs; i<xend; i++) {
662: 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; */
664: sfctemp1 = (double)Ts;
665: sfctemp1 = (double)X[j][i].Ts;
666: calcfluxs(sfctemp1,airtemp,emma,fract,Tc,&fsfc1); /* calculates surface net radiative flux */
667: sensibleflux(sfctemp1,airtemp,wind,&sheat); /* calculate sensible heat flux */
668: latentflux(sfctemp1,dewtemp,wind,pressure1,&latentheat); /* calculates latent heat flux */
669: calc_gflux(sfctemp1,deep_grnd_temp,&groundflux); /* calculates flux from earth below surface soil layer by conduction */
670: calcfluxa(sfctemp1,airtemp,emma,&Ra); /* Calculates the change in downward radiative flux */
671: fsfc1 = fsfc1 + latentheat + sheat + groundflux; /* adds radiative, sensible heat, latent heat, and ground heat flux yielding net flux */
673: /* convective coefficients for upwinding */
674: u_abs = PetscAbsScalar(u);
675: u_plus = .5*(u + u_abs); /* u if u>0; 0 if u<0 */
676: u_minus = .5*(u - u_abs); /* u if u <0; 0 if u>0 */
678: v_abs = PetscAbsScalar(v);
679: v_plus = .5*(v + v_abs); /* v if v>0; 0 if v<0 */
680: v_minus = .5*(v - v_abs); /* v if v <0; 0 if v>0 */
682: /* Solve governing equations */
683: /* P = p*Rd*Ts; */
685: /* du/dt -> time change of east-west component of the wind */
686: 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) */
687: - v_plus*(u - X[j-1][i].u)*dhy - v_minus*(X[j+1][i].u - u)*dhy /* - v(du/dy) */
688: -(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)] */
689: /* -(1/p)*(X[j][i+1].P - X[j][i-1].P)*dhx */
690: + f*v;
692: /* dv/dt -> time change of north-south component of the wind */
693: 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) */
694: - v_plus*(v - X[j-1][i].v)*dhy - v_minus*(X[j+1][i].v - v)*dhy /* - v(dv/dy) */
695: -(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)] */
696: /* -(1/p)*(X[j+1][i].P - X[j-1][i].P)*dhy */
697: -f*u;
699: /* dT/dt -> time change of temperature */
700: Frhs[j][i].Ts = (fsfc1/(csoil*dzlay)) /* Fnet/(Cp*dz) diabatic change in T */
701: -u_plus*(Ts - X[j][i-1].Ts)*dhx - u_minus*(X[j][i+1].Ts - Ts)*dhx /* - u*(dTs/dx) advection x */
702: -v_plus*(Ts - X[j-1][i].Ts)*dhy - v_minus*(X[j+1][i].Ts - Ts)*dhy /* - v*(dTs/dy) advection y */
703: + diffconst*((X[j][i+1].Ts - 2*Ts + X[j][i-1].Ts)*dhx*dhx /* + D(Ts_xx + Ts_yy) diffusion */
704: + (X[j+1][i].Ts - 2*Ts + X[j-1][i].Ts)*dhy*dhy);
706: /* dp/dt -> time change of */
707: 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) */
708: -v_plus*(p - X[j-1][i].p)*dhy - v_minus*(X[j+1][i].p - p)*dhy; /* - v*(dp/dy) */
710: Frhs[j][i].Ta = Ra/Cp; /* dTa/dt time change of air temperature */
711: }
712: }
714: /* Restore vectors */
715: DMDAVecRestoreArray(da,localT,&X);
716: DMDAVecRestoreArray(da,F,&Frhs);
717: DMRestoreLocalVector(da,&localT);
718: return(0);
719: }
723: PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal time,Vec T,void *ctx)
724: {
726: PetscScalar *array;
727: MonitorCtx *user = (MonitorCtx*)ctx;
728: PetscViewer viewer = user->drawviewer;
729: PetscMPIInt rank;
730: PetscReal norm;
733: MPI_Comm_rank(PetscObjectComm((PetscObject)ts),&rank);
734: VecNorm(T,NORM_INFINITY,&norm);
736: if (step%user->interval == 0) {
737: VecGetArray(T,&array);
738: 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]);
739: VecRestoreArray(T,&array);
740: }
742: if (user->drawcontours) {
743: VecView(T,viewer);
744: }
745: return(0);
746: }