Actual source code: ex5.c
petsc-3.4.3 2013-10-15
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;
175: PetscInitialize(&argc,&argv,(char*)0,help);
176: MPI_Comm_size(PETSC_COMM_WORLD,&size);
177: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
179: /* Inputs */
180: readinput(&put);
182: sfctemp = put.Ts;
183: dewtemp = put.Td;
184: cloudTemp = put.Tc;
185: airtemp = put.Ta;
186: pwat = put.pwt;
188: if (!rank) PetscPrintf(PETSC_COMM_SELF,"Initial Temperature = %g\n",sfctemp); /* input surface temperature */
190: deep_grnd_temp = sfctemp - 10; /* set underlying ground layer temperature */
191: emma = emission(pwat); /* accounts for radiative effects of water vapor */
193: /* Converts from Fahrenheit to Celsuis */
194: sfctemp = fahr_to_cel(sfctemp);
195: airtemp = fahr_to_cel(airtemp);
196: dewtemp = fahr_to_cel(dewtemp);
197: cloudTemp = fahr_to_cel(cloudTemp);
198: deep_grnd_temp = fahr_to_cel(deep_grnd_temp);
200: /* Converts from Celsius to Kelvin */
201: sfctemp += 273;
202: airtemp += 273;
203: dewtemp += 273;
204: cloudTemp += 273;
205: deep_grnd_temp += 273;
207: /* Calculates initial relative humidity */
208: x = calcmixingr(dewtemp,pressure1);
209: mixratio = calcmixingr(sfctemp,pressure1);
210: rh = (x/mixratio)*100;
212: if (!rank) printf("Initial RH = %.1f percent\n\n",rh); /* prints initial relative humidity */
214: time = 3600*put.time; /* sets amount of timesteps to run model */
216: /* Configure PETSc TS solver */
217: /*------------------------------------------*/
219: /* Create grid */
220: DMDACreate2d(PETSC_COMM_WORLD,DMDA_BOUNDARY_PERIODIC,DMDA_BOUNDARY_PERIODIC,DMDA_STENCIL_STAR,-20,-20,
221: PETSC_DECIDE,PETSC_DECIDE,dof,1,NULL,NULL,&da);
222: DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0);
224: /* Define output window for each variable of interest */
225: DMDASetFieldName(da,0,"Ts");
226: DMDASetFieldName(da,1,"Ta");
227: DMDASetFieldName(da,2,"u");
228: DMDASetFieldName(da,3,"v");
229: DMDASetFieldName(da,4,"p");
231: /* set values for appctx */
232: user.da = da;
233: user.Ts = sfctemp;
234: user.fract = put.fr; /* fraction of sky covered by clouds */
235: user.dewtemp = dewtemp; /* dew point temperature (mositure in air) */
236: user.csoil = 2000000; /* heat constant for layer */
237: user.dzlay = 0.08; /* thickness of top soil layer */
238: user.emma = emma; /* emission parameter */
239: user.wind = put.wnd; /* wind spped */
240: user.pressure1 = pressure1; /* sea level pressure */
241: user.airtemp = airtemp; /* temperature of air near boundar layer inversion */
242: user.Tc = cloudTemp; /* temperature at base of lowest cloud layer */
243: user.init = put.init; /* user chosen initiation scenario */
244: user.lat = 70*0.0174532; /* converts latitude degrees to latitude in radians */
245: user.deep_grnd_temp = deep_grnd_temp; /* temp in lowest ground layer */
247: /* set values for MonitorCtx */
248: usermonitor.drawcontours = PETSC_FALSE;
249: PetscOptionsHasName(NULL,"-drawcontours",&usermonitor.drawcontours);
250: if (usermonitor.drawcontours) {
251: PetscReal bounds[] = {1000.0,-1000., -1000.,-1000., 1000.,-1000., 1000.,-1000., 1000,-1000, 100700,100800};
252: PetscViewerDrawOpen(PETSC_COMM_WORLD,0,0,0,0,300,300,&usermonitor.drawviewer);
253: PetscViewerDrawSetBounds(usermonitor.drawviewer,dof,bounds);
254: }
255: usermonitor.interval = 1;
256: PetscOptionsGetInt(NULL,"-monitor_interval",&usermonitor.interval,NULL);
258: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
259: Extract global vectors from DA;
260: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
261: DMCreateGlobalVector(da,&T);
262: VecDuplicate(T,&rhs); /* r: vector to put the computed right hand side */
264: TSCreate(PETSC_COMM_WORLD,&ts);
265: TSSetProblemType(ts,TS_NONLINEAR);
266: TSSetType(ts,TSBEULER);
267: TSSetRHSFunction(ts,rhs,RhsFunc,&user);
269: /* Set Jacobian evaluation routine - use coloring to compute finite difference Jacobian efficiently */
270: PetscBool use_coloring = PETSC_TRUE;
271: MatFDColoring matfdcoloring = 0;
272: DMCreateMatrix(da,MATAIJ,&J);
273: TSGetSNES(ts,&snes);
274: if (use_coloring) {
275: ISColoring iscoloring;
276: DMCreateColoring(da,IS_COLORING_GLOBAL,MATAIJ,&iscoloring);
277: MatFDColoringCreate(J,iscoloring,&matfdcoloring);
278: MatFDColoringSetFromOptions(matfdcoloring);
279: ISColoringDestroy(&iscoloring);
280: MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode (*)(void))SNESTSFormFunction,ts);
281: SNESSetJacobian(snes,J,J,SNESComputeJacobianDefaultColor,matfdcoloring);
282: } else {
283: SNESSetJacobian(snes,J,J,SNESComputeJacobianDefault,NULL);
284: }
286: /* Define what to print for ts_monitor option */
287: PetscBool monitor_off = PETSC_FALSE;
288: PetscOptionsHasName(NULL,"-monitor_off",&monitor_off);
289: if (!monitor_off) {
290: TSMonitorSet(ts,Monitor,&usermonitor,NULL);
291: }
292: FormInitialSolution(da,T,&user);
293: dt = TIMESTEP; /* initial time step */
294: ftime = TIMESTEP*time;
295: if (!rank) printf("time %d, ftime %g hour, TIMESTEP %g\n",time,ftime/3600,dt);
297: TSSetInitialTimeStep(ts,0.0,dt);
298: TSSetDuration(ts,time,ftime);
299: TSSetSolution(ts,T);
300: TSSetDM(ts,da);
302: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
303: Set runtime options
304: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
305: TSSetFromOptions(ts);
307: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
308: Solve nonlinear system
309: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
310: TSSolve(ts,T);
311: TSGetSolveTime(ts,&ftime);
312: TSGetTimeStepNumber(ts,&steps);
313: if (!rank) PetscPrintf(PETSC_COMM_WORLD,"Solution T after %g hours %d steps\n",ftime/3600,steps);
316: if (matfdcoloring) {MatFDColoringDestroy(&matfdcoloring);}
317: if (usermonitor.drawcontours) {
318: PetscViewerDestroy(&usermonitor.drawviewer);
319: }
320: MatDestroy(&J);
321: VecDestroy(&T);
322: VecDestroy(&rhs);
323: TSDestroy(&ts);
324: DMDestroy(&da);
326: PetscFinalize();
327: return 0;
328: }
329: /*****************************end main program********************************/
330: /*****************************************************************************/
331: /*****************************************************************************/
332: /*****************************************************************************/
335: PetscErrorCode calcfluxs(PetscScalar sfctemp, PetscScalar airtemp, PetscScalar emma, PetscScalar fract, PetscScalar cloudTemp, PetscScalar *flux)
336: {
338: *flux = SIG*((EMMSFC*emma*pow(airtemp,4)) + (EMMSFC*fract*(1 - emma)*pow(cloudTemp,4)) - (EMMSFC*pow(sfctemp,4))); /* calculates flux using Stefan-Boltzmann relation */
339: return(0);
340: }
344: PetscErrorCode calcfluxa(PetscScalar sfctemp, PetscScalar airtemp, PetscScalar emma, PetscScalar *flux) /* this function is not currently called upon */
345: {
346: PetscScalar emm = 0.001;
349: *flux = SIG*(-emm*(pow(airtemp,4))); /* calculates flux usinge Stefan-Boltzmann relation */
350: return(0);
351: }
354: PetscErrorCode sensibleflux(PetscScalar sfctemp, PetscScalar airtemp, PetscScalar wind, PetscScalar *sheat)
355: {
356: PetscScalar density = 1; /* air density */
357: PetscScalar Cp = 1005; /* heat capicity for dry air */
358: PetscScalar wndmix; /* temperature change from wind mixing: wind*Ch */
361: wndmix = 0.0025 + 0.0042*wind; /* regression equation valid for neutral and stable BL */
362: *sheat = density*Cp*wndmix*(airtemp - sfctemp); /* calculates sensible heat flux */
363: return(0);
364: }
368: PetscErrorCode latentflux(PetscScalar sfctemp, PetscScalar dewtemp, PetscScalar wind, PetscScalar pressure1, PetscScalar *latentheat)
369: {
370: PetscScalar density = 1; /* density of dry air */
371: PetscScalar q; /* actual specific humitity */
372: PetscScalar qs; /* saturation specific humidity */
373: PetscScalar wndmix; /* temperature change from wind mixing: wind*Ch */
374: PetscScalar beta = .4; /* moisture availability */
375: PetscScalar mr; /* mixing ratio */
376: PetscScalar lhcnst; /* latent heat of vaporization constant = 2501000 J/kg at 0c */
377: /* latent heat of saturation const = 2834000 J/kg */
378: /* latent heat of fusion const = 333700 J/kg */
381: wind = mph2mpers(wind); /* converts wind from mph to meters per second */
382: wndmix = 0.0025 + 0.0042*wind; /* regression equation valid for neutral BL */
383: lhcnst = Lconst(sfctemp); /* calculates latent heat of evaporation */
384: mr = calcmixingr(sfctemp,pressure1); /* calculates saturation mixing ratio */
385: qs = calc_q(mr); /* calculates saturation specific humidty */
386: mr = calcmixingr(dewtemp,pressure1); /* calculates mixing ratio */
387: q = calc_q(mr); /* calculates specific humidty */
389: *latentheat = density*wndmix*beta*lhcnst*(q - qs); /* calculates latent heat flux */
390: return(0);
391: }
395: PetscErrorCode potential_temperature(PetscScalar temp, PetscScalar pressure1, PetscScalar pressure2, PetscScalar sfctemp, PetscScalar *pottemp)
396: {
397: PetscScalar kdry; /* poisson constant for dry atmosphere */
398: PetscScalar kmoist; /* poisson constant for moist atmosphere */
399: PetscScalar pavg; /* average atmospheric pressure */
400: PetscScalar mixratio; /* mixing ratio */
403: mixratio = calcmixingr(sfctemp,pressure1);
405: /* initialize poisson constant */
406: kdry = 0.2854;
407: kmoist = 0.2854*(1 - 0.24*mixratio);
409: pavg = ((0.7*pressure1)+pressure2)/2; /* calculates simple average press */
410: *pottemp = temp*(pow((pressure1/pavg),kdry)); /* calculates potential temperature */
411: return(0);
412: }
413: extern PetscScalar calcmixingr(PetscScalar dtemp, PetscScalar pressure1)
414: {
415: PetscScalar e; /* vapor pressure */
416: PetscScalar mixratio; /* mixing ratio */
418: dtemp = dtemp - 273; /* converts from Kelvin to Celsuis */
419: e = 6.11*(pow(10,((7.5*dtemp)/(237.7+dtemp)))); /* converts from dew point temp to vapor pressure */
420: e = e*100; /* converts from hPa to Pa */
421: mixratio = (0.622*e)/(pressure1 - e); /* computes mixing ratio */
422: mixratio = mixratio*1; /* convert to g/Kg */
424: return mixratio;
425: }
426: extern PetscScalar calc_q(PetscScalar rv)
427: {
428: PetscScalar specific_humidity; /* define specific humidity variable */
429: specific_humidity = rv/(1 + rv); /* calculates specific humidity */
430: return specific_humidity;
431: }
435: PetscErrorCode calc_gflux(PetscScalar sfctemp, PetscScalar deep_grnd_temp, PetscScalar* Gflux)
436: {
437: PetscScalar k; /* thermal conductivity parameter */
438: PetscScalar n = 0.38; /* value of soil porosity */
439: PetscScalar dz = 1; /* depth of layer between soil surface and deep soil layer */
440: PetscScalar unit_soil_weight = 2700; /* unit soil weight in kg/m^3 */
443: k = ((0.135*(1-n)*unit_soil_weight) + 64.7)/(unit_soil_weight - (0.947*(1-n)*unit_soil_weight)); /* dry soil conductivity */
444: *Gflux = (k*(deep_grnd_temp - sfctemp)/dz); /* calculates flux from deep ground layer */
445: return(0);
446: }
449: extern PetscScalar emission(PetscScalar pwat)
450: {
451: PetscScalar emma;
453: emma = 0.725 + 0.17*log10(pwat);
455: return emma;
456: }
457: extern PetscScalar cloud(PetscScalar fract)
458: {
459: PetscScalar emma = 0;
461: /* modifies radiative balance depending on cloud cover */
462: if (fract >= 0.9) emma = 1;
463: else if (0.9 > fract && fract >= 0.8) emma = 0.9;
464: else if (0.8 > fract && fract >= 0.7) emma = 0.85;
465: else if (0.7 > fract && fract >= 0.6) emma = 0.75;
466: else if (0.6 > fract && fract >= 0.5) emma = 0.65;
467: else if (0.4 > fract && fract >= 0.3) emma = emma*1.086956;
468: return emma;
469: }
470: extern PetscScalar Lconst(PetscScalar sfctemp)
471: {
472: PetscScalar Lheat;
473: sfctemp -=273; /* converts from kelvin to celsius */
474: Lheat = 4186.8*(597.31 - 0.5625*sfctemp); /* calculates latent heat constant */
475: return Lheat;
476: }
477: extern PetscScalar mph2mpers(PetscScalar wind)
478: {
479: wind = ((wind*1.6*1000)/3600); /* converts wind from mph to meters per second */
480: return wind;
481: }
482: extern PetscScalar fahr_to_cel(PetscScalar temp)
483: {
484: temp = (5*(temp-32))/9; /* converts from farhrenheit to celsuis */
485: return temp;
486: }
487: extern PetscScalar cel_to_fahr(PetscScalar temp)
488: {
489: temp = ((temp*9)/5) + 32; /* converts from celsuis to farhrenheit */
490: return temp;
491: }
492: void readinput(struct in *put)
493: {
494: int i;
495: char x;
496: FILE *ifp;
498: ifp = fopen("ex5_control.txt", "r");
500: for (i=0; i<110; i++) fscanf(ifp, "%c", &x);
501: fscanf(ifp, "%lf", &put->Ts);
503: for (i=0; i<43; i++) fscanf(ifp, "%c", &x);
504: fscanf(ifp, "%lf", &put->Td);
506: for (i=0; i<43; i++) fscanf(ifp, "%c", &x);
507: fscanf(ifp, "%lf", &put->Ta);
509: for (i=0; i<43; i++) fscanf(ifp, "%c", &x);
510: fscanf(ifp, "%lf", &put->Tc);
512: for (i=0; i<43; i++) fscanf(ifp, "%c", &x);
513: fscanf(ifp, "%lf", &put->fr);
515: for (i=0; i<43; i++) fscanf(ifp, "%c", &x);
516: fscanf(ifp, "%lf", &put->wnd);
518: for (i=0; i<43; i++) fscanf(ifp, "%c", &x);
519: fscanf(ifp, "%lf", &put->pwt);
521: for (i=0; i<43; i++) fscanf(ifp, "%c", &x);
522: fscanf(ifp, "%lf", &put->wndDir);
524: for (i=0; i<43; i++) fscanf(ifp, "%c", &x);
525: fscanf(ifp, "%lf", &put->time);
527: for (i=0; i<63; i++) fscanf(ifp, "%c", &x);
528: fscanf(ifp, "%lf", &put->init);
529: }
531: /* ------------------------------------------------------------------- */
534: PetscErrorCode FormInitialSolution(DM da,Vec Xglobal,void *ctx)
535: {
537: AppCtx *user = (AppCtx*)ctx; /* user-defined application context */
538: PetscInt i,j,xs,ys,xm,ym,Mx,My;
539: Field **X;
540: PetscScalar deltT;
541: PetscReal hx,hy;
542: FILE *ifp;
543: FILE *ofp;
546: ofp = fopen("swing", "w");
547: ifp = fopen("grid.in", "r");
548: deltT = 0.8;
550: DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
551: PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
553: hx = 1/(PetscReal)(Mx-1);
554: hy = 1/(PetscReal)(My-1);
556: /* Get pointers to vector data */
557: DMDAVecGetArray(da,Xglobal,&X);
559: /* Get local grid boundaries */
560: DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);
562: /* Compute function over the locally owned part of the grid */
564: if (user->init == 1) {
565: for (j=ys; j<ys+ym; j++) {
566: for (i=xs; i<xs+xm; i++) {
567: X[j][i].Ts = user->Ts - i*0.0001;
568: X[j][i].Ta = X[j][i].Ts - 5;
569: X[j][i].u = 0;
570: X[j][i].v = 0;
571: X[j][i].p = 1.25;
572: if ((j == 5 || j == 6) && (i == 4 || i == 5)) X[j][i].p += 0.00001;
573: if ((j == 5 || j == 6) && (i == 12 || i == 13)) X[j][i].p += 0.00001;
574: }
575: }
576: } else {
577: for (j=ys; j<ys+ym; j++) {
578: for (i=xs; i<xs+xm; i++) {
579: X[j][i].Ts = user->Ts;
580: X[j][i].Ta = X[j][i].Ts - 5;
581: X[j][i].u = 0;
582: X[j][i].v = 0;
583: X[j][i].p = 1.25;
584: }
585: }
586: }
588: /* Restore vectors */
589: DMDAVecRestoreArray(da,Xglobal,&X);
590: return(0);
591: }
595: /*
596: RhsFunc - Evaluates nonlinear function F(u).
598: Input Parameters:
599: . ts - the TS context
600: . t - current time
601: . Xglobal - input vector
602: . F - output vector
603: . ptr - optional user-defined context, as set by SNESSetFunction()
605: Output Parameter:
606: . F - rhs function vector
607: */
608: PetscErrorCode RhsFunc(TS ts,PetscReal t,Vec Xglobal,Vec F,void *ctx)
609: {
610: AppCtx *user = (AppCtx*)ctx; /* user-defined application context */
611: DM da = user->da;
613: PetscInt i,j,Mx,My,xs,ys,xm,ym;
614: PetscReal dhx,dhy;
615: Vec localT;
616: Field **X,**Frhs; /* structures that contain variables of interest and left hand side of governing equations respectively */
617: PetscScalar csoil = user->csoil; /* heat constant for layer */
618: PetscScalar dzlay = user->dzlay; /* thickness of top soil layer */
619: PetscScalar emma = user->emma; /* emission parameter */
620: PetscScalar wind = user->wind; /* wind speed */
621: PetscScalar dewtemp = user->dewtemp; /* dew point temperature (moisture in air) */
622: PetscScalar pressure1 = user->pressure1; /* sea level pressure */
623: PetscScalar airtemp = user->airtemp; /* temperature of air near boundary layer inversion */
624: PetscScalar fract = user->fract; /* fraction of the sky covered by clouds */
625: PetscScalar Tc = user->Tc; /* temperature at base of lowest cloud layer */
626: PetscScalar lat = user->lat; /* latitude */
627: PetscScalar Cp = 1005.7; /* specific heat of air at constant pressure */
628: PetscScalar Rd = 287.058; /* gas constant for dry air */
629: PetscScalar diffconst = 1000; /* diffusion coefficient */
630: PetscScalar f = 2*0.0000727*sin(lat); /* coriolis force */
631: PetscScalar deep_grnd_temp = user->deep_grnd_temp; /* temp in lowest ground layer */
632: PetscScalar Ts,u,v,p,P;
633: PetscScalar u_abs,u_plus,u_minus,v_abs,v_plus,v_minus;
635: PetscScalar sfctemp1,fsfc1,Ra;
636: PetscScalar sheat; /* sensible heat flux */
637: PetscScalar latentheat; /* latent heat flux */
638: PetscScalar groundflux; /* flux from conduction of deep ground layer in contact with top soil */
639: PetscInt xend,yend;
642: DMGetLocalVector(da,&localT);
643: DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
644: PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
646: dhx = (PetscReal)(Mx-1)/(5000*(Mx-1)); /* dhx = 1/dx; assume 2D space domain: [0.0, 1.e5] x [0.0, 1.e5] */
647: dhy = (PetscReal)(My-1)/(5000*(Mx-1)); /* dhy = 1/dy; */
650: /*
651: Scatter ghost points to local vector,using the 2-step process
652: DAGlobalToLocalBegin(),DAGlobalToLocalEnd().
653: By placing code between these two statements, computations can be
654: done while messages are in transition.
655: */
656: DMGlobalToLocalBegin(da,Xglobal,INSERT_VALUES,localT);
657: DMGlobalToLocalEnd(da,Xglobal,INSERT_VALUES,localT);
659: /* Get pointers to vector data */
660: DMDAVecGetArray(da,localT,&X);
661: DMDAVecGetArray(da,F,&Frhs);
663: /* Get local grid boundaries */
664: DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);
666: /* Compute function over the locally owned part of the grid */
667: /* the interior points */
668: xend=xs+xm; yend=ys+ym;
669: for (j=ys; j<yend; j++) {
670: for (i=xs; i<xend; i++) {
671: 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; */
673: sfctemp1 = (double)Ts;
674: sfctemp1 = (double)X[j][i].Ts;
675: calcfluxs(sfctemp1,airtemp,emma,fract,Tc,&fsfc1); /* calculates surface net radiative flux */
676: sensibleflux(sfctemp1,airtemp,wind,&sheat); /* calculate sensible heat flux */
677: latentflux(sfctemp1,dewtemp,wind,pressure1,&latentheat); /* calculates latent heat flux */
678: calc_gflux(sfctemp1,deep_grnd_temp,&groundflux); /* calculates flux from earth below surface soil layer by conduction */
679: calcfluxa(sfctemp1,airtemp,emma,&Ra); /* Calculates the change in downward radiative flux */
680: fsfc1 = fsfc1 + latentheat + sheat + groundflux; /* adds radiative, sensible heat, latent heat, and ground heat flux yielding net flux */
682: /* convective coefficients for upwinding */
683: u_abs = PetscAbsScalar(u);
684: u_plus = .5*(u + u_abs); /* u if u>0; 0 if u<0 */
685: u_minus = .5*(u - u_abs); /* u if u <0; 0 if u>0 */
687: v_abs = PetscAbsScalar(v);
688: v_plus = .5*(v + v_abs); /* v if v>0; 0 if v<0 */
689: v_minus = .5*(v - v_abs); /* v if v <0; 0 if v>0 */
691: /* Solve governing equations */
692: P = p*Rd*Ts;
694: /* du/dt -> time change of east-west component of the wind */
695: 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) */
696: - v_plus*(u - X[j-1][i].u)*dhy - v_minus*(X[j+1][i].u - u)*dhy /* - v(du/dy) */
697: -(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)] */
698: /* -(1/p)*(X[j][i+1].P - X[j][i-1].P)*dhx */
699: + f*v;
701: /* dv/dt -> time change of north-south component of the wind */
702: 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) */
703: - v_plus*(v - X[j-1][i].v)*dhy - v_minus*(X[j+1][i].v - v)*dhy /* - v(dv/dy) */
704: -(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)] */
705: /* -(1/p)*(X[j+1][i].P - X[j-1][i].P)*dhy */
706: -f*u;
708: /* dT/dt -> time change of temperature */
709: Frhs[j][i].Ts = (fsfc1/(csoil*dzlay)) /* Fnet/(Cp*dz) diabatic change in T */
710: -u_plus*(Ts - X[j][i-1].Ts)*dhx - u_minus*(X[j][i+1].Ts - Ts)*dhx /* - u*(dTs/dx) advection x */
711: -v_plus*(Ts - X[j-1][i].Ts)*dhy - v_minus*(X[j+1][i].Ts - Ts)*dhy /* - v*(dTs/dy) advection y */
712: + diffconst*((X[j][i+1].Ts - 2*Ts + X[j][i-1].Ts)*dhx*dhx /* + D(Ts_xx + Ts_yy) diffusion */
713: + (X[j+1][i].Ts - 2*Ts + X[j-1][i].Ts)*dhy*dhy);
715: /* dp/dt -> time change of */
716: 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) */
717: -v_plus*(p - X[j-1][i].p)*dhy - v_minus*(X[j+1][i].p - p)*dhy; /* - v*(dp/dy) */
719: Frhs[j][i].Ta = Ra/Cp; /* dTa/dt time change of air temperature */
720: }
721: }
723: /* Restore vectors */
724: DMDAVecRestoreArray(da,localT,&X);
725: DMDAVecRestoreArray(da,F,&Frhs);
726: DMRestoreLocalVector(da,&localT);
727: return(0);
728: }
732: PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal time,Vec T,void *ctx)
733: {
735: PetscScalar *array;
736: MonitorCtx *user = (MonitorCtx*)ctx;
737: PetscViewer viewer = user->drawviewer;
738: PetscMPIInt rank;
739: PetscReal norm;
742: MPI_Comm_rank(PetscObjectComm((PetscObject)ts),&rank);
743: VecNorm(T,NORM_INFINITY,&norm);
745: if (step%user->interval == 0) {
746: VecGetArray(T,&array);
747: 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]);
748: VecRestoreArray(T,&array);
749: }
751: if (user->drawcontours) {
752: VecView(T,viewer);
753: }
754: return(0);
755: }