xref: /petsc/src/ts/tests/ex5.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown static char help[] = "Nonlinear, time-dependent. Developed from radiative_surface_balance.c \n";
2*c4762a1bSJed Brown /*
3*c4762a1bSJed Brown   Contributed by Steve Froehlich, Illinois Institute of Technology
4*c4762a1bSJed Brown 
5*c4762a1bSJed Brown    Usage:
6*c4762a1bSJed Brown     mpiexec -n <np> ./ex5 [options]
7*c4762a1bSJed Brown     ./ex5 -help  [view petsc options]
8*c4762a1bSJed Brown     ./ex5 -ts_type sundials -ts_view
9*c4762a1bSJed Brown     ./ex5 -da_grid_x 20 -da_grid_y 20 -log_view
10*c4762a1bSJed Brown     ./ex5 -da_grid_x 20 -da_grid_y 20 -ts_type rosw -ts_atol 1.e-6 -ts_rtol 1.e-6
11*c4762a1bSJed Brown     ./ex5 -drawcontours -draw_pause 0.1 -draw_fields 0,1,2,3,4
12*c4762a1bSJed Brown */
13*c4762a1bSJed Brown 
14*c4762a1bSJed Brown /*
15*c4762a1bSJed Brown    -----------------------------------------------------------------------
16*c4762a1bSJed Brown 
17*c4762a1bSJed Brown    Governing equations:
18*c4762a1bSJed Brown 
19*c4762a1bSJed Brown         R      = s*(Ea*Ta^4 - Es*Ts^4)
20*c4762a1bSJed Brown         SH     = p*Cp*Ch*wind*(Ta - Ts)
21*c4762a1bSJed Brown         LH     = p*L*Ch*wind*B(q(Ta) - q(Ts))
22*c4762a1bSJed Brown         G      = k*(Tgnd - Ts)/dz
23*c4762a1bSJed Brown 
24*c4762a1bSJed Brown         Fnet   = R + SH + LH + G
25*c4762a1bSJed Brown 
26*c4762a1bSJed Brown         du/dt  = -u*(du/dx) - v*(du/dy) - 2*omeg*sin(lat)*v - (1/p)*(dP/dx)
27*c4762a1bSJed Brown         dv/dt  = -u*(dv/dx) - v*(dv/dy) + 2*omeg*sin(lat)*u - (1/p)*(dP/dy)
28*c4762a1bSJed Brown         dTs/dt = Fnet/(Cp*dz) - Div([u*Ts, v*Ts]) + D*Lap(Ts)
29*c4762a1bSJed Brown                = Fnet/(Cs*dz) - u*(dTs/dx) - v*(dTs/dy) + D*(Ts_xx + Ts_yy)
30*c4762a1bSJed Brown         dp/dt  = -Div([u*p,v*p])
31*c4762a1bSJed Brown                = - u*dp/dx - v*dp/dy
32*c4762a1bSJed Brown         dTa/dt = Fnet/Cp
33*c4762a1bSJed Brown 
34*c4762a1bSJed Brown    Equation of State:
35*c4762a1bSJed Brown 
36*c4762a1bSJed Brown         P = p*R*Ts
37*c4762a1bSJed Brown 
38*c4762a1bSJed Brown    -----------------------------------------------------------------------
39*c4762a1bSJed Brown 
40*c4762a1bSJed Brown    Program considers the evolution of a two dimensional atmosphere from
41*c4762a1bSJed Brown    sunset to sunrise. There are two components:
42*c4762a1bSJed Brown                 1. Surface energy balance model to compute diabatic dT (Fnet)
43*c4762a1bSJed Brown                 2. Dynamical model using simplified primitive equations
44*c4762a1bSJed Brown 
45*c4762a1bSJed Brown    Program is to be initiated at sunset and run to sunrise.
46*c4762a1bSJed Brown 
47*c4762a1bSJed Brown    Inputs are:
48*c4762a1bSJed Brown                 Surface temperature
49*c4762a1bSJed Brown                 Dew point temperature
50*c4762a1bSJed Brown                 Air temperature
51*c4762a1bSJed Brown                 Temperature at cloud base (if clouds are present)
52*c4762a1bSJed Brown                 Fraction of sky covered by clouds
53*c4762a1bSJed Brown                 Wind speed
54*c4762a1bSJed Brown                 Precipitable water in centimeters
55*c4762a1bSJed Brown                 Wind direction
56*c4762a1bSJed Brown 
57*c4762a1bSJed Brown    Inputs are are read in from the text file ex5_control.txt. To change an
58*c4762a1bSJed Brown    input value use ex5_control.txt.
59*c4762a1bSJed Brown 
60*c4762a1bSJed Brown    Solvers:
61*c4762a1bSJed Brown             Backward Euler = default solver
62*c4762a1bSJed Brown             Sundials = fastest and most accurate, requires Sundials libraries
63*c4762a1bSJed Brown 
64*c4762a1bSJed Brown    This model is under development and should be used only as an example
65*c4762a1bSJed Brown    and not as a predictive weather model.
66*c4762a1bSJed Brown */
67*c4762a1bSJed Brown 
68*c4762a1bSJed Brown #include <petscts.h>
69*c4762a1bSJed Brown #include <petscdm.h>
70*c4762a1bSJed Brown #include <petscdmda.h>
71*c4762a1bSJed Brown 
72*c4762a1bSJed Brown /* stefan-boltzmann constant */
73*c4762a1bSJed Brown #define SIG 0.000000056703
74*c4762a1bSJed Brown /* absorption-emission constant for surface */
75*c4762a1bSJed Brown #define EMMSFC   1
76*c4762a1bSJed Brown /* amount of time (seconds) that passes before new flux is calculated */
77*c4762a1bSJed Brown #define TIMESTEP 1
78*c4762a1bSJed Brown 
79*c4762a1bSJed Brown /* variables of interest to be solved at each grid point */
80*c4762a1bSJed Brown typedef struct {
81*c4762a1bSJed Brown   PetscScalar Ts,Ta; /* surface and air temperature */
82*c4762a1bSJed Brown   PetscScalar u,v;   /* wind speed */
83*c4762a1bSJed Brown   PetscScalar p;     /* density */
84*c4762a1bSJed Brown } Field;
85*c4762a1bSJed Brown 
86*c4762a1bSJed Brown /* User defined variables. Used in solving for variables of interest */
87*c4762a1bSJed Brown typedef struct {
88*c4762a1bSJed Brown   DM          da;        /* grid */
89*c4762a1bSJed Brown   PetscScalar csoil;     /* heat constant for layer */
90*c4762a1bSJed Brown   PetscScalar dzlay;     /* thickness of top soil layer */
91*c4762a1bSJed Brown   PetscScalar emma;      /* emission parameter */
92*c4762a1bSJed Brown   PetscScalar wind;      /* wind speed */
93*c4762a1bSJed Brown   PetscScalar dewtemp;   /* dew point temperature (moisture in air) */
94*c4762a1bSJed Brown   PetscScalar pressure1; /* sea level pressure */
95*c4762a1bSJed Brown   PetscScalar airtemp;   /* temperature of air near boundary layer inversion */
96*c4762a1bSJed Brown   PetscScalar Ts;        /* temperature at the surface */
97*c4762a1bSJed Brown   PetscScalar fract;     /* fraction of sky covered by clouds */
98*c4762a1bSJed Brown   PetscScalar Tc;        /* temperature at base of lowest cloud layer */
99*c4762a1bSJed Brown   PetscScalar lat;       /* Latitude in degrees */
100*c4762a1bSJed Brown   PetscScalar init;      /* initialization scenario */
101*c4762a1bSJed Brown   PetscScalar deep_grnd_temp; /* temperature of ground under top soil surface layer */
102*c4762a1bSJed Brown } AppCtx;
103*c4762a1bSJed Brown 
104*c4762a1bSJed Brown /* Struct for visualization */
105*c4762a1bSJed Brown typedef struct {
106*c4762a1bSJed Brown   PetscBool   drawcontours;   /* flag - 1 indicates drawing contours */
107*c4762a1bSJed Brown   PetscViewer drawviewer;
108*c4762a1bSJed Brown   PetscInt    interval;
109*c4762a1bSJed Brown } MonitorCtx;
110*c4762a1bSJed Brown 
111*c4762a1bSJed Brown 
112*c4762a1bSJed Brown /* Inputs read in from text file */
113*c4762a1bSJed Brown struct in {
114*c4762a1bSJed Brown   PetscScalar Ts;     /* surface temperature  */
115*c4762a1bSJed Brown   PetscScalar Td;     /* dewpoint temperature */
116*c4762a1bSJed Brown   PetscScalar Tc;     /* temperature of cloud base */
117*c4762a1bSJed Brown   PetscScalar fr;     /* fraction of sky covered by clouds */
118*c4762a1bSJed Brown   PetscScalar wnd;    /* wind speed */
119*c4762a1bSJed Brown   PetscScalar Ta;     /* air temperature */
120*c4762a1bSJed Brown   PetscScalar pwt;    /* precipitable water */
121*c4762a1bSJed Brown   PetscScalar wndDir; /* wind direction */
122*c4762a1bSJed Brown   PetscScalar lat;    /* latitude */
123*c4762a1bSJed Brown   PetscReal   time;   /* time in hours */
124*c4762a1bSJed Brown   PetscScalar init;
125*c4762a1bSJed Brown };
126*c4762a1bSJed Brown 
127*c4762a1bSJed Brown /* functions */
128*c4762a1bSJed Brown extern PetscScalar emission(PetscScalar);                           /* sets emission/absorption constant depending on water vapor content */
129*c4762a1bSJed Brown extern PetscScalar calc_q(PetscScalar);                             /* calculates specific humidity */
130*c4762a1bSJed Brown extern PetscScalar mph2mpers(PetscScalar);                          /* converts miles per hour to meters per second */
131*c4762a1bSJed Brown 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*c4762a1bSJed Brown extern PetscScalar fahr_to_cel(PetscScalar);                        /* converts Fahrenheit to Celsius */
133*c4762a1bSJed Brown extern PetscScalar cel_to_fahr(PetscScalar);                        /* converts Celsius to Fahrenheit */
134*c4762a1bSJed Brown extern PetscScalar calcmixingr(PetscScalar, PetscScalar);           /* calculates mixing ratio */
135*c4762a1bSJed Brown extern PetscScalar cloud(PetscScalar);                              /* cloud radiative parameterization */
136*c4762a1bSJed Brown extern PetscErrorCode FormInitialSolution(DM,Vec,void*);            /* Specifies initial conditions for the system of equations (PETSc defined function) */
137*c4762a1bSJed Brown extern PetscErrorCode RhsFunc(TS,PetscReal,Vec,Vec,void*);          /* Specifies the user defined functions                     (PETSc defined function) */
138*c4762a1bSJed Brown extern PetscErrorCode Monitor(TS,PetscInt,PetscReal,Vec,void*);     /* Specifies output and visualization tools                 (PETSc defined function) */
139*c4762a1bSJed Brown extern void readinput(struct in *put);                              /* reads input from text file */
140*c4762a1bSJed Brown extern PetscErrorCode calcfluxs(PetscScalar, PetscScalar, PetscScalar, PetscScalar, PetscScalar, PetscScalar*); /* calculates upward IR from surface */
141*c4762a1bSJed Brown extern PetscErrorCode calcfluxa(PetscScalar, PetscScalar, PetscScalar, PetscScalar*);                           /* calculates downward IR from atmosphere */
142*c4762a1bSJed Brown extern PetscErrorCode sensibleflux(PetscScalar, PetscScalar, PetscScalar, PetscScalar*);                        /* calculates sensible heat flux */
143*c4762a1bSJed Brown extern PetscErrorCode potential_temperature(PetscScalar, PetscScalar, PetscScalar, PetscScalar, PetscScalar*);  /* calculates potential temperature */
144*c4762a1bSJed Brown extern PetscErrorCode latentflux(PetscScalar, PetscScalar, PetscScalar, PetscScalar, PetscScalar*);             /* calculates latent heat flux */
145*c4762a1bSJed Brown extern PetscErrorCode calc_gflux(PetscScalar, PetscScalar, PetscScalar*);                                       /* calculates flux between top soil layer and underlying earth */
146*c4762a1bSJed Brown 
147*c4762a1bSJed Brown int main(int argc,char **argv)
148*c4762a1bSJed Brown {
149*c4762a1bSJed Brown   PetscErrorCode ierr;
150*c4762a1bSJed Brown   int            time;           /* amount of loops */
151*c4762a1bSJed Brown   struct in      put;
152*c4762a1bSJed Brown   PetscScalar    rh;             /* relative humidity */
153*c4762a1bSJed Brown   PetscScalar    x;              /* memory varialbe for relative humidity calculation */
154*c4762a1bSJed Brown   PetscScalar    deep_grnd_temp; /* temperature of ground under top soil surface layer */
155*c4762a1bSJed Brown   PetscScalar    emma;           /* absorption-emission constant for air */
156*c4762a1bSJed Brown   PetscScalar    pressure1 = 101300; /* surface pressure */
157*c4762a1bSJed Brown   PetscScalar    mixratio;       /* mixing ratio */
158*c4762a1bSJed Brown   PetscScalar    airtemp;        /* temperature of air near boundary layer inversion */
159*c4762a1bSJed Brown   PetscScalar    dewtemp;        /* dew point temperature */
160*c4762a1bSJed Brown   PetscScalar    sfctemp;        /* temperature at surface */
161*c4762a1bSJed Brown   PetscScalar    pwat;           /* total column precipitable water */
162*c4762a1bSJed Brown   PetscScalar    cloudTemp;      /* temperature at base of cloud */
163*c4762a1bSJed Brown   AppCtx         user;           /*  user-defined work context */
164*c4762a1bSJed Brown   MonitorCtx     usermonitor;    /* user-defined monitor context */
165*c4762a1bSJed Brown   PetscMPIInt    rank,size;
166*c4762a1bSJed Brown   TS             ts;
167*c4762a1bSJed Brown   SNES           snes;
168*c4762a1bSJed Brown   DM             da;
169*c4762a1bSJed Brown   Vec            T,rhs;          /* solution vector */
170*c4762a1bSJed Brown   Mat            J;              /* Jacobian matrix */
171*c4762a1bSJed Brown   PetscReal      ftime,dt;
172*c4762a1bSJed Brown   PetscInt       steps,dof = 5;
173*c4762a1bSJed Brown   PetscBool      use_coloring  = PETSC_TRUE;
174*c4762a1bSJed Brown   MatFDColoring  matfdcoloring = 0;
175*c4762a1bSJed Brown   PetscBool      monitor_off = PETSC_FALSE;
176*c4762a1bSJed Brown 
177*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
178*c4762a1bSJed Brown   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
179*c4762a1bSJed Brown   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
180*c4762a1bSJed Brown 
181*c4762a1bSJed Brown   /* Inputs */
182*c4762a1bSJed Brown   readinput(&put);
183*c4762a1bSJed Brown 
184*c4762a1bSJed Brown   sfctemp   = put.Ts;
185*c4762a1bSJed Brown   dewtemp   = put.Td;
186*c4762a1bSJed Brown   cloudTemp = put.Tc;
187*c4762a1bSJed Brown   airtemp   = put.Ta;
188*c4762a1bSJed Brown   pwat      = put.pwt;
189*c4762a1bSJed Brown 
190*c4762a1bSJed Brown   if (!rank) PetscPrintf(PETSC_COMM_SELF,"Initial Temperature = %g\n",(double)sfctemp); /* input surface temperature */
191*c4762a1bSJed Brown 
192*c4762a1bSJed Brown   deep_grnd_temp = sfctemp - 10;   /* set underlying ground layer temperature */
193*c4762a1bSJed Brown   emma           = emission(pwat); /* accounts for radiative effects of water vapor */
194*c4762a1bSJed Brown 
195*c4762a1bSJed Brown   /* Converts from Fahrenheit to Celsuis */
196*c4762a1bSJed Brown   sfctemp        = fahr_to_cel(sfctemp);
197*c4762a1bSJed Brown   airtemp        = fahr_to_cel(airtemp);
198*c4762a1bSJed Brown   dewtemp        = fahr_to_cel(dewtemp);
199*c4762a1bSJed Brown   cloudTemp      = fahr_to_cel(cloudTemp);
200*c4762a1bSJed Brown   deep_grnd_temp = fahr_to_cel(deep_grnd_temp);
201*c4762a1bSJed Brown 
202*c4762a1bSJed Brown   /* Converts from Celsius to Kelvin */
203*c4762a1bSJed Brown   sfctemp        += 273;
204*c4762a1bSJed Brown   airtemp        += 273;
205*c4762a1bSJed Brown   dewtemp        += 273;
206*c4762a1bSJed Brown   cloudTemp      += 273;
207*c4762a1bSJed Brown   deep_grnd_temp += 273;
208*c4762a1bSJed Brown 
209*c4762a1bSJed Brown   /* Calculates initial relative humidity */
210*c4762a1bSJed Brown   x        = calcmixingr(dewtemp,pressure1);
211*c4762a1bSJed Brown   mixratio = calcmixingr(sfctemp,pressure1);
212*c4762a1bSJed Brown   rh       = (x/mixratio)*100;
213*c4762a1bSJed Brown 
214*c4762a1bSJed Brown   if (!rank) printf("Initial RH = %.1f percent\n\n",(double)rh);   /* prints initial relative humidity */
215*c4762a1bSJed Brown 
216*c4762a1bSJed Brown   time = 3600*put.time;                         /* sets amount of timesteps to run model */
217*c4762a1bSJed Brown 
218*c4762a1bSJed Brown   /* Configure PETSc TS solver */
219*c4762a1bSJed Brown   /*------------------------------------------*/
220*c4762a1bSJed Brown 
221*c4762a1bSJed Brown   /* Create grid */
222*c4762a1bSJed Brown   ierr = DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_PERIODIC,DM_BOUNDARY_PERIODIC,DMDA_STENCIL_STAR,20,20,PETSC_DECIDE,PETSC_DECIDE,dof,1,NULL,NULL,&da);CHKERRQ(ierr);
223*c4762a1bSJed Brown   ierr = DMSetFromOptions(da);CHKERRQ(ierr);
224*c4762a1bSJed Brown   ierr = DMSetUp(da);CHKERRQ(ierr);
225*c4762a1bSJed Brown   ierr = DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0);CHKERRQ(ierr);
226*c4762a1bSJed Brown 
227*c4762a1bSJed Brown   /* Define output window for each variable of interest */
228*c4762a1bSJed Brown   ierr = DMDASetFieldName(da,0,"Ts");CHKERRQ(ierr);
229*c4762a1bSJed Brown   ierr = DMDASetFieldName(da,1,"Ta");CHKERRQ(ierr);
230*c4762a1bSJed Brown   ierr = DMDASetFieldName(da,2,"u");CHKERRQ(ierr);
231*c4762a1bSJed Brown   ierr = DMDASetFieldName(da,3,"v");CHKERRQ(ierr);
232*c4762a1bSJed Brown   ierr = DMDASetFieldName(da,4,"p");CHKERRQ(ierr);
233*c4762a1bSJed Brown 
234*c4762a1bSJed Brown   /* set values for appctx */
235*c4762a1bSJed Brown   user.da             = da;
236*c4762a1bSJed Brown   user.Ts             = sfctemp;
237*c4762a1bSJed Brown   user.fract          = put.fr;          /* fraction of sky covered by clouds */
238*c4762a1bSJed Brown   user.dewtemp        = dewtemp;         /* dew point temperature (mositure in air) */
239*c4762a1bSJed Brown   user.csoil          = 2000000;         /* heat constant for layer */
240*c4762a1bSJed Brown   user.dzlay          = 0.08;            /* thickness of top soil layer */
241*c4762a1bSJed Brown   user.emma           = emma;            /* emission parameter */
242*c4762a1bSJed Brown   user.wind           = put.wnd;         /* wind spped */
243*c4762a1bSJed Brown   user.pressure1      = pressure1;       /* sea level pressure */
244*c4762a1bSJed Brown   user.airtemp        = airtemp;         /* temperature of air near boundar layer inversion */
245*c4762a1bSJed Brown   user.Tc             = cloudTemp;       /* temperature at base of lowest cloud layer */
246*c4762a1bSJed Brown   user.init           = put.init;        /* user chosen initiation scenario */
247*c4762a1bSJed Brown   user.lat            = 70*0.0174532;    /* converts latitude degrees to latitude in radians */
248*c4762a1bSJed Brown   user.deep_grnd_temp = deep_grnd_temp;  /* temp in lowest ground layer */
249*c4762a1bSJed Brown 
250*c4762a1bSJed Brown   /* set values for MonitorCtx */
251*c4762a1bSJed Brown   usermonitor.drawcontours = PETSC_FALSE;
252*c4762a1bSJed Brown   ierr = PetscOptionsHasName(NULL,NULL,"-drawcontours",&usermonitor.drawcontours);CHKERRQ(ierr);
253*c4762a1bSJed Brown   if (usermonitor.drawcontours) {
254*c4762a1bSJed Brown     PetscReal bounds[] = {1000.0,-1000.,  -1000.,-1000.,  1000.,-1000.,  1000.,-1000.,  1000,-1000, 100700,100800};
255*c4762a1bSJed Brown     ierr = PetscViewerDrawOpen(PETSC_COMM_WORLD,0,0,0,0,300,300,&usermonitor.drawviewer);CHKERRQ(ierr);
256*c4762a1bSJed Brown     ierr = PetscViewerDrawSetBounds(usermonitor.drawviewer,dof,bounds);CHKERRQ(ierr);
257*c4762a1bSJed Brown   }
258*c4762a1bSJed Brown   usermonitor.interval = 1;
259*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-monitor_interval",&usermonitor.interval,NULL);CHKERRQ(ierr);
260*c4762a1bSJed Brown 
261*c4762a1bSJed Brown   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
262*c4762a1bSJed Brown      Extract global vectors from DA;
263*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
264*c4762a1bSJed Brown   ierr = DMCreateGlobalVector(da,&T);CHKERRQ(ierr);
265*c4762a1bSJed Brown   ierr = VecDuplicate(T,&rhs);CHKERRQ(ierr); /* r: vector to put the computed right hand side */
266*c4762a1bSJed Brown 
267*c4762a1bSJed Brown   ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr);
268*c4762a1bSJed Brown   ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr);
269*c4762a1bSJed Brown   ierr = TSSetType(ts,TSBEULER);CHKERRQ(ierr);
270*c4762a1bSJed Brown   ierr = TSSetRHSFunction(ts,rhs,RhsFunc,&user);CHKERRQ(ierr);
271*c4762a1bSJed Brown 
272*c4762a1bSJed Brown   /* Set Jacobian evaluation routine - use coloring to compute finite difference Jacobian efficiently */
273*c4762a1bSJed Brown   ierr = DMSetMatType(da,MATAIJ);CHKERRQ(ierr);
274*c4762a1bSJed Brown   ierr = DMCreateMatrix(da,&J);CHKERRQ(ierr);
275*c4762a1bSJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
276*c4762a1bSJed Brown   if (use_coloring) {
277*c4762a1bSJed Brown     ISColoring iscoloring;
278*c4762a1bSJed Brown     ierr = DMCreateColoring(da,IS_COLORING_GLOBAL,&iscoloring);CHKERRQ(ierr);
279*c4762a1bSJed Brown     ierr = MatFDColoringCreate(J,iscoloring,&matfdcoloring);CHKERRQ(ierr);
280*c4762a1bSJed Brown     ierr = MatFDColoringSetFromOptions(matfdcoloring);CHKERRQ(ierr);
281*c4762a1bSJed Brown     ierr = MatFDColoringSetUp(J,iscoloring,matfdcoloring);CHKERRQ(ierr);
282*c4762a1bSJed Brown     ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr);
283*c4762a1bSJed Brown     ierr = MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode (*)(void))SNESTSFormFunction,ts);CHKERRQ(ierr);
284*c4762a1bSJed Brown     ierr = SNESSetJacobian(snes,J,J,SNESComputeJacobianDefaultColor,matfdcoloring);CHKERRQ(ierr);
285*c4762a1bSJed Brown   } else {
286*c4762a1bSJed Brown     ierr = SNESSetJacobian(snes,J,J,SNESComputeJacobianDefault,NULL);CHKERRQ(ierr);
287*c4762a1bSJed Brown   }
288*c4762a1bSJed Brown 
289*c4762a1bSJed Brown   /* Define what to print for ts_monitor option */
290*c4762a1bSJed Brown   ierr = PetscOptionsHasName(NULL,NULL,"-monitor_off",&monitor_off);CHKERRQ(ierr);
291*c4762a1bSJed Brown   if (!monitor_off) {
292*c4762a1bSJed Brown     ierr = TSMonitorSet(ts,Monitor,&usermonitor,NULL);CHKERRQ(ierr);
293*c4762a1bSJed Brown   }
294*c4762a1bSJed Brown   ierr  = FormInitialSolution(da,T,&user);CHKERRQ(ierr);
295*c4762a1bSJed Brown   dt    = TIMESTEP; /* initial time step */
296*c4762a1bSJed Brown   ftime = TIMESTEP*time;
297*c4762a1bSJed Brown   if (!rank) printf("time %d, ftime %g hour, TIMESTEP %g\n",time,(double)(ftime/3600),(double)dt);
298*c4762a1bSJed Brown 
299*c4762a1bSJed Brown   ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr);
300*c4762a1bSJed Brown   ierr = TSSetMaxSteps(ts,time);CHKERRQ(ierr);
301*c4762a1bSJed Brown   ierr = TSSetMaxTime(ts,ftime);CHKERRQ(ierr);
302*c4762a1bSJed Brown   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr);
303*c4762a1bSJed Brown   ierr = TSSetSolution(ts,T);CHKERRQ(ierr);
304*c4762a1bSJed Brown   ierr = TSSetDM(ts,da);CHKERRQ(ierr);
305*c4762a1bSJed Brown 
306*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
307*c4762a1bSJed Brown      Set runtime options
308*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
309*c4762a1bSJed Brown   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
310*c4762a1bSJed Brown 
311*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
312*c4762a1bSJed Brown      Solve nonlinear system
313*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
314*c4762a1bSJed Brown   ierr = TSSolve(ts,T);CHKERRQ(ierr);
315*c4762a1bSJed Brown   ierr = TSGetSolveTime(ts,&ftime);CHKERRQ(ierr);
316*c4762a1bSJed Brown   ierr = TSGetStepNumber(ts,&steps);CHKERRQ(ierr);
317*c4762a1bSJed Brown   if (!rank) PetscPrintf(PETSC_COMM_WORLD,"Solution T after %g hours %d steps\n",(double)(ftime/3600),steps);
318*c4762a1bSJed Brown 
319*c4762a1bSJed Brown 
320*c4762a1bSJed Brown   if (matfdcoloring) {ierr = MatFDColoringDestroy(&matfdcoloring);CHKERRQ(ierr);}
321*c4762a1bSJed Brown   if (usermonitor.drawcontours) {
322*c4762a1bSJed Brown     ierr = PetscViewerDestroy(&usermonitor.drawviewer);CHKERRQ(ierr);
323*c4762a1bSJed Brown   }
324*c4762a1bSJed Brown   ierr = MatDestroy(&J);CHKERRQ(ierr);
325*c4762a1bSJed Brown   ierr = VecDestroy(&T);CHKERRQ(ierr);
326*c4762a1bSJed Brown   ierr = VecDestroy(&rhs);CHKERRQ(ierr);
327*c4762a1bSJed Brown   ierr = TSDestroy(&ts);CHKERRQ(ierr);
328*c4762a1bSJed Brown   ierr = DMDestroy(&da);CHKERRQ(ierr);
329*c4762a1bSJed Brown 
330*c4762a1bSJed Brown   ierr = PetscFinalize();
331*c4762a1bSJed Brown   return ierr;
332*c4762a1bSJed Brown }
333*c4762a1bSJed Brown /*****************************end main program********************************/
334*c4762a1bSJed Brown /*****************************************************************************/
335*c4762a1bSJed Brown /*****************************************************************************/
336*c4762a1bSJed Brown /*****************************************************************************/
337*c4762a1bSJed Brown PetscErrorCode calcfluxs(PetscScalar sfctemp, PetscScalar airtemp, PetscScalar emma, PetscScalar fract, PetscScalar cloudTemp, PetscScalar *flux)
338*c4762a1bSJed Brown {
339*c4762a1bSJed Brown   PetscFunctionBeginUser;
340*c4762a1bSJed Brown   *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*c4762a1bSJed Brown   PetscFunctionReturn(0);
342*c4762a1bSJed Brown }
343*c4762a1bSJed Brown 
344*c4762a1bSJed Brown PetscErrorCode calcfluxa(PetscScalar sfctemp, PetscScalar airtemp, PetscScalar emma, PetscScalar *flux)   /* this function is not currently called upon */
345*c4762a1bSJed Brown {
346*c4762a1bSJed Brown   PetscScalar emm = 0.001;
347*c4762a1bSJed Brown 
348*c4762a1bSJed Brown   PetscFunctionBeginUser;
349*c4762a1bSJed Brown   *flux = SIG*(-emm*(PetscPowScalarInt(airtemp,4)));      /* calculates flux usinge Stefan-Boltzmann relation */
350*c4762a1bSJed Brown   PetscFunctionReturn(0);
351*c4762a1bSJed Brown }
352*c4762a1bSJed Brown PetscErrorCode sensibleflux(PetscScalar sfctemp, PetscScalar airtemp, PetscScalar wind, PetscScalar *sheat)
353*c4762a1bSJed Brown {
354*c4762a1bSJed Brown   PetscScalar density = 1;    /* air density */
355*c4762a1bSJed Brown   PetscScalar Cp      = 1005; /* heat capicity for dry air */
356*c4762a1bSJed Brown   PetscScalar wndmix;         /* temperature change from wind mixing: wind*Ch */
357*c4762a1bSJed Brown 
358*c4762a1bSJed Brown   PetscFunctionBeginUser;
359*c4762a1bSJed Brown   wndmix = 0.0025 + 0.0042*wind;                               /* regression equation valid for neutral and stable BL */
360*c4762a1bSJed Brown   *sheat = density*Cp*wndmix*(airtemp - sfctemp);              /* calculates sensible heat flux */
361*c4762a1bSJed Brown   PetscFunctionReturn(0);
362*c4762a1bSJed Brown }
363*c4762a1bSJed Brown 
364*c4762a1bSJed Brown PetscErrorCode latentflux(PetscScalar sfctemp, PetscScalar dewtemp, PetscScalar wind, PetscScalar pressure1, PetscScalar *latentheat)
365*c4762a1bSJed Brown {
366*c4762a1bSJed Brown   PetscScalar density = 1;   /* density of dry air */
367*c4762a1bSJed Brown   PetscScalar q;             /* actual specific humitity */
368*c4762a1bSJed Brown   PetscScalar qs;            /* saturation specific humidity */
369*c4762a1bSJed Brown   PetscScalar wndmix;        /* temperature change from wind mixing: wind*Ch */
370*c4762a1bSJed Brown   PetscScalar beta = .4;     /* moisture availability */
371*c4762a1bSJed Brown   PetscScalar mr;            /* mixing ratio */
372*c4762a1bSJed Brown   PetscScalar lhcnst;        /* latent heat of vaporization constant = 2501000 J/kg at 0c */
373*c4762a1bSJed Brown                              /* latent heat of saturation const = 2834000 J/kg */
374*c4762a1bSJed Brown                              /* latent heat of fusion const = 333700 J/kg */
375*c4762a1bSJed Brown 
376*c4762a1bSJed Brown   PetscFunctionBeginUser;
377*c4762a1bSJed Brown   wind   = mph2mpers(wind);                /* converts wind from mph to meters per second */
378*c4762a1bSJed Brown   wndmix = 0.0025 + 0.0042*wind;           /* regression equation valid for neutral BL */
379*c4762a1bSJed Brown   lhcnst = Lconst(sfctemp);                /* calculates latent heat of evaporation */
380*c4762a1bSJed Brown   mr     = calcmixingr(sfctemp,pressure1); /* calculates saturation mixing ratio */
381*c4762a1bSJed Brown   qs     = calc_q(mr);                     /* calculates saturation specific humidty */
382*c4762a1bSJed Brown   mr     = calcmixingr(dewtemp,pressure1); /* calculates mixing ratio */
383*c4762a1bSJed Brown   q      = calc_q(mr);                     /* calculates specific humidty */
384*c4762a1bSJed Brown 
385*c4762a1bSJed Brown   *latentheat = density*wndmix*beta*lhcnst*(q - qs); /* calculates latent heat flux */
386*c4762a1bSJed Brown   PetscFunctionReturn(0);
387*c4762a1bSJed Brown }
388*c4762a1bSJed Brown 
389*c4762a1bSJed Brown PetscErrorCode potential_temperature(PetscScalar temp, PetscScalar pressure1, PetscScalar pressure2, PetscScalar sfctemp, PetscScalar *pottemp)
390*c4762a1bSJed Brown {
391*c4762a1bSJed Brown   PetscScalar kdry;     /* poisson constant for dry atmosphere */
392*c4762a1bSJed Brown   PetscScalar pavg;     /* average atmospheric pressure */
393*c4762a1bSJed Brown   /* PetscScalar mixratio; mixing ratio */
394*c4762a1bSJed Brown   /* PetscScalar kmoist;   poisson constant for moist atmosphere */
395*c4762a1bSJed Brown 
396*c4762a1bSJed Brown   PetscFunctionBeginUser;
397*c4762a1bSJed Brown   /* mixratio = calcmixingr(sfctemp,pressure1); */
398*c4762a1bSJed Brown 
399*c4762a1bSJed Brown /* initialize poisson constant */
400*c4762a1bSJed Brown   kdry   = 0.2854;
401*c4762a1bSJed Brown   /* kmoist = 0.2854*(1 - 0.24*mixratio); */
402*c4762a1bSJed Brown 
403*c4762a1bSJed Brown   pavg     = ((0.7*pressure1)+pressure2)/2;     /* calculates simple average press */
404*c4762a1bSJed Brown   *pottemp = temp*(PetscPowScalar((pressure1/pavg),kdry)); /* calculates potential temperature */
405*c4762a1bSJed Brown   PetscFunctionReturn(0);
406*c4762a1bSJed Brown }
407*c4762a1bSJed Brown extern PetscScalar calcmixingr(PetscScalar dtemp, PetscScalar pressure1)
408*c4762a1bSJed Brown {
409*c4762a1bSJed Brown   PetscScalar e;        /* vapor pressure */
410*c4762a1bSJed Brown   PetscScalar mixratio; /* mixing ratio */
411*c4762a1bSJed Brown 
412*c4762a1bSJed Brown   dtemp    = dtemp - 273;                                /* converts from Kelvin to Celsuis */
413*c4762a1bSJed Brown   e        = 6.11*(PetscPowScalar(10,((7.5*dtemp)/(237.7+dtemp)))); /* converts from dew point temp to vapor pressure */
414*c4762a1bSJed Brown   e        = e*100;                                      /* converts from hPa to Pa */
415*c4762a1bSJed Brown   mixratio = (0.622*e)/(pressure1 - e);                  /* computes mixing ratio */
416*c4762a1bSJed Brown   mixratio = mixratio*1;                                 /* convert to g/Kg */
417*c4762a1bSJed Brown 
418*c4762a1bSJed Brown   return mixratio;
419*c4762a1bSJed Brown }
420*c4762a1bSJed Brown extern PetscScalar calc_q(PetscScalar rv)
421*c4762a1bSJed Brown {
422*c4762a1bSJed Brown   PetscScalar specific_humidity;        /* define specific humidity variable */
423*c4762a1bSJed Brown   specific_humidity = rv/(1 + rv);      /* calculates specific humidity */
424*c4762a1bSJed Brown   return specific_humidity;
425*c4762a1bSJed Brown }
426*c4762a1bSJed Brown 
427*c4762a1bSJed Brown PetscErrorCode calc_gflux(PetscScalar sfctemp, PetscScalar deep_grnd_temp, PetscScalar* Gflux)
428*c4762a1bSJed Brown {
429*c4762a1bSJed Brown   PetscScalar k;                       /* thermal conductivity parameter */
430*c4762a1bSJed Brown   PetscScalar n                = 0.38; /* value of soil porosity */
431*c4762a1bSJed Brown   PetscScalar dz               = 1;    /* depth of layer between soil surface and deep soil layer */
432*c4762a1bSJed Brown   PetscScalar unit_soil_weight = 2700; /* unit soil weight in kg/m^3 */
433*c4762a1bSJed Brown 
434*c4762a1bSJed Brown   PetscFunctionBeginUser;
435*c4762a1bSJed Brown   k      = ((0.135*(1-n)*unit_soil_weight) + 64.7)/(unit_soil_weight - (0.947*(1-n)*unit_soil_weight)); /* dry soil conductivity */
436*c4762a1bSJed Brown   *Gflux = (k*(deep_grnd_temp - sfctemp)/dz);   /* calculates flux from deep ground layer */
437*c4762a1bSJed Brown   PetscFunctionReturn(0);
438*c4762a1bSJed Brown }
439*c4762a1bSJed Brown extern PetscScalar emission(PetscScalar pwat)
440*c4762a1bSJed Brown {
441*c4762a1bSJed Brown   PetscScalar emma;
442*c4762a1bSJed Brown 
443*c4762a1bSJed Brown   emma = 0.725 + 0.17*PetscLog10Real(PetscRealPart(pwat));
444*c4762a1bSJed Brown 
445*c4762a1bSJed Brown   return emma;
446*c4762a1bSJed Brown }
447*c4762a1bSJed Brown extern PetscScalar cloud(PetscScalar fract)
448*c4762a1bSJed Brown {
449*c4762a1bSJed Brown   PetscScalar emma = 0;
450*c4762a1bSJed Brown 
451*c4762a1bSJed Brown   /* modifies radiative balance depending on cloud cover */
452*c4762a1bSJed Brown   if (fract >= 0.9)                     emma = 1;
453*c4762a1bSJed Brown   else if (0.9 > fract && fract >= 0.8) emma = 0.9;
454*c4762a1bSJed Brown   else if (0.8 > fract && fract >= 0.7) emma = 0.85;
455*c4762a1bSJed Brown   else if (0.7 > fract && fract >= 0.6) emma = 0.75;
456*c4762a1bSJed Brown   else if (0.6 > fract && fract >= 0.5) emma = 0.65;
457*c4762a1bSJed Brown   else if (0.4 > fract && fract >= 0.3) emma = emma*1.086956;
458*c4762a1bSJed Brown   return emma;
459*c4762a1bSJed Brown }
460*c4762a1bSJed Brown extern PetscScalar Lconst(PetscScalar sfctemp)
461*c4762a1bSJed Brown {
462*c4762a1bSJed Brown   PetscScalar Lheat;
463*c4762a1bSJed Brown   sfctemp -=273;                               /* converts from kelvin to celsius */
464*c4762a1bSJed Brown   Lheat    = 4186.8*(597.31 - 0.5625*sfctemp); /* calculates latent heat constant */
465*c4762a1bSJed Brown   return Lheat;
466*c4762a1bSJed Brown }
467*c4762a1bSJed Brown extern PetscScalar mph2mpers(PetscScalar wind)
468*c4762a1bSJed Brown {
469*c4762a1bSJed Brown   wind = ((wind*1.6*1000)/3600);                 /* converts wind from mph to meters per second */
470*c4762a1bSJed Brown   return wind;
471*c4762a1bSJed Brown }
472*c4762a1bSJed Brown extern PetscScalar fahr_to_cel(PetscScalar temp)
473*c4762a1bSJed Brown {
474*c4762a1bSJed Brown   temp = (5*(temp-32))/9; /* converts from farhrenheit to celsuis */
475*c4762a1bSJed Brown   return temp;
476*c4762a1bSJed Brown }
477*c4762a1bSJed Brown extern PetscScalar cel_to_fahr(PetscScalar temp)
478*c4762a1bSJed Brown {
479*c4762a1bSJed Brown   temp = ((temp*9)/5) + 32; /* converts from celsuis to farhrenheit */
480*c4762a1bSJed Brown   return temp;
481*c4762a1bSJed Brown }
482*c4762a1bSJed Brown void readinput(struct in *put)
483*c4762a1bSJed Brown {
484*c4762a1bSJed Brown   int    i;
485*c4762a1bSJed Brown   char   x;
486*c4762a1bSJed Brown   FILE   *ifp;
487*c4762a1bSJed Brown   double tmp;
488*c4762a1bSJed Brown 
489*c4762a1bSJed Brown   ifp = fopen("ex5_control.txt", "r");
490*c4762a1bSJed Brown 
491*c4762a1bSJed Brown   for (i=0; i<110; i++) { if (fscanf(ifp, "%c", &x) != 1) abort();}
492*c4762a1bSJed Brown   if (fscanf(ifp, "%lf", &tmp) != 1) abort();
493*c4762a1bSJed Brown   put->Ts = tmp;
494*c4762a1bSJed Brown 
495*c4762a1bSJed Brown   for (i=0; i<43; i++) { if (fscanf(ifp, "%c", &x) != 1) abort();}
496*c4762a1bSJed Brown   if (fscanf(ifp, "%lf", &tmp) != 1) abort();
497*c4762a1bSJed Brown   put->Td = tmp;
498*c4762a1bSJed Brown 
499*c4762a1bSJed Brown   for (i=0; i<43; i++) { if (fscanf(ifp, "%c", &x) != 1) abort();}
500*c4762a1bSJed Brown   if (fscanf(ifp, "%lf", &tmp) != 1) abort();
501*c4762a1bSJed Brown   put->Ta = tmp;
502*c4762a1bSJed Brown 
503*c4762a1bSJed Brown   for (i=0; i<43; i++) { if (fscanf(ifp, "%c", &x) != 1) abort();}
504*c4762a1bSJed Brown   if (fscanf(ifp, "%lf", &tmp)!= 1) abort();
505*c4762a1bSJed Brown   put->Tc = tmp;
506*c4762a1bSJed Brown 
507*c4762a1bSJed Brown   for (i=0; i<43; i++) { if (fscanf(ifp, "%c", &x) != 1) abort();}
508*c4762a1bSJed Brown   if (fscanf(ifp, "%lf", &tmp) != 1) abort();
509*c4762a1bSJed Brown   put->fr = tmp;
510*c4762a1bSJed Brown 
511*c4762a1bSJed Brown   for (i=0; i<43; i++) {if (fscanf(ifp, "%c", &x) != 1) abort();}
512*c4762a1bSJed Brown   if (fscanf(ifp, "%lf", &tmp) != 1) abort();
513*c4762a1bSJed Brown   put->wnd = tmp;
514*c4762a1bSJed Brown 
515*c4762a1bSJed Brown   for (i=0; i<43; i++) {if (fscanf(ifp, "%c", &x) != 1) abort();}
516*c4762a1bSJed Brown   if (fscanf(ifp, "%lf", &tmp) != 1) abort();
517*c4762a1bSJed Brown   put->pwt = tmp;
518*c4762a1bSJed Brown 
519*c4762a1bSJed Brown   for (i=0; i<43; i++) {if (fscanf(ifp, "%c", &x) != 1) abort();}
520*c4762a1bSJed Brown   if (fscanf(ifp, "%lf", &tmp) != 1) abort();
521*c4762a1bSJed Brown   put->wndDir = tmp;
522*c4762a1bSJed Brown 
523*c4762a1bSJed Brown   for (i=0; i<43; i++) {if (fscanf(ifp, "%c", &x) != 1) abort();}
524*c4762a1bSJed Brown   if (fscanf(ifp, "%lf", &tmp) != 1) abort();
525*c4762a1bSJed Brown   put->time = tmp;
526*c4762a1bSJed Brown 
527*c4762a1bSJed Brown   for (i=0; i<63; i++) {if (fscanf(ifp, "%c", &x) != 1) abort();}
528*c4762a1bSJed Brown   if (fscanf(ifp, "%lf", &tmp) != 1) abort();
529*c4762a1bSJed Brown   put->init = tmp;
530*c4762a1bSJed Brown 
531*c4762a1bSJed Brown }
532*c4762a1bSJed Brown 
533*c4762a1bSJed Brown /* ------------------------------------------------------------------- */
534*c4762a1bSJed Brown PetscErrorCode FormInitialSolution(DM da,Vec Xglobal,void *ctx)
535*c4762a1bSJed Brown {
536*c4762a1bSJed Brown   PetscErrorCode ierr;
537*c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ctx;       /* user-defined application context */
538*c4762a1bSJed Brown   PetscInt       i,j,xs,ys,xm,ym,Mx,My;
539*c4762a1bSJed Brown   Field          **X;
540*c4762a1bSJed Brown 
541*c4762a1bSJed Brown   PetscFunctionBeginUser;
542*c4762a1bSJed Brown   ierr = DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
543*c4762a1bSJed Brown                      PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);CHKERRQ(ierr);
544*c4762a1bSJed Brown 
545*c4762a1bSJed Brown   /* Get pointers to vector data */
546*c4762a1bSJed Brown   ierr = DMDAVecGetArray(da,Xglobal,&X);CHKERRQ(ierr);
547*c4762a1bSJed Brown 
548*c4762a1bSJed Brown   /* Get local grid boundaries */
549*c4762a1bSJed Brown   ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
550*c4762a1bSJed Brown 
551*c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
552*c4762a1bSJed Brown 
553*c4762a1bSJed Brown   if (user->init == 1) {
554*c4762a1bSJed Brown     for (j=ys; j<ys+ym; j++) {
555*c4762a1bSJed Brown       for (i=xs; i<xs+xm; i++) {
556*c4762a1bSJed Brown         X[j][i].Ts = user->Ts - i*0.0001;
557*c4762a1bSJed Brown         X[j][i].Ta = X[j][i].Ts - 5;
558*c4762a1bSJed Brown         X[j][i].u  = 0;
559*c4762a1bSJed Brown         X[j][i].v  = 0;
560*c4762a1bSJed Brown         X[j][i].p  = 1.25;
561*c4762a1bSJed Brown         if ((j == 5 || j == 6) && (i == 4 || i == 5))   X[j][i].p += 0.00001;
562*c4762a1bSJed Brown         if ((j == 5 || j == 6) && (i == 12 || i == 13)) X[j][i].p += 0.00001;
563*c4762a1bSJed Brown       }
564*c4762a1bSJed Brown     }
565*c4762a1bSJed Brown   } else {
566*c4762a1bSJed Brown     for (j=ys; j<ys+ym; j++) {
567*c4762a1bSJed Brown       for (i=xs; i<xs+xm; i++) {
568*c4762a1bSJed Brown         X[j][i].Ts = user->Ts;
569*c4762a1bSJed Brown         X[j][i].Ta = X[j][i].Ts - 5;
570*c4762a1bSJed Brown         X[j][i].u  = 0;
571*c4762a1bSJed Brown         X[j][i].v  = 0;
572*c4762a1bSJed Brown         X[j][i].p  = 1.25;
573*c4762a1bSJed Brown       }
574*c4762a1bSJed Brown     }
575*c4762a1bSJed Brown   }
576*c4762a1bSJed Brown 
577*c4762a1bSJed Brown   /* Restore vectors */
578*c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da,Xglobal,&X);CHKERRQ(ierr);
579*c4762a1bSJed Brown   PetscFunctionReturn(0);
580*c4762a1bSJed Brown }
581*c4762a1bSJed Brown 
582*c4762a1bSJed Brown /*
583*c4762a1bSJed Brown    RhsFunc - Evaluates nonlinear function F(u).
584*c4762a1bSJed Brown 
585*c4762a1bSJed Brown    Input Parameters:
586*c4762a1bSJed Brown .  ts - the TS context
587*c4762a1bSJed Brown .  t - current time
588*c4762a1bSJed Brown .  Xglobal - input vector
589*c4762a1bSJed Brown .  F - output vector
590*c4762a1bSJed Brown .  ptr - optional user-defined context, as set by SNESSetFunction()
591*c4762a1bSJed Brown 
592*c4762a1bSJed Brown    Output Parameter:
593*c4762a1bSJed Brown .  F - rhs function vector
594*c4762a1bSJed Brown  */
595*c4762a1bSJed Brown PetscErrorCode RhsFunc(TS ts,PetscReal t,Vec Xglobal,Vec F,void *ctx)
596*c4762a1bSJed Brown {
597*c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ctx;       /* user-defined application context */
598*c4762a1bSJed Brown   DM             da    = user->da;
599*c4762a1bSJed Brown   PetscErrorCode ierr;
600*c4762a1bSJed Brown   PetscInt       i,j,Mx,My,xs,ys,xm,ym;
601*c4762a1bSJed Brown   PetscReal      dhx,dhy;
602*c4762a1bSJed Brown   Vec            localT;
603*c4762a1bSJed Brown   Field          **X,**Frhs;                            /* structures that contain variables of interest and left hand side of governing equations respectively */
604*c4762a1bSJed Brown   PetscScalar    csoil          = user->csoil;          /* heat constant for layer */
605*c4762a1bSJed Brown   PetscScalar    dzlay          = user->dzlay;          /* thickness of top soil layer */
606*c4762a1bSJed Brown   PetscScalar    emma           = user->emma;           /* emission parameter */
607*c4762a1bSJed Brown   PetscScalar    wind           = user->wind;           /* wind speed */
608*c4762a1bSJed Brown   PetscScalar    dewtemp        = user->dewtemp;        /* dew point temperature (moisture in air) */
609*c4762a1bSJed Brown   PetscScalar    pressure1      = user->pressure1;      /* sea level pressure */
610*c4762a1bSJed Brown   PetscScalar    airtemp        = user->airtemp;        /* temperature of air near boundary layer inversion */
611*c4762a1bSJed Brown   PetscScalar    fract          = user->fract;          /* fraction of the sky covered by clouds */
612*c4762a1bSJed Brown   PetscScalar    Tc             = user->Tc;             /* temperature at base of lowest cloud layer */
613*c4762a1bSJed Brown   PetscScalar    lat            = user->lat;            /* latitude */
614*c4762a1bSJed Brown   PetscScalar    Cp             = 1005.7;               /* specific heat of air at constant pressure */
615*c4762a1bSJed Brown   PetscScalar    Rd             = 287.058;              /* gas constant for dry air */
616*c4762a1bSJed Brown   PetscScalar    diffconst      = 1000;                 /* diffusion coefficient */
617*c4762a1bSJed Brown   PetscScalar    f              = 2*0.0000727*PetscSinScalar(lat); /* coriolis force */
618*c4762a1bSJed Brown   PetscScalar    deep_grnd_temp = user->deep_grnd_temp; /* temp in lowest ground layer */
619*c4762a1bSJed Brown   PetscScalar    Ts,u,v,p;
620*c4762a1bSJed Brown   PetscScalar    u_abs,u_plus,u_minus,v_abs,v_plus,v_minus;
621*c4762a1bSJed Brown 
622*c4762a1bSJed Brown   PetscScalar sfctemp1,fsfc1,Ra;
623*c4762a1bSJed Brown   PetscScalar sheat;                   /* sensible heat flux */
624*c4762a1bSJed Brown   PetscScalar latentheat;              /* latent heat flux */
625*c4762a1bSJed Brown   PetscScalar groundflux;              /* flux from conduction of deep ground layer in contact with top soil */
626*c4762a1bSJed Brown   PetscInt    xend,yend;
627*c4762a1bSJed Brown 
628*c4762a1bSJed Brown   PetscFunctionBeginUser;
629*c4762a1bSJed Brown   ierr = DMGetLocalVector(da,&localT);CHKERRQ(ierr);
630*c4762a1bSJed Brown   ierr = DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);CHKERRQ(ierr);
631*c4762a1bSJed Brown 
632*c4762a1bSJed Brown   dhx = (PetscReal)(Mx-1)/(5000*(Mx-1));  /* dhx = 1/dx; assume 2D space domain: [0.0, 1.e5] x [0.0, 1.e5] */
633*c4762a1bSJed Brown   dhy = (PetscReal)(My-1)/(5000*(Mx-1));  /* dhy = 1/dy; */
634*c4762a1bSJed Brown 
635*c4762a1bSJed Brown 
636*c4762a1bSJed Brown   /*
637*c4762a1bSJed Brown      Scatter ghost points to local vector,using the 2-step process
638*c4762a1bSJed Brown         DAGlobalToLocalBegin(),DAGlobalToLocalEnd().
639*c4762a1bSJed Brown      By placing code between these two statements, computations can be
640*c4762a1bSJed Brown      done while messages are in transition.
641*c4762a1bSJed Brown   */
642*c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(da,Xglobal,INSERT_VALUES,localT);CHKERRQ(ierr);
643*c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(da,Xglobal,INSERT_VALUES,localT);CHKERRQ(ierr);
644*c4762a1bSJed Brown 
645*c4762a1bSJed Brown   /* Get pointers to vector data */
646*c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(da,localT,&X);CHKERRQ(ierr);
647*c4762a1bSJed Brown   ierr = DMDAVecGetArray(da,F,&Frhs);CHKERRQ(ierr);
648*c4762a1bSJed Brown 
649*c4762a1bSJed Brown   /* Get local grid boundaries */
650*c4762a1bSJed Brown   ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
651*c4762a1bSJed Brown 
652*c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
653*c4762a1bSJed Brown   /* the interior points */
654*c4762a1bSJed Brown   xend=xs+xm; yend=ys+ym;
655*c4762a1bSJed Brown   for (j=ys; j<yend; j++) {
656*c4762a1bSJed Brown     for (i=xs; i<xend; i++) {
657*c4762a1bSJed Brown       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; */
658*c4762a1bSJed Brown 
659*c4762a1bSJed Brown       sfctemp1 = (double)Ts;
660*c4762a1bSJed Brown       ierr     = calcfluxs(sfctemp1,airtemp,emma,fract,Tc,&fsfc1);CHKERRQ(ierr);        /* calculates surface net radiative flux */
661*c4762a1bSJed Brown       ierr     = sensibleflux(sfctemp1,airtemp,wind,&sheat);CHKERRQ(ierr);              /* calculate sensible heat flux */
662*c4762a1bSJed Brown       ierr     = latentflux(sfctemp1,dewtemp,wind,pressure1,&latentheat);CHKERRQ(ierr); /* calculates latent heat flux */
663*c4762a1bSJed Brown       ierr     = calc_gflux(sfctemp1,deep_grnd_temp,&groundflux);CHKERRQ(ierr);         /* calculates flux from earth below surface soil layer by conduction */
664*c4762a1bSJed Brown       ierr     = calcfluxa(sfctemp1,airtemp,emma,&Ra);CHKERRQ(ierr);                    /* Calculates the change in downward radiative flux */
665*c4762a1bSJed Brown       fsfc1    = fsfc1 + latentheat + sheat + groundflux;                               /* adds radiative, sensible heat, latent heat, and ground heat flux yielding net flux */
666*c4762a1bSJed Brown 
667*c4762a1bSJed Brown       /* convective coefficients for upwinding */
668*c4762a1bSJed Brown       u_abs   = PetscAbsScalar(u);
669*c4762a1bSJed Brown       u_plus  = .5*(u + u_abs); /* u if u>0; 0 if u<0 */
670*c4762a1bSJed Brown       u_minus = .5*(u - u_abs); /* u if u <0; 0 if u>0 */
671*c4762a1bSJed Brown 
672*c4762a1bSJed Brown       v_abs   = PetscAbsScalar(v);
673*c4762a1bSJed Brown       v_plus  = .5*(v + v_abs); /* v if v>0; 0 if v<0 */
674*c4762a1bSJed Brown       v_minus = .5*(v - v_abs); /* v if v <0; 0 if v>0 */
675*c4762a1bSJed Brown 
676*c4762a1bSJed Brown       /* Solve governing equations */
677*c4762a1bSJed Brown       /* P = p*Rd*Ts; */
678*c4762a1bSJed Brown 
679*c4762a1bSJed Brown       /* du/dt -> time change of east-west component of the wind */
680*c4762a1bSJed Brown       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) */
681*c4762a1bSJed Brown                      - v_plus*(u - X[j-1][i].u)*dhy - v_minus*(X[j+1][i].u - u)*dhy       /* - v(du/dy) */
682*c4762a1bSJed Brown                      -(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)] */
683*c4762a1bSJed Brown /*                     -(1/p)*(X[j][i+1].P - X[j][i-1].P)*dhx */
684*c4762a1bSJed Brown                      + f*v;
685*c4762a1bSJed Brown 
686*c4762a1bSJed Brown       /* dv/dt -> time change of north-south component of the wind */
687*c4762a1bSJed Brown       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) */
688*c4762a1bSJed Brown                      - v_plus*(v - X[j-1][i].v)*dhy - v_minus*(X[j+1][i].v - v)*dhy       /* - v(dv/dy) */
689*c4762a1bSJed Brown                      -(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)] */
690*c4762a1bSJed Brown /*                     -(1/p)*(X[j+1][i].P - X[j-1][i].P)*dhy */
691*c4762a1bSJed Brown                      -f*u;
692*c4762a1bSJed Brown 
693*c4762a1bSJed Brown       /* dT/dt -> time change of temperature */
694*c4762a1bSJed Brown       Frhs[j][i].Ts = (fsfc1/(csoil*dzlay))                                            /* Fnet/(Cp*dz)  diabatic change in T */
695*c4762a1bSJed Brown                       -u_plus*(Ts - X[j][i-1].Ts)*dhx - u_minus*(X[j][i+1].Ts - Ts)*dhx  /* - u*(dTs/dx)  advection x */
696*c4762a1bSJed Brown                       -v_plus*(Ts - X[j-1][i].Ts)*dhy - v_minus*(X[j+1][i].Ts - Ts)*dhy  /* - v*(dTs/dy)  advection y */
697*c4762a1bSJed Brown                       + diffconst*((X[j][i+1].Ts - 2*Ts + X[j][i-1].Ts)*dhx*dhx               /* + D(Ts_xx + Ts_yy)  diffusion */
698*c4762a1bSJed Brown                                    + (X[j+1][i].Ts - 2*Ts + X[j-1][i].Ts)*dhy*dhy);
699*c4762a1bSJed Brown 
700*c4762a1bSJed Brown       /* dp/dt -> time change of */
701*c4762a1bSJed Brown       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) */
702*c4762a1bSJed Brown                      -v_plus*(p - X[j-1][i].p)*dhy - v_minus*(X[j+1][i].p - p)*dhy;    /* - v*(dp/dy) */
703*c4762a1bSJed Brown 
704*c4762a1bSJed Brown       Frhs[j][i].Ta = Ra/Cp;  /* dTa/dt time change of air temperature */
705*c4762a1bSJed Brown     }
706*c4762a1bSJed Brown   }
707*c4762a1bSJed Brown 
708*c4762a1bSJed Brown   /* Restore vectors */
709*c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(da,localT,&X);CHKERRQ(ierr);
710*c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da,F,&Frhs);CHKERRQ(ierr);
711*c4762a1bSJed Brown   ierr = DMRestoreLocalVector(da,&localT);CHKERRQ(ierr);
712*c4762a1bSJed Brown   PetscFunctionReturn(0);
713*c4762a1bSJed Brown }
714*c4762a1bSJed Brown 
715*c4762a1bSJed Brown PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal time,Vec T,void *ctx)
716*c4762a1bSJed Brown {
717*c4762a1bSJed Brown   PetscErrorCode    ierr;
718*c4762a1bSJed Brown   const PetscScalar *array;
719*c4762a1bSJed Brown   MonitorCtx        *user  = (MonitorCtx*)ctx;
720*c4762a1bSJed Brown   PetscViewer       viewer = user->drawviewer;
721*c4762a1bSJed Brown   PetscMPIInt       rank;
722*c4762a1bSJed Brown   PetscReal         norm;
723*c4762a1bSJed Brown 
724*c4762a1bSJed Brown   PetscFunctionBeginUser;
725*c4762a1bSJed Brown   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)ts),&rank);CHKERRQ(ierr);
726*c4762a1bSJed Brown   ierr = VecNorm(T,NORM_INFINITY,&norm);CHKERRQ(ierr);
727*c4762a1bSJed Brown 
728*c4762a1bSJed Brown   if (step%user->interval == 0) {
729*c4762a1bSJed Brown     ierr = VecGetArrayRead(T,&array);CHKERRQ(ierr);
730*c4762a1bSJed Brown     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]);
731*c4762a1bSJed Brown     ierr = VecRestoreArrayRead(T,&array);CHKERRQ(ierr);
732*c4762a1bSJed Brown   }
733*c4762a1bSJed Brown 
734*c4762a1bSJed Brown   if (user->drawcontours) {
735*c4762a1bSJed Brown     ierr = VecView(T,viewer);CHKERRQ(ierr);
736*c4762a1bSJed Brown   }
737*c4762a1bSJed Brown   PetscFunctionReturn(0);
738*c4762a1bSJed Brown }
739*c4762a1bSJed Brown 
740*c4762a1bSJed Brown 
741*c4762a1bSJed Brown 
742*c4762a1bSJed Brown /*TEST
743*c4762a1bSJed Brown 
744*c4762a1bSJed Brown    build:
745*c4762a1bSJed Brown       requires: !complex !single
746*c4762a1bSJed Brown 
747*c4762a1bSJed Brown    test:
748*c4762a1bSJed Brown       args: -ts_max_steps 130 -monitor_interval 60
749*c4762a1bSJed Brown       output_file: output/ex5.out
750*c4762a1bSJed Brown       requires: !complex !single
751*c4762a1bSJed Brown       localrunfiles: ex5_control.txt
752*c4762a1bSJed Brown 
753*c4762a1bSJed Brown    test:
754*c4762a1bSJed Brown       suffix: 2
755*c4762a1bSJed Brown       nsize: 4
756*c4762a1bSJed Brown       args: -ts_max_steps 130 -monitor_interval 60
757*c4762a1bSJed Brown       output_file: output/ex5.out
758*c4762a1bSJed Brown       localrunfiles: ex5_control.txt
759*c4762a1bSJed Brown       requires: !complex !single
760*c4762a1bSJed Brown 
761*c4762a1bSJed Brown TEST*/
762