Actual source code: ex5.c
petsc-3.4.2 2013-07-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: } MonitorCtx;
110: /* Inputs read in from text file */
111: struct in {
112: PetscScalar Ts; /* surface temperature */
113: PetscScalar Td; /* dewpoint temperature */
114: PetscScalar Tc; /* temperature of cloud base */
115: PetscScalar fr; /* fraction of sky covered by clouds */
116: PetscScalar wnd; /* wind speed */
117: PetscScalar Ta; /* air temperature */
118: PetscScalar pwt; /* precipitable water */
119: PetscScalar wndDir; /* wind direction */
120: PetscScalar lat; /* latitude */
121: PetscReal time; /* time in hours */
122: PetscScalar init;
123: };
125: /* functions */
126: extern PetscScalar emission(PetscScalar); /* sets emission/absorption constant depending on water vapor content */
127: extern PetscScalar calc_q(PetscScalar); /* calculates specific humidity */
128: extern PetscScalar mph2mpers(PetscScalar); /* converts miles per hour to meters per second */
129: 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. */
130: extern PetscScalar fahr_to_cel(PetscScalar); /* converts Fahrenheit to Celsius */
131: extern PetscScalar cel_to_fahr(PetscScalar); /* converts Celsius to Fahrenheit */
132: extern PetscScalar calcmixingr(PetscScalar, PetscScalar); /* calculates mixing ratio */
133: extern PetscScalar cloud(PetscScalar); /* cloud radiative parameterization */
134: extern PetscErrorCode FormInitialSolution(DM,Vec,void*); /* Specifies initial conditions for the system of equations (PETSc defined function) */
135: extern PetscErrorCode RhsFunc(TS,PetscReal,Vec,Vec,void*); /* Specifies the user defined functions (PETSc defined function) */
136: extern PetscErrorCode Monitor(TS,PetscInt,PetscReal,Vec,void*); /* Specifies output and visualization tools (PETSc defined function) */
137: extern void readinput(struct in *put); /* reads input from text file */
138: extern PetscErrorCode calcfluxs(PetscScalar, PetscScalar, PetscScalar, PetscScalar, PetscScalar, PetscScalar*); /* calculates upward IR from surface */
139: extern PetscErrorCode calcfluxa(PetscScalar, PetscScalar, PetscScalar, PetscScalar*); /* calculates downward IR from atmosphere */
140: extern PetscErrorCode sensibleflux(PetscScalar, PetscScalar, PetscScalar, PetscScalar*); /* calculates sensible heat flux */
141: extern PetscErrorCode potential_temperature(PetscScalar, PetscScalar, PetscScalar, PetscScalar, PetscScalar*); /* calculates potential temperature */
142: extern PetscErrorCode latentflux(PetscScalar, PetscScalar, PetscScalar, PetscScalar, PetscScalar*); /* calculates latent heat flux */
143: extern PetscErrorCode calc_gflux(PetscScalar, PetscScalar, PetscScalar*); /* calculates flux between top soil layer and underlying earth */
147: int main(int argc,char **argv)
148: {
150: int time; /* amount of loops */
151: struct in put;
152: PetscScalar rh; /* relative humidity */
153: PetscScalar x; /* memory varialbe for relative humidity calculation */
154: PetscScalar deep_grnd_temp; /* temperature of ground under top soil surface layer */
155: PetscScalar emma; /* absorption-emission constant for air */
156: PetscScalar pressure1 = 101300; /* surface pressure */
157: PetscScalar mixratio; /* mixing ratio */
158: PetscScalar airtemp; /* temperature of air near boundary layer inversion */
159: PetscScalar dewtemp; /* dew point temperature */
160: PetscScalar sfctemp; /* temperature at surface */
161: PetscScalar pwat; /* total column precipitable water */
162: PetscScalar cloudTemp; /* temperature at base of cloud */
163: AppCtx user; /* user-defined work context */
164: MonitorCtx usermonitor; /* user-defined monitor context */
165: PetscMPIInt rank,size;
166: TS ts;
167: SNES snes;
168: DM da;
169: Vec T,rhs; /* solution vector */
170: Mat J; /* Jacobian matrix */
171: PetscReal ftime,dt;
172: PetscInt steps,dof = 5;
174: PetscInitialize(&argc,&argv,(char*)0,help);
175: MPI_Comm_size(PETSC_COMM_WORLD,&size);
176: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
178: /* Inputs */
179: readinput(&put);
181: sfctemp = put.Ts;
182: dewtemp = put.Td;
183: cloudTemp = put.Tc;
184: airtemp = put.Ta;
185: pwat = put.pwt;
187: if (!rank) PetscPrintf(PETSC_COMM_SELF,"Initial Temperature = %g\n",sfctemp); /* input surface temperature */
189: deep_grnd_temp = sfctemp - 10; /* set underlying ground layer temperature */
190: emma = emission(pwat); /* accounts for radiative effects of water vapor */
192: /* Converts from Fahrenheit to Celsuis */
193: sfctemp = fahr_to_cel(sfctemp);
194: airtemp = fahr_to_cel(airtemp);
195: dewtemp = fahr_to_cel(dewtemp);
196: cloudTemp = fahr_to_cel(cloudTemp);
197: deep_grnd_temp = fahr_to_cel(deep_grnd_temp);
199: /* Converts from Celsius to Kelvin */
200: sfctemp += 273;
201: airtemp += 273;
202: dewtemp += 273;
203: cloudTemp += 273;
204: deep_grnd_temp += 273;
206: /* Calculates initial relative humidity */
207: x = calcmixingr(dewtemp,pressure1);
208: mixratio = calcmixingr(sfctemp,pressure1);
209: rh = (x/mixratio)*100;
211: if (!rank) printf("Initial RH = %.1f percent\n\n",rh); /* prints initial relative humidity */
213: time = 3600*put.time; /* sets amount of timesteps to run model */
215: /* Configure PETSc TS solver */
216: /*------------------------------------------*/
218: /* Create grid */
219: DMDACreate2d(PETSC_COMM_WORLD,DMDA_BOUNDARY_PERIODIC,DMDA_BOUNDARY_PERIODIC,DMDA_STENCIL_STAR,-20,-20,
220: PETSC_DECIDE,PETSC_DECIDE,dof,1,NULL,NULL,&da);
221: DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0);
223: /* Define output window for each variable of interest */
224: DMDASetFieldName(da,0,"Ts");
225: DMDASetFieldName(da,1,"Ta");
226: DMDASetFieldName(da,2,"u");
227: DMDASetFieldName(da,3,"v");
228: DMDASetFieldName(da,4,"p");
230: /* set values for appctx */
231: user.da = da;
232: user.Ts = sfctemp;
233: user.fract = put.fr; /* fraction of sky covered by clouds */
234: user.dewtemp = dewtemp; /* dew point temperature (mositure in air) */
235: user.csoil = 2000000; /* heat constant for layer */
236: user.dzlay = 0.08; /* thickness of top soil layer */
237: user.emma = emma; /* emission parameter */
238: user.wind = put.wnd; /* wind spped */
239: user.pressure1 = pressure1; /* sea level pressure */
240: user.airtemp = airtemp; /* temperature of air near boundar layer inversion */
241: user.Tc = cloudTemp; /* temperature at base of lowest cloud layer */
242: user.init = put.init; /* user chosen initiation scenario */
243: user.lat = 70*0.0174532; /* converts latitude degrees to latitude in radians */
244: user.deep_grnd_temp = deep_grnd_temp; /* temp in lowest ground layer */
246: /* set values for MonitorCtx */
247: usermonitor.drawcontours = PETSC_FALSE;
248: PetscOptionsHasName(NULL,"-drawcontours",&usermonitor.drawcontours);
249: if (usermonitor.drawcontours) {
250: PetscReal bounds[] = {1000.0,-1000., -1000.,-1000., 1000.,-1000., 1000.,-1000., 1000,-1000, 100700,100800};
251: PetscViewerDrawOpen(PETSC_COMM_WORLD,0,0,0,0,300,300,&usermonitor.drawviewer);
252: PetscViewerDrawSetBounds(usermonitor.drawviewer,dof,bounds);
253: }
255: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
256: Extract global vectors from DA;
257: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
258: DMCreateGlobalVector(da,&T);
259: VecDuplicate(T,&rhs); /* r: vector to put the computed right hand side */
261: TSCreate(PETSC_COMM_WORLD,&ts);
262: TSSetProblemType(ts,TS_NONLINEAR);
263: TSSetType(ts,TSBEULER);
264: TSSetRHSFunction(ts,rhs,RhsFunc,&user);
266: /* Set Jacobian evaluation routine - use coloring to compute finite difference Jacobian efficiently */
267: PetscBool use_coloring = PETSC_TRUE;
268: MatFDColoring matfdcoloring = 0;
269: DMCreateMatrix(da,MATAIJ,&J);
270: TSGetSNES(ts,&snes);
271: if (use_coloring) {
272: ISColoring iscoloring;
273: DMCreateColoring(da,IS_COLORING_GLOBAL,MATAIJ,&iscoloring);
274: MatFDColoringCreate(J,iscoloring,&matfdcoloring);
275: MatFDColoringSetFromOptions(matfdcoloring);
276: ISColoringDestroy(&iscoloring);
277: MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode (*)(void))SNESTSFormFunction,ts);
278: SNESSetJacobian(snes,J,J,SNESComputeJacobianDefaultColor,matfdcoloring);
279: } else {
280: SNESSetJacobian(snes,J,J,SNESComputeJacobianDefault,NULL);
281: }
283: /* Define what to print for ts_monitor option */
284: PetscBool monitor_off = PETSC_FALSE;
285: PetscOptionsHasName(NULL,"-monitor_off",&monitor_off);
286: if (!monitor_off) {
287: TSMonitorSet(ts,Monitor,&usermonitor,NULL);
288: }
289: FormInitialSolution(da,T,&user);
290: dt = TIMESTEP; /* initial time step */
291: ftime = TIMESTEP*time;
292: if (!rank) printf("time %d, ftime %g hour, TIMESTEP %g\n",time,ftime/3600,dt);
294: TSSetInitialTimeStep(ts,0.0,dt);
295: TSSetDuration(ts,time,ftime);
296: TSSetSolution(ts,T);
297: TSSetDM(ts,da);
299: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
300: Set runtime options
301: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
302: TSSetFromOptions(ts);
304: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
305: Solve nonlinear system
306: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
307: TSSolve(ts,T);
308: TSGetSolveTime(ts,&ftime);
309: TSGetTimeStepNumber(ts,&steps);
310: if (!rank) PetscPrintf(PETSC_COMM_WORLD,"Solution T after %g hours %d steps\n",ftime/3600,steps);
313: if (matfdcoloring) {MatFDColoringDestroy(&matfdcoloring);}
314: if (usermonitor.drawcontours) {
315: PetscViewerDestroy(&usermonitor.drawviewer);
316: }
317: MatDestroy(&J);
318: VecDestroy(&T);
319: VecDestroy(&rhs);
320: TSDestroy(&ts);
321: DMDestroy(&da);
323: PetscFinalize();
324: return 0;
325: }
326: /*****************************end main program********************************/
327: /*****************************************************************************/
328: /*****************************************************************************/
329: /*****************************************************************************/
332: PetscErrorCode calcfluxs(PetscScalar sfctemp, PetscScalar airtemp, PetscScalar emma, PetscScalar fract, PetscScalar cloudTemp, PetscScalar *flux)
333: {
335: *flux = SIG*((EMMSFC*emma*pow(airtemp,4)) + (EMMSFC*fract*(1 - emma)*pow(cloudTemp,4)) - (EMMSFC*pow(sfctemp,4))); /* calculates flux using Stefan-Boltzmann relation */
336: return(0);
337: }
341: PetscErrorCode calcfluxa(PetscScalar sfctemp, PetscScalar airtemp, PetscScalar emma, PetscScalar *flux) /* this function is not currently called upon */
342: {
343: PetscScalar emm = 0.001;
346: *flux = SIG*(-emm*(pow(airtemp,4))); /* calculates flux usinge Stefan-Boltzmann relation */
347: return(0);
348: }
351: PetscErrorCode sensibleflux(PetscScalar sfctemp, PetscScalar airtemp, PetscScalar wind, PetscScalar *sheat)
352: {
353: PetscScalar density = 1; /* air density */
354: PetscScalar Cp = 1005; /* heat capicity for dry air */
355: PetscScalar wndmix; /* temperature change from wind mixing: wind*Ch */
358: wndmix = 0.0025 + 0.0042*wind; /* regression equation valid for neutral and stable BL */
359: *sheat = density*Cp*wndmix*(airtemp - sfctemp); /* calculates sensible heat flux */
360: return(0);
361: }
365: PetscErrorCode latentflux(PetscScalar sfctemp, PetscScalar dewtemp, PetscScalar wind, PetscScalar pressure1, PetscScalar *latentheat)
366: {
367: PetscScalar density = 1; /* density of dry air */
368: PetscScalar q; /* actual specific humitity */
369: PetscScalar qs; /* saturation specific humidity */
370: PetscScalar wndmix; /* temperature change from wind mixing: wind*Ch */
371: PetscScalar beta = .4; /* moisture availability */
372: PetscScalar mr; /* mixing ratio */
373: PetscScalar lhcnst; /* latent heat of vaporization constant = 2501000 J/kg at 0c */
374: /* latent heat of saturation const = 2834000 J/kg */
375: /* latent heat of fusion const = 333700 J/kg */
378: wind = mph2mpers(wind); /* converts wind from mph to meters per second */
379: wndmix = 0.0025 + 0.0042*wind; /* regression equation valid for neutral BL */
380: lhcnst = Lconst(sfctemp); /* calculates latent heat of evaporation */
381: mr = calcmixingr(sfctemp,pressure1); /* calculates saturation mixing ratio */
382: qs = calc_q(mr); /* calculates saturation specific humidty */
383: mr = calcmixingr(dewtemp,pressure1); /* calculates mixing ratio */
384: q = calc_q(mr); /* calculates specific humidty */
386: *latentheat = density*wndmix*beta*lhcnst*(q - qs); /* calculates latent heat flux */
387: return(0);
388: }
392: PetscErrorCode potential_temperature(PetscScalar temp, PetscScalar pressure1, PetscScalar pressure2, PetscScalar sfctemp, PetscScalar *pottemp)
393: {
394: PetscScalar kdry; /* poisson constant for dry atmosphere */
395: PetscScalar kmoist; /* poisson constant for moist atmosphere */
396: PetscScalar pavg; /* average atmospheric pressure */
397: PetscScalar mixratio; /* mixing ratio */
400: mixratio = calcmixingr(sfctemp,pressure1);
402: /* initialize poisson constant */
403: kdry = 0.2854;
404: kmoist = 0.2854*(1 - 0.24*mixratio);
406: pavg = ((0.7*pressure1)+pressure2)/2; /* calculates simple average press */
407: *pottemp = temp*(pow((pressure1/pavg),kdry)); /* calculates potential temperature */
408: return(0);
409: }
410: extern PetscScalar calcmixingr(PetscScalar dtemp, PetscScalar pressure1)
411: {
412: PetscScalar e; /* vapor pressure */
413: PetscScalar mixratio; /* mixing ratio */
415: dtemp = dtemp - 273; /* converts from Kelvin to Celsuis */
416: e = 6.11*(pow(10,((7.5*dtemp)/(237.7+dtemp)))); /* converts from dew point temp to vapor pressure */
417: e = e*100; /* converts from hPa to Pa */
418: mixratio = (0.622*e)/(pressure1 - e); /* computes mixing ratio */
419: mixratio = mixratio*1; /* convert to g/Kg */
421: return mixratio;
422: }
423: extern PetscScalar calc_q(PetscScalar rv)
424: {
425: PetscScalar specific_humidity; /* define specific humidity variable */
426: specific_humidity = rv/(1 + rv); /* calculates specific humidity */
427: return specific_humidity;
428: }
432: PetscErrorCode calc_gflux(PetscScalar sfctemp, PetscScalar deep_grnd_temp, PetscScalar* Gflux)
433: {
434: PetscScalar k; /* thermal conductivity parameter */
435: PetscScalar n = 0.38; /* value of soil porosity */
436: PetscScalar dz = 1; /* depth of layer between soil surface and deep soil layer */
437: PetscScalar unit_soil_weight = 2700; /* unit soil weight in kg/m^3 */
440: k = ((0.135*(1-n)*unit_soil_weight) + 64.7)/(unit_soil_weight - (0.947*(1-n)*unit_soil_weight)); /* dry soil conductivity */
441: *Gflux = (k*(deep_grnd_temp - sfctemp)/dz); /* calculates flux from deep ground layer */
442: return(0);
443: }
446: extern PetscScalar emission(PetscScalar pwat)
447: {
448: PetscScalar emma;
450: emma = 0.725 + 0.17*log10(pwat);
452: return emma;
453: }
454: extern PetscScalar cloud(PetscScalar fract)
455: {
456: PetscScalar emma = 0;
458: /* modifies radiative balance depending on cloud cover */
459: if (fract >= 0.9) emma = 1;
460: else if (0.9 > fract && fract >= 0.8) emma = 0.9;
461: else if (0.8 > fract && fract >= 0.7) emma = 0.85;
462: else if (0.7 > fract && fract >= 0.6) emma = 0.75;
463: else if (0.6 > fract && fract >= 0.5) emma = 0.65;
464: else if (0.4 > fract && fract >= 0.3) emma = emma*1.086956;
465: return emma;
466: }
467: extern PetscScalar Lconst(PetscScalar sfctemp)
468: {
469: PetscScalar Lheat;
470: sfctemp -=273; /* converts from kelvin to celsius */
471: Lheat = 4186.8*(597.31 - 0.5625*sfctemp); /* calculates latent heat constant */
472: return Lheat;
473: }
474: extern PetscScalar mph2mpers(PetscScalar wind)
475: {
476: wind = ((wind*1.6*1000)/3600); /* converts wind from mph to meters per second */
477: return wind;
478: }
479: extern PetscScalar fahr_to_cel(PetscScalar temp)
480: {
481: temp = (5*(temp-32))/9; /* converts from farhrenheit to celsuis */
482: return temp;
483: }
484: extern PetscScalar cel_to_fahr(PetscScalar temp)
485: {
486: temp = ((temp*9)/5) + 32; /* converts from celsuis to farhrenheit */
487: return temp;
488: }
489: void readinput(struct in *put)
490: {
491: int i;
492: char x;
493: FILE *ifp;
495: ifp = fopen("ex5_control.txt", "r");
497: for (i=0; i<110; i++) fscanf(ifp, "%c", &x);
498: fscanf(ifp, "%lf", &put->Ts);
500: for (i=0; i<43; i++) fscanf(ifp, "%c", &x);
501: fscanf(ifp, "%lf", &put->Td);
503: for (i=0; i<43; i++) fscanf(ifp, "%c", &x);
504: fscanf(ifp, "%lf", &put->Ta);
506: for (i=0; i<43; i++) fscanf(ifp, "%c", &x);
507: fscanf(ifp, "%lf", &put->Tc);
509: for (i=0; i<43; i++) fscanf(ifp, "%c", &x);
510: fscanf(ifp, "%lf", &put->fr);
512: for (i=0; i<43; i++) fscanf(ifp, "%c", &x);
513: fscanf(ifp, "%lf", &put->wnd);
515: for (i=0; i<43; i++) fscanf(ifp, "%c", &x);
516: fscanf(ifp, "%lf", &put->pwt);
518: for (i=0; i<43; i++) fscanf(ifp, "%c", &x);
519: fscanf(ifp, "%lf", &put->wndDir);
521: for (i=0; i<43; i++) fscanf(ifp, "%c", &x);
522: fscanf(ifp, "%lf", &put->time);
524: for (i=0; i<63; i++) fscanf(ifp, "%c", &x);
525: fscanf(ifp, "%lf", &put->init);
526: }
528: /* ------------------------------------------------------------------- */
531: PetscErrorCode FormInitialSolution(DM da,Vec Xglobal,void *ctx)
532: {
534: AppCtx *user = (AppCtx*)ctx; /* user-defined application context */
535: PetscInt i,j,xs,ys,xm,ym,Mx,My;
536: Field **X;
537: PetscScalar deltT;
538: PetscReal hx,hy;
539: FILE *ifp;
540: FILE *ofp;
543: ofp = fopen("swing", "w");
544: ifp = fopen("grid.in", "r");
545: deltT = 0.8;
547: DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
548: PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
550: hx = 1/(PetscReal)(Mx-1);
551: hy = 1/(PetscReal)(My-1);
553: /* Get pointers to vector data */
554: DMDAVecGetArray(da,Xglobal,&X);
556: /* Get local grid boundaries */
557: DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);
559: /* Compute function over the locally owned part of the grid */
561: if (user->init == 1) {
562: for (j=ys; j<ys+ym; j++) {
563: for (i=xs; i<xs+xm; i++) {
564: X[j][i].Ts = user->Ts - i*0.0001;
565: X[j][i].Ta = X[j][i].Ts - 5;
566: X[j][i].u = 0;
567: X[j][i].v = 0;
568: X[j][i].p = 1.25;
569: if ((j == 5 || j == 6) && (i == 4 || i == 5)) X[j][i].p += 0.00001;
570: if ((j == 5 || j == 6) && (i == 12 || i == 13)) X[j][i].p += 0.00001;
571: }
572: }
573: } else {
574: for (j=ys; j<ys+ym; j++) {
575: for (i=xs; i<xs+xm; i++) {
576: X[j][i].Ts = user->Ts;
577: X[j][i].Ta = X[j][i].Ts - 5;
578: X[j][i].u = 0;
579: X[j][i].v = 0;
580: X[j][i].p = 1.25;
581: }
582: }
583: }
585: /* Restore vectors */
586: DMDAVecRestoreArray(da,Xglobal,&X);
587: return(0);
588: }
592: /*
593: RhsFunc - Evaluates nonlinear function F(u).
595: Input Parameters:
596: . ts - the TS context
597: . t - current time
598: . Xglobal - input vector
599: . F - output vector
600: . ptr - optional user-defined context, as set by SNESSetFunction()
602: Output Parameter:
603: . F - rhs function vector
604: */
605: PetscErrorCode RhsFunc(TS ts,PetscReal t,Vec Xglobal,Vec F,void *ctx)
606: {
607: AppCtx *user = (AppCtx*)ctx; /* user-defined application context */
608: DM da = user->da;
610: PetscInt i,j,Mx,My,xs,ys,xm,ym;
611: PetscReal dhx,dhy;
612: Vec localT;
613: Field **X,**Frhs; /* structures that contain variables of interest and left hand side of governing equations respectively */
614: PetscScalar csoil = user->csoil; /* heat constant for layer */
615: PetscScalar dzlay = user->dzlay; /* thickness of top soil layer */
616: PetscScalar emma = user->emma; /* emission parameter */
617: PetscScalar wind = user->wind; /* wind speed */
618: PetscScalar dewtemp = user->dewtemp; /* dew point temperature (moisture in air) */
619: PetscScalar pressure1 = user->pressure1; /* sea level pressure */
620: PetscScalar airtemp = user->airtemp; /* temperature of air near boundary layer inversion */
621: PetscScalar fract = user->fract; /* fraction of the sky covered by clouds */
622: PetscScalar Tc = user->Tc; /* temperature at base of lowest cloud layer */
623: PetscScalar lat = user->lat; /* latitude */
624: PetscScalar Cp = 1005.7; /* specific heat of air at constant pressure */
625: PetscScalar Rd = 287.058; /* gas constant for dry air */
626: PetscScalar diffconst = 1000; /* diffusion coefficient */
627: PetscScalar f = 2*0.0000727*sin(lat); /* coriolis force */
628: PetscScalar deep_grnd_temp = user->deep_grnd_temp; /* temp in lowest ground layer */
629: PetscScalar Ts,u,v,p,P;
630: PetscScalar u_abs,u_plus,u_minus,v_abs,v_plus,v_minus;
632: PetscScalar sfctemp1,fsfc1,Ra;
633: PetscScalar sheat; /* sensible heat flux */
634: PetscScalar latentheat; /* latent heat flux */
635: PetscScalar groundflux; /* flux from conduction of deep ground layer in contact with top soil */
636: PetscInt xend,yend;
639: DMGetLocalVector(da,&localT);
640: DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
641: PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
643: dhx = (PetscReal)(Mx-1)/(5000*(Mx-1)); /* dhx = 1/dx; assume 2D space domain: [0.0, 1.e5] x [0.0, 1.e5] */
644: dhy = (PetscReal)(My-1)/(5000*(Mx-1)); /* dhy = 1/dy; */
647: /*
648: Scatter ghost points to local vector,using the 2-step process
649: DAGlobalToLocalBegin(),DAGlobalToLocalEnd().
650: By placing code between these two statements, computations can be
651: done while messages are in transition.
652: */
653: DMGlobalToLocalBegin(da,Xglobal,INSERT_VALUES,localT);
654: DMGlobalToLocalEnd(da,Xglobal,INSERT_VALUES,localT);
656: /* Get pointers to vector data */
657: DMDAVecGetArray(da,localT,&X);
658: DMDAVecGetArray(da,F,&Frhs);
660: /* Get local grid boundaries */
661: DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);
663: /* Compute function over the locally owned part of the grid */
664: /* the interior points */
665: xend=xs+xm; yend=ys+ym;
666: for (j=ys; j<yend; j++) {
667: for (i=xs; i<xend; i++) {
668: 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; */
670: sfctemp1 = (double)Ts;
671: sfctemp1 = (double)X[j][i].Ts;
672: calcfluxs(sfctemp1,airtemp,emma,fract,Tc,&fsfc1); /* calculates surface net radiative flux */
673: sensibleflux(sfctemp1,airtemp,wind,&sheat); /* calculate sensible heat flux */
674: latentflux(sfctemp1,dewtemp,wind,pressure1,&latentheat); /* calculates latent heat flux */
675: calc_gflux(sfctemp1,deep_grnd_temp,&groundflux); /* calculates flux from earth below surface soil layer by conduction */
676: calcfluxa(sfctemp1,airtemp,emma,&Ra); /* Calculates the change in downward radiative flux */
677: fsfc1 = fsfc1 + latentheat + sheat + groundflux; /* adds radiative, sensible heat, latent heat, and ground heat flux yielding net flux */
679: /* convective coefficients for upwinding */
680: u_abs = PetscAbsScalar(u);
681: u_plus = .5*(u + u_abs); /* u if u>0; 0 if u<0 */
682: u_minus = .5*(u - u_abs); /* u if u <0; 0 if u>0 */
684: v_abs = PetscAbsScalar(v);
685: v_plus = .5*(v + v_abs); /* v if v>0; 0 if v<0 */
686: v_minus = .5*(v - v_abs); /* v if v <0; 0 if v>0 */
688: /* Solve governing equations */
689: P = p*Rd*Ts;
691: /* du/dt -> time change of east-west component of the wind */
692: 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) */
693: - v_plus*(u - X[j-1][i].u)*dhy - v_minus*(X[j+1][i].u - u)*dhy /* - v(du/dy) */
694: -(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)] */
695: /* -(1/p)*(X[j][i+1].P - X[j][i-1].P)*dhx */
696: + f*v;
698: /* dv/dt -> time change of north-south component of the wind */
699: 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) */
700: - v_plus*(v - X[j-1][i].v)*dhy - v_minus*(X[j+1][i].v - v)*dhy /* - v(dv/dy) */
701: -(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)] */
702: /* -(1/p)*(X[j+1][i].P - X[j-1][i].P)*dhy */
703: -f*u;
705: /* dT/dt -> time change of temperature */
706: Frhs[j][i].Ts = (fsfc1/(csoil*dzlay)) /* Fnet/(Cp*dz) diabatic change in T */
707: -u_plus*(Ts - X[j][i-1].Ts)*dhx - u_minus*(X[j][i+1].Ts - Ts)*dhx /* - u*(dTs/dx) advection x */
708: -v_plus*(Ts - X[j-1][i].Ts)*dhy - v_minus*(X[j+1][i].Ts - Ts)*dhy /* - v*(dTs/dy) advection y */
709: + diffconst*((X[j][i+1].Ts - 2*Ts + X[j][i-1].Ts)*dhx*dhx /* + D(Ts_xx + Ts_yy) diffusion */
710: + (X[j+1][i].Ts - 2*Ts + X[j-1][i].Ts)*dhy*dhy);
712: /* dp/dt -> time change of */
713: 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) */
714: -v_plus*(p - X[j-1][i].p)*dhy - v_minus*(X[j+1][i].p - p)*dhy; /* - v*(dp/dy) */
716: Frhs[j][i].Ta = Ra/Cp; /* dTa/dt time change of air temperature */
717: }
718: }
720: /* Restore vectors */
721: DMDAVecRestoreArray(da,localT,&X);
722: DMDAVecRestoreArray(da,F,&Frhs);
723: DMRestoreLocalVector(da,&localT);
724: return(0);
725: }
729: PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal time,Vec T,void *ctx)
730: {
732: PetscScalar *array;
733: MonitorCtx *user = (MonitorCtx*)ctx;
734: PetscViewer viewer = user->drawviewer;
735: PetscMPIInt rank;
736: PetscReal norm;
739: MPI_Comm_rank(PetscObjectComm((PetscObject)ts),&rank);
740: VecNorm(T,NORM_INFINITY,&norm);
742: VecGetArray(T,&array);
743: if (!rank) printf("step %4d, time %8.1f, %6.4f, %6.4f, %6.4f, %6.4f, %6.4f, %6.4f\n",step,time,(((array[0]-273)*9)/5 + 32),(((array[1]-273)*9)/5 + 32),array[2],array[3],array[4],array[5]);
744: VecRestoreArray(T,&array);
746: if (user->drawcontours) {
747: VecView(T,viewer);
748: }
749: return(0);
750: }