surface_th_aux.F90 coverage: 83.33 %func 84.06 %block
1) module Surface_TH_Aux_module
2)
3) use PFLOTRAN_Constants_module
4)
5) implicit none
6)
7) private
8)
9) #include "petsc/finclude/petscsys.h"
10)
11) type, public :: Surface_TH_auxvar_type
12) PetscReal :: h ! enthalpy -- not currently used
13) PetscReal :: u ! internal energy -- not currently used
14) PetscReal :: pc ! pressure change -- not currently used
15) PetscReal :: Cw ! Specific heat capacity of surface water
16) PetscReal :: Ci ! Specific heat capacity of surface ice
17) PetscReal :: Cwi ! Weighted average of Cw and Ci
18) ! RTM: Note that I believe we can simple make Cw and Ci Fortran
19) ! parameters, but I'm keeping things as is for now.
20) PetscReal :: k_therm ! Thermal conductivity of surface water
21) PetscReal :: unfrozen_fraction ! Proportion of unfrozen surface water.
22) PetscReal :: den_water_kg ! Density [kg/m^3] of liquid water ONLY.
23) ! Note that we currently use the den_kg(1) field of the global aux var type
24) ! to store the density of the liquid/ice mixture. We also need to track
25) ! the liquid density because it is required in in the advective term of
26) ! the temperature equation.
27) end type Surface_TH_auxvar_type
28)
29) type, public :: Surface_TH_type
30) PetscInt :: n_zero_rows
31) PetscInt, pointer :: zero_rows_local(:), zero_rows_local_ghosted(:)
32) PetscBool :: auxvars_up_to_date
33) PetscBool :: inactive_cells_exist
34) PetscInt :: num_aux, num_aux_bc, num_aux_ss
35) type(Surface_TH_auxvar_type), pointer :: auxvars(:)
36) type(Surface_TH_auxvar_type), pointer :: auxvars_bc(:)
37) type(Surface_TH_auxvar_type), pointer :: auxvars_ss(:)
38) end type Surface_TH_type
39)
40) public :: SurfaceTHAuxCreate, &
41) SurfaceTHAuxDestroy, &
42) SurfaceTHAuxVarCompute, &
43) SurfaceTHAuxVarComputeUnfrozen, &
44) SurfaceTHAuxVarInit, &
45) SurfaceTHAuxVarCopy
46)
47) contains
48)
49) ! ************************************************************************** !
50)
51) function SurfaceTHAuxCreate(option)
52) !
53) ! This routine creates an empty object
54) !
55) ! Author: Gautam Bisht, LBNL
56) ! Date: 03/07/13
57) !
58)
59) use Option_module
60)
61) implicit none
62)
63) type(option_type) :: option
64) type(Surface_TH_type), pointer :: SurfaceTHAuxCreate
65)
66) type(Surface_TH_type), pointer :: aux
67)
68) allocate(aux)
69) aux%auxvars_up_to_date = PETSC_FALSE
70) aux%inactive_cells_exist = PETSC_FALSE
71) aux%num_aux = 0
72) aux%num_aux_bc = 0
73) aux%num_aux_ss = 0
74) nullify(aux%auxvars)
75) nullify(aux%auxvars_bc)
76) nullify(aux%auxvars_ss)
77) aux%n_zero_rows = 0
78) nullify(aux%zero_rows_local)
79) nullify(aux%zero_rows_local_ghosted)
80)
81) SurfaceTHAuxCreate => aux
82)
83) end function SurfaceTHAuxCreate
84)
85) ! ************************************************************************** !
86)
87) subroutine SurfaceTHAuxVarInit(auxvar,option)
88) !
89) ! This routine initilizes an object.
90) !
91) ! Author: Gautam Bisht, LBNL
92) ! Date: 03/07/13
93) !
94)
95) use Option_module
96) use PFLOTRAN_Constants_module, only : UNINITIALIZED_DOUBLE
97)
98) implicit none
99)
100) type(Surface_TH_auxvar_type) :: auxvar
101) type(option_type) :: option
102)
103) PetscReal :: uninit_value
104) uninit_value = UNINITIALIZED_DOUBLE
105)
106) auxvar%h = uninit_value
107) auxvar%u = uninit_value
108) auxvar%pc = uninit_value
109) auxvar%Cw = 4.188d3 ! [J/kg/K]
110) !auxvar%Ci = 2.050d3 ! [J/kg/K]
111) auxvar%Ci = 4.188d3 ! [J/kg/K]
112) auxvar%Cwi = 4.188d3 ! [J/kg/K]
113) auxvar%k_therm = 0.57d0 ! [J/s/m/K]
114) auxvar%unfrozen_fraction = 1.d0
115) auxvar%den_water_kg = 1.d3 ! [kg/m^3]
116)
117) end subroutine SurfaceTHAuxVarInit
118)
119) ! ************************************************************************** !
120)
121) subroutine SurfaceTHAuxVarCopy(auxvar,auxvar2,option)
122) !
123) ! This routine makes a copy of an object.
124) !
125) ! Author: Gautam Bisht, LBNL
126) ! Date: 03/07/13
127) !
128)
129) use Option_module
130)
131) implicit none
132)
133) type(Surface_TH_auxvar_type) :: auxvar, auxvar2
134) type(option_type) :: option
135)
136) auxvar2%h = auxvar%h
137) auxvar2%u = auxvar%u
138) auxvar2%pc = auxvar%pc
139) auxvar2%Cw = auxvar%Cw
140) auxvar2%Ci = auxvar%Ci
141) auxvar2%Cwi = auxvar%Cwi
142) auxvar2%k_therm = auxvar%k_therm
143) auxvar2%unfrozen_fraction = auxvar%unfrozen_fraction
144) auxvar2%den_water_kg = auxvar%den_water_kg
145)
146) end subroutine SurfaceTHAuxVarCopy
147)
148) ! ************************************************************************** !
149)
150) subroutine SurfaceTHAuxVarCompute(xx,auxvar,global_auxvar, &
151) option)
152) !
153) ! This routine computes values for auxiliary variables.
154) !
155) ! Author: Gautam Bisht, LBNL
156) ! Date: 03/07/13
157) !
158)
159) use Option_module
160) use Surface_Global_Aux_module
161)
162) use EOS_Water_module
163) use PFLOTRAN_Constants_module, only : DUMMY_VALUE,MIN_SURFACE_WATER_HEIGHT
164)
165) implicit none
166)
167) type(option_type) :: option
168) PetscReal :: xx(option%nflowdof)
169) type(Surface_TH_auxvar_type) :: auxvar
170) type(surface_global_auxvar_type) :: global_auxvar
171) PetscReal :: por, perm
172)
173) PetscErrorCode :: ierr
174) PetscReal :: pw,dw_kg,dw_mol,hw,sat_pressure,visl
175) PetscReal :: di_kg, dwi_kg
176) ! Densities of ice and water-ice mixture, respectively.
177) PetscReal :: unfrozen_fraction
178) PetscReal :: kr, ds_dp, dkr_dp
179) PetscReal :: dvis_dt, dvis_dp, dvis_dpsat
180) PetscReal :: dw_dp, dw_dt, hw_dp, hw_dt
181) PetscReal :: dpw_dp
182) PetscReal :: dpsat_dt
183) PetscReal :: k_therm_w, k_therm_i
184)
185) global_auxvar%den_kg(1) = 0.d0
186)
187) auxvar%h = 0.d0
188) auxvar%u = 0.d0
189) kr = 0.d0
190)
191) global_auxvar%head(1) = xx(1)
192) !global_auxvar%temp = xx(2)
193) ! RTM: Why is the above commented out? Is one of these internal
194) ! energy instead of temperature?
195)
196)
197) !*************** Liquid phase properties **************************
198)
199) pw = option%reference_pressure
200) ds_dp = 0.d0
201) dkr_dp = 0.d0
202)
203) if (global_auxvar%head(1) < MIN_SURFACE_WATER_HEIGHT) then
204) global_auxvar%is_dry = PETSC_TRUE
205) call EOSWaterDensity(0.0d0,pw,dw_kg,dw_mol,ierr)
206) call EOSWaterEnthalpy(0.0d0,pw,hw,ierr)
207) else
208) global_auxvar%is_dry = PETSC_FALSE
209) call EOSWaterDensity(global_auxvar%temp,pw,dw_kg,dw_mol,ierr)
210) call EOSWaterEnthalpy(global_auxvar%temp,pw,hw,ierr)
211) endif
212)
213) ! J/kmol -> whatever units
214) hw = hw * option%scale
215)
216) global_auxvar%den_kg(1) = dw_kg
217) !di_kg = 917.d0 ![kg/m^3]
218) di_kg = dw_kg
219)
220) k_therm_w = 0.57d0 ! [J/s/m/K]
221) k_therm_i = 2.18d0 ! [J/s/m/K]
222) ! RTM: Same warning for thermal conductivities--these should be computed.
223)
224) ! RTM: These are being set but we are not using them now. We should get rid
225) ! of them if we settle on not using an enthalpy formulation for the energy
226) ! equation.
227) auxvar%h = hw
228) auxvar%u = auxvar%h - pw / dw_mol * option%scale
229)
230) ! Compute unfrozen fraction, and then compute the weighted averages of
231) ! density, specific heat capacity, thermal conductivity
232) unfrozen_fraction = SurfaceTHAuxVarComputeUnfrozen(global_auxvar%temp)
233) auxvar%unfrozen_fraction = unfrozen_fraction
234) global_auxvar%den_kg(1) = unfrozen_fraction * dw_kg + (1.d0 - unfrozen_fraction) * di_kg
235) auxvar%Cwi = unfrozen_fraction * auxvar%Cw + (1.d0 - unfrozen_fraction) * auxvar%Ci
236) auxvar%k_therm = unfrozen_fraction * k_therm_w + (1.d0 - unfrozen_fraction) * k_therm_i
237)
238)
239) end subroutine SurfaceTHAuxVarCompute
240)
241) ! ************************************************************************** !
242)
243) function SurfaceTHAuxVarComputeUnfrozen(temp)
244) !
245) ! This returns the unfrozen fraction, given a temperature.
246) ! RTM: This is currently a public function, but is only called within
247) ! Surface_TH_Aux_module, so we may want to make it private.
248) ! We may also want to move this into its own module if we add some
249) ! more sophisticated ways of computing this quantity.
250) !
251)
252) implicit none
253)
254) PetscReal :: SurfaceTHAuxVarComputeUnfrozen
255)
256) PetscReal :: temp
257)
258) ! For now, we simply say that everything is unfrozen if the temperature is
259) ! above the nominal freezing point, and everything is frozen otherwise.
260) ! We may want to turn this step function into a steep sigmoidal one for
261) ! numerical reasons.
262) ! At some point, we may also want to consider a more detailed model of the
263) ! freezing and thawing process that might partition the frozen and unfrozen
264) ! phases according to a more physical mechanism (maybe using state transition
265) ! theory?).
266) if (temp > 0.d0) then
267) SurfaceTHAuxVarComputeUnfrozen = 1.d0
268) else
269) SurfaceTHAuxVarComputeUnfrozen = 0.d0
270) endif
271)
272) end function SurfaceTHAuxVarComputeUnfrozen
273)
274) ! ************************************************************************** !
275)
276) subroutine SurfaceTHAuxDestroy(aux)
277) !
278) ! This routine deallocates an object.
279) !
280) ! Author: Gautam Bisht, LBNL
281) ! Date: 03/07/13
282) !
283)
284) implicit none
285)
286) type(Surface_TH_type), pointer :: aux
287) PetscInt :: iaux
288)
289) if (.not.associated(aux)) return
290)
291) if (associated(aux%auxvars)) deallocate(aux%auxvars)
292) nullify(aux%auxvars)
293) if (associated(aux%auxvars_bc)) deallocate(aux%auxvars_bc)
294) nullify(aux%auxvars_bc)
295) if (associated(aux%auxvars_ss)) deallocate(aux%auxvars_ss)
296) nullify(aux%auxvars_ss)
297) if (associated(aux%zero_rows_local)) deallocate(aux%zero_rows_local)
298) nullify(aux%zero_rows_local)
299) if (associated(aux%zero_rows_local_ghosted)) deallocate(aux%zero_rows_local_ghosted)
300) nullify(aux%zero_rows_local_ghosted)
301)
302) deallocate(aux)
303) nullify(aux)
304)
305) end subroutine SurfaceTHAuxDestroy
306)
307) end module Surface_TH_Aux_module