richards_aux.F90 coverage: 100.00 %func 91.82 %block
1) module Richards_Aux_module
2)
3) #ifdef BUFFER_MATRIX
4) use Matrix_Buffer_module
5) #endif
6)
7) use PFLOTRAN_Constants_module
8)
9) implicit none
10)
11) private
12)
13) #include "petsc/finclude/petscsys.h"
14)
15) PetscReal, public :: richards_itol_scaled_res = 1.d-5
16) PetscReal, public :: richards_itol_rel_update = UNINITIALIZED_DOUBLE
17)
18) type, public :: richards_auxvar_type
19)
20) PetscReal :: pc
21) #ifdef USE_ANISOTROPIC_MOBILITY
22) PetscReal :: kvr_x
23) PetscReal :: kvr_y
24) PetscReal :: kvr_z
25) PetscReal :: dkvr_x_dp
26) PetscReal :: dkvr_y_dp
27) PetscReal :: dkvr_z_dp
28) #else
29) PetscReal :: kvr
30) PetscReal :: dkvr_dp
31) #endif
32) PetscReal :: dsat_dp
33) PetscReal :: dden_dp
34)
35) ! OLD-VAR-NAMES = NEW-VAR
36) ! ------------------------------------------------
37) ! P_min = vars_for_sflow(1)
38) ! P_max = vars_for_sflow(2)
39) ! coeff_for_cubic_approx = vars_for_sflow(3:6)
40) ! range_for_linear_approx = vars_for_sflow(7:10)
41) ! bcflux_default_scheme = vars_for_sflow(11)
42) PetscReal, pointer :: vars_for_sflow(:)
43)
44) end type richards_auxvar_type
45)
46) type, public :: richards_type
47) PetscInt :: n_zero_rows
48) PetscInt, pointer :: zero_rows_local(:), zero_rows_local_ghosted(:)
49)
50) PetscBool :: auxvars_up_to_date
51) PetscBool :: auxvars_cell_pressures_up_to_date
52) PetscBool :: inactive_cells_exist
53) PetscInt :: num_aux, num_aux_bc, num_aux_ss
54) type(richards_auxvar_type), pointer :: auxvars(:)
55) type(richards_auxvar_type), pointer :: auxvars_bc(:)
56) type(richards_auxvar_type), pointer :: auxvars_ss(:)
57) #ifdef BUFFER_MATRIX
58) type(matrix_buffer_type), pointer :: matrix_buffer
59) #endif
60) end type richards_type
61)
62) public :: RichardsAuxCreate, RichardsAuxDestroy, &
63) RichardsAuxVarCompute, RichardsAuxVarInit, &
64) RichardsAuxVarCopy
65)
66) contains
67)
68) ! ************************************************************************** !
69)
70) function RichardsAuxCreate()
71) !
72) ! Allocate and initialize auxiliary object
73) !
74) ! Author: Glenn Hammond
75) ! Date: 02/14/08
76) !
77)
78) use Option_module
79)
80) implicit none
81)
82) type(richards_type), pointer :: RichardsAuxCreate
83)
84) type(richards_type), pointer :: aux
85)
86) allocate(aux)
87) aux%auxvars_up_to_date = PETSC_FALSE
88) aux%auxvars_cell_pressures_up_to_date = PETSC_FALSE
89) aux%inactive_cells_exist = PETSC_FALSE
90) aux%num_aux = 0
91) aux%num_aux_bc = 0
92) aux%num_aux_ss = 0
93) nullify(aux%auxvars)
94) nullify(aux%auxvars_bc)
95) nullify(aux%auxvars_ss)
96) aux%n_zero_rows = 0
97) nullify(aux%zero_rows_local)
98) nullify(aux%zero_rows_local_ghosted)
99) #ifdef BUFFER_MATRIX
100) nullify(aux%matrix_buffer)
101) #endif
102)
103) RichardsAuxCreate => aux
104)
105) end function RichardsAuxCreate
106)
107) ! ************************************************************************** !
108)
109) subroutine RichardsAuxVarInit(auxvar,option)
110) !
111) ! Initialize auxiliary object
112) !
113) ! Author: Glenn Hammond
114) ! Date: 02/14/08
115) !
116)
117) use Option_module
118)
119) implicit none
120)
121) type(richards_auxvar_type) :: auxvar
122) type(option_type) :: option
123)
124) auxvar%pc = 0.d0
125)
126) #ifdef USE_ANISOTROPIC_MOBILITY
127) auxvar%kvr_x = 0.d0
128) auxvar%kvr_y = 0.d0
129) auxvar%kvr_z = 0.d0
130) auxvar%dkvr_x_dp = 0.d0
131) auxvar%dkvr_y_dp = 0.d0
132) auxvar%dkvr_z_dp = 0.d0
133) #else
134) auxvar%kvr = 0.d0
135) auxvar%dkvr_dp = 0.d0
136) #endif
137)
138) auxvar%dsat_dp = 0.d0
139) auxvar%dden_dp = 0.d0
140)
141) if (option%surf_flow_on) then
142) allocate(auxvar%vars_for_sflow(11))
143) auxvar%vars_for_sflow(:) = 0.d0
144) else
145) nullify(auxvar%vars_for_sflow)
146) endif
147)
148)
149) end subroutine RichardsAuxVarInit
150)
151) ! ************************************************************************** !
152)
153) subroutine RichardsAuxVarCopy(auxvar,auxvar2,option)
154) !
155) ! Copies an auxiliary variable
156) !
157) ! Author: Glenn Hammond
158) ! Date: 12/13/07
159) !
160)
161) use Option_module
162)
163) implicit none
164)
165) type(richards_auxvar_type) :: auxvar, auxvar2
166) type(option_type) :: option
167)
168) auxvar2%pc = auxvar%pc
169)
170) #ifdef USE_ANISOTROPIC_MOBILITY
171) auxvar2%kvr_x = auxvar%kvr_x
172) auxvar2%kvr_y = auxvar%kvr_y
173) auxvar2%kvr_z = auxvar%kvr_z
174) auxvar2%dkvr_x_dp = auxvar%dkvr_x_dp
175) auxvar2%dkvr_y_dp = auxvar%dkvr_y_dp
176) auxvar2%dkvr_z_dp = auxvar%dkvr_z_dp
177) #else
178) auxvar2%kvr = auxvar%kvr
179) auxvar2%dkvr_dp = auxvar%dkvr_dp
180) #endif
181)
182) auxvar2%dsat_dp = auxvar%dsat_dp
183) auxvar2%dden_dp = auxvar%dden_dp
184)
185) if (option%surf_flow_on) &
186) auxvar2%vars_for_sflow(:) = auxvar%vars_for_sflow(:)
187)
188) end subroutine RichardsAuxVarCopy
189)
190) ! ************************************************************************** !
191)
192) subroutine RichardsAuxVarCompute(x,auxvar,global_auxvar,material_auxvar, &
193) characteristic_curves,option)
194) !
195) ! Computes auxiliary variables for each grid cell
196) !
197) ! Author: Glenn Hammond
198) ! Date: 02/22/08
199) !
200)
201) use Option_module
202) use Global_Aux_module
203)
204) use EOS_Water_module
205) use Characteristic_Curves_module
206) use Material_Aux_class
207)
208) implicit none
209)
210) type(option_type) :: option
211) class(characteristic_curves_type) :: characteristic_curves
212) PetscReal :: x(option%nflowdof)
213) type(richards_auxvar_type) :: auxvar
214) type(global_auxvar_type) :: global_auxvar
215) class(material_auxvar_type) :: material_auxvar
216)
217) PetscInt :: i
218) PetscBool :: saturated
219) PetscErrorCode :: ierr
220) PetscReal :: pw,dw_kg,dw_mol,hw,sat_pressure,visl
221) PetscReal :: kr, ds_dp, dkr_dp
222) PetscReal :: dvis_dt, dvis_dp
223) PetscReal :: dw_dp, dw_dt, hw_dp, hw_dt
224) PetscReal :: pert, pw_pert, dw_kg_pert
225) PetscReal :: fs, ani_A, ani_B, ani_C, ani_n, ani_coef
226) PetscReal :: dkr_sat
227) PetscReal :: aux(1)
228) PetscReal, parameter :: tol = 1.d-3
229)
230) global_auxvar%sat = 0.d0
231) global_auxvar%den = 0.d0
232) global_auxvar%den_kg = 0.d0
233)
234) #ifdef USE_ANISOTROPIC_MOBILITY
235) auxvar%kvr_x = 0.d0
236) auxvar%kvr_y = 0.d0
237) auxvar%kvr_z = 0.d0
238) #else
239) auxvar%kvr = 0.d0
240) #endif
241)
242) kr = 0.d0
243)
244) global_auxvar%pres = x(1)
245) global_auxvar%temp = option%reference_temperature
246)
247) auxvar%pc = option%reference_pressure - global_auxvar%pres(1)
248)
249) !*************** Liquid phase properties **************************
250) pw = option%reference_pressure
251) ds_dp = 0.d0
252) dkr_dp = 0.d0
253) if (auxvar%pc > 0.d0) then
254) saturated = PETSC_FALSE
255) call characteristic_curves%saturation_function% &
256) Saturation(auxvar%pc, &
257) global_auxvar%sat(1), &
258) ds_dp,option)
259) ! if ds_dp is 0, we consider the cell saturated.
260) if (ds_dp < 1.d-40) then
261) saturated = PETSC_TRUE
262) else
263) call characteristic_curves%liq_rel_perm_function% &
264) RelativePermeability(global_auxvar%sat(1), &
265) kr,dkr_sat,option)
266) dkr_dp = ds_dp * dkr_sat
267) endif
268) else
269) saturated = PETSC_TRUE
270) endif
271)
272) ! the purpose for splitting this condition from the 'else' statement
273) ! above is due to SaturationFunctionCompute switching a cell to
274) ! saturated to prevent unstable (potentially infinite) derivatives when
275) ! capillary pressure is very small
276) if (saturated) then
277) auxvar%pc = 0.d0
278) global_auxvar%sat = 1.d0
279) kr = 1.d0
280) pw = global_auxvar%pres(1)
281) endif
282)
283) if (.not.option%flow%density_depends_on_salinity) then
284) call EOSWaterDensity(global_auxvar%temp,pw,dw_kg,dw_mol, &
285) dw_dp,dw_dt,ierr)
286) ! may need to compute dpsat_dt to pass to VISW
287) call EOSWaterSaturationPressure(global_auxvar%temp,sat_pressure,ierr)
288) !geh: 0.d0 passed in for derivative of pressure w/respect to temp
289) call EOSWaterViscosity(global_auxvar%temp,pw,sat_pressure,0.d0, &
290) visl,dvis_dt,dvis_dp,ierr)
291) else
292) aux(1) = global_auxvar%m_nacl(1)
293) call EOSWaterDensityExt(global_auxvar%temp,pw,aux, &
294) dw_kg,dw_mol,dw_dp,dw_dt,ierr)
295) call EOSWaterViscosityExt(global_auxvar%temp,pw,sat_pressure,0.d0,aux, &
296) visl,dvis_dt,dvis_dp,ierr)
297) endif
298) if (.not.saturated) then !kludge since pw is constant in the unsat zone
299) dvis_dp = 0.d0
300) dw_dp = 0.d0
301) hw_dp = 0.d0
302) endif
303)
304) global_auxvar%den = dw_mol
305) global_auxvar%den_kg = dw_kg
306) auxvar%dsat_dp = ds_dp
307) auxvar%dden_dp = dw_dp
308) auxvar%kvr = kr/visl
309) auxvar%dkvr_dp = dkr_dp/visl - kr/(visl*visl)*dvis_dp
310)
311) end subroutine RichardsAuxVarCompute
312)
313) ! ************************************************************************** !
314)
315) subroutine AuxVarDestroy(auxvar)
316) !
317) ! Deallocates a richards auxiliary object
318) !
319) ! Author: Glenn Hammond
320) ! Date: 02/14/08
321) !
322)
323) implicit none
324)
325) type(richards_auxvar_type) :: auxvar
326)
327) end subroutine AuxVarDestroy
328)
329) ! ************************************************************************** !
330)
331) subroutine RichardsAuxDestroy(aux)
332) !
333) ! Deallocates a richards auxiliary object
334) !
335) ! Author: Glenn Hammond
336) ! Date: 02/14/08
337) !
338) use Utility_module, only : DeallocateArray
339)
340) implicit none
341)
342) type(richards_type), pointer :: aux
343) PetscInt :: iaux
344)
345) if (.not.associated(aux)) return
346)
347) if (associated(aux%auxvars)) then
348) do iaux = 1, aux%num_aux
349) call AuxVarDestroy(aux%auxvars(iaux))
350) enddo
351) deallocate(aux%auxvars)
352) endif
353) nullify(aux%auxvars)
354) if (associated(aux%auxvars_bc)) then
355) do iaux = 1, aux%num_aux_bc
356) call AuxVarDestroy(aux%auxvars_bc(iaux))
357) enddo
358) deallocate(aux%auxvars_bc)
359) endif
360) nullify(aux%auxvars_bc)
361) if (associated(aux%auxvars_ss)) then
362) do iaux = 1, aux%num_aux_ss
363) call AuxVarDestroy(aux%auxvars_ss(iaux))
364) enddo
365) deallocate(aux%auxvars_ss)
366) endif
367) nullify(aux%auxvars_ss)
368)
369) call DeallocateArray(aux%zero_rows_local)
370) call DeallocateArray(aux%zero_rows_local_ghosted)
371)
372) #ifdef BUFFER_MATRIX
373) if (associated(aux%matrix_buffer)) then
374) call MatrixBufferDestroy(aux%matrix_buffer)
375) endif
376) nullify(aux%matrix_buffer)
377) #endif
378)
379) deallocate(aux)
380) nullify(aux)
381)
382) end subroutine RichardsAuxDestroy
383)
384) end module Richards_Aux_module