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