srcsink_sandbox_downreg.F90 coverage: 100.00 %func 77.78 %block
1) module SrcSink_Sandbox_Downreg_class
2)
3) ! Source/sink Sandbox for downregulating source or sink terms to avoid
4) ! overpressurization (too big a pressure, + or -)
5) ! the source is regulated by multiplying it by
6) ! (1 - x^2)^2
7) ! with
8) ! x = (pressure - pref - pressure_max - delta/2)/delta
9) ! for pressure in between of x*delta and (x+1)*delta
10) ! the sink is regulated by multiplying it by
11) ! 1 - (1 - x^2)^2
12) ! with
13) ! x = (pressure - pref - pressure_min - delta/2)/delta
14) ! for pressure in between of x*delta and (x+1)*delta
15) ! for transpiration sink, the physical interpretation can be that the soil water
16) ! extration rate is decreased from q to 0 as pressure decreases from
17) ! pressure_min - delta/2 to pressure_min + delta/2
18)
19) ! extended from srcsink_sandbox_mass_rate
20) ! currently work in RICHARDS mode
21)
22) use PFLOTRAN_Constants_module
23) use SrcSink_Sandbox_Base_class
24) use Dataset_Base_class
25)
26) implicit none
27)
28) private
29)
30) #include "petsc/finclude/petscsys.h"
31)
32) type, public, &
33) extends(srcsink_sandbox_base_type) :: srcsink_sandbox_downreg_type
34) PetscReal :: pressure_min ! Pa
35) PetscReal :: pressure_max ! Pa
36) PetscReal :: pressure_delta ! Pa
37) class(dataset_base_type), pointer :: dataset
38) contains
39) procedure, public :: ReadInput => DownregRead
40) procedure, public :: Setup => DownregSetup
41) procedure, public :: Update => DownregUpdate
42) procedure, public :: Evaluate => DownregSrcSink
43) procedure, public :: Destroy => DownregDestroy
44) end type srcsink_sandbox_downreg_type
45)
46) public :: DownregCreate
47)
48) contains
49)
50) ! ************************************************************************** !
51)
52) function DownregCreate()
53) !
54) ! Allocates mass rate src/sink object.
55) !
56) ! Author: Guoping Tang
57) ! Date: 06/03/14
58)
59) implicit none
60)
61) class(srcsink_sandbox_downreg_type), pointer :: DownregCreate
62)
63) allocate(DownregCreate)
64) call SSSandboxBaseInit(DownregCreate)
65) DownregCreate%pressure_min = -1.0d9
66) DownregCreate%pressure_max = 1.0d9
67) DownregCreate%pressure_delta = 1.0d4
68) nullify(DownregCreate%dataset)
69)
70) end function DownregCreate
71)
72) ! ************************************************************************** !
73)
74) subroutine DownregRead(this,input,option)
75) !
76) ! Reads input deck for mass rate src/sink parameters
77) !
78) ! Author: Guoping Tang
79) ! Date: 06/03/14
80)
81) use Option_module
82) use String_module
83) use Input_Aux_module
84) use Units_module, only : UnitsConvertToInternal
85) use Condition_module
86) use Dataset_module
87) use Dataset_Ascii_class
88) use Time_Storage_module
89)
90) implicit none
91)
92) class(srcsink_sandbox_downreg_type) :: this
93) type(input_type), pointer :: input
94) type(option_type) :: option
95) class(dataset_ascii_type), pointer :: dataset_ascii
96)
97) PetscInt :: i
98) character(len=MAXWORDLENGTH) :: word
99) character(len=MAXSTRINGLENGTH) :: string
100) character(len=MAXWORDLENGTH) :: units, internal_units
101) type(time_storage_type), pointer :: null_time_storage
102) PetscBool :: found
103)
104) nullify(null_time_storage)
105)
106) dataset_ascii => DatasetAsciiCreate()
107) call DatasetAsciiInit(dataset_ascii)
108) dataset_ascii%array_width = option%nflowdof
109) dataset_ascii%data_type = DATASET_REAL
110) this%dataset => dataset_ascii
111) nullify(dataset_ascii)
112)
113) do
114) call InputReadPflotranString(input,option)
115) if (InputError(input)) exit
116) if (InputCheckExit(input,option)) exit
117)
118) call InputReadWord(input,option,word,PETSC_TRUE)
119) call InputErrorMsg(input,option,'keyword', &
120) 'SOURCE_SINK_SANDBOX,DOWNREG')
121) call StringToUpper(word)
122)
123) call SSSandboxBaseSelectCase(this,input,option,word,found)
124) if (found) cycle
125)
126) select case(trim(word))
127)
128) case('RATE')
129) internal_units = 'unitless/sec'
130) call ConditionReadValues(input,option,word,this%dataset, &
131) units,internal_units)
132) case('POSITIVE_REG_PRESSURE')
133) call InputReadDouble(input,option,this%pressure_max)
134) call InputErrorMsg(input,option,'maximum pressure', &
135) 'SOURCE_SINK_SANDBOX,DOWNREG')
136) call InputReadWord(input,option,word,PETSC_TRUE)
137) if (input%ierr == 0) then
138) internal_units = 'Pa'
139) this%pressure_max = this%pressure_max * &
140) UnitsConvertToInternal(word,internal_units,option)
141) endif
142) case('NEGATIVE_REG_PRESSURE')
143) call InputReadDouble(input,option,this%pressure_min)
144) call InputErrorMsg(input,option,'minimum pressure', &
145) 'SOURCE_SINK_SANDBOX,DOWNREG')
146) call InputReadWord(input,option,word,PETSC_TRUE)
147) if (input%ierr == 0) then
148) internal_units = 'Pa'
149) this%pressure_min = this%pressure_min * &
150) UnitsConvertToInternal(word,internal_units,option)
151) endif
152) case('DELTA_REG_PRESSURE')
153) call InputReadDouble(input,option,this%pressure_delta)
154) call InputErrorMsg(input,option,'delta pressure', &
155) 'SOURCE_SINK_SANDBOX,DOWNREG')
156) call InputReadWord(input,option,word,PETSC_TRUE)
157) if (input%ierr == 0) then
158) internal_units = 'Pa'
159) this%pressure_delta = this%pressure_delta * &
160) UnitsConvertToInternal(word,internal_units,option)
161) endif
162) if (this%pressure_delta <= 1.0d-10) then
163) option%io_buffer = 'SRCSINK_SANDBOX,DOWNREG,DELTA_REG_PRESSURE' // &
164) ': the pressure delta is too close to 0 Pa.'
165) call printErrMsg(option)
166) endif
167) case default
168) call InputKeywordUnrecognized(word,'SRCSINK_SANDBOX,DOWNREG',option)
169) end select
170) enddo
171)
172) if (associated(this%dataset%time_storage)) then
173) ! for now, forcing to step, which makes sense for src/sinks.
174) this%dataset%time_storage%time_interpolation_method = INTERPOLATION_STEP
175) endif
176) call DatasetVerify(this%dataset,null_time_storage,option)
177)
178) end subroutine DownregRead
179)
180) ! ************************************************************************** !
181)
182) subroutine DownregSetup(this,grid,option)
183) !
184) ! Sets up the mass rate src/sink
185) !
186) ! Author: Guoping Tang
187) ! Date: 06/03/14
188)
189) use Option_module
190) use Grid_module
191) use General_Aux_module, only : general_fmw_com => fmw_comp
192)
193) implicit none
194)
195) class(srcsink_sandbox_downreg_type) :: this
196) type(grid_type) :: grid
197) type(option_type) :: option
198)
199) call SSSandboxBaseSetup(this,grid,option)
200)
201) end subroutine DownregSetup
202)
203) ! ************************************************************************** !
204)
205) subroutine DownregUpdate(this,time,option)
206)
207) use Option_module
208) use Dataset_module
209)
210) implicit none
211)
212) class(srcsink_sandbox_downreg_type) :: this
213) PetscReal :: time
214) type(option_type) :: option
215)
216) call DatasetUpdate(this%dataset,time,option)
217)
218) end subroutine DownregUpdate
219)
220) ! ************************************************************************** !
221)
222) subroutine DownregSrcSink(this,Residual,Jacobian,compute_derivative, &
223) material_auxvar,aux_real,option)
224) !
225) ! Evaluates src/sink storing residual and/or Jacobian
226) !
227) ! Author: Guoping Tang
228) ! Date: 06/04/14
229) ! Glenn suggested to use reference pressure for RICHARDS mode
230) ! Gautam suggested to use Nathan's smooth function (6/8/2014)
231)
232) use Option_module
233) use Reaction_Aux_module
234) use Material_Aux_class
235)
236) implicit none
237)
238) class(srcsink_sandbox_downreg_type) :: this
239) type(option_type) :: option
240) PetscBool :: compute_derivative
241) PetscReal :: Residual(option%nflowdof)
242) PetscReal :: Jacobian(option%nflowdof,option%nflowdof)
243) class(material_auxvar_type) :: material_auxvar
244) PetscReal :: aux_real(:)
245) PetscReal :: pressure
246) PetscReal :: pressure_lower, pressure_upper, x
247) PetscReal :: rate_regulator
248) PetscReal :: drate_regulator
249) PetscReal :: temp_real
250) PetscReal :: rate
251)
252) PetscInt :: idof
253)
254) do idof = 1, option%nflowdof
255) if (option%iflowmode == RICHARDS_MODE .and. idof == ONE_INTEGER) then
256) ! regulate liquid pressure in Richards mode
257) pressure = aux_real(3)
258) rate = this%dataset%rarray(1)
259) rate = rate / FMWH2O ! from kg/s to kmol/s (for regression tests)
260) ! rate = rate / aux_real(2) ! from m^3/s to kmol/s (later on, we wll assume m^3/s)
261) if (rate > 0.0d0) then
262) ! source
263) pressure_lower = this%pressure_max - this%pressure_delta/2.0d0 &
264) - option%reference_pressure
265) pressure_upper = pressure_lower + this%pressure_delta
266) if (pressure <= pressure_lower) then
267) rate_regulator = 1.0d0
268) drate_regulator = 0.0d0
269) elseif (pressure >= pressure_upper) then
270) rate_regulator = 0.0d0
271) drate_regulator = 0.0d0
272) else
273) x = (pressure - pressure_lower)/this%pressure_delta
274) rate_regulator = (1.0d0 - x * x) * (1.0d0 - x * x)
275) drate_regulator = 2.0d0 * (1.0d0 - x * x) * (-2.0d0) * x / &
276) this%pressure_delta
277) endif
278) Residual(idof) = rate * rate_regulator
279) else
280) ! sink
281) pressure_lower = this%pressure_min - this%pressure_delta/2.0d0 &
282) - option%reference_pressure
283) pressure_upper = pressure_lower + this%pressure_delta
284) if (pressure <= pressure_lower) then
285) rate_regulator = 0.0d0
286) drate_regulator = 0.0d0
287) elseif (pressure >= pressure_upper) then
288) rate_regulator = 1.0d0
289) drate_regulator = 0.0d0
290) else
291) x = (pressure - pressure_lower)/this%pressure_delta
292) rate_regulator = 1.0d0 - (1.0d0 - x * x) * (1.0d0 - x * x)
293) drate_regulator = (-2.0d0) * (1.0d0 - x * x) * (-2.0d0) * x / &
294) this%pressure_delta
295) endif
296) Residual(idof) = rate * rate_regulator
297) endif
298) else
299) option%io_buffer = 'srcsink_sandbox_downreg is implemented ' // &
300) 'only for RICHARDS mode.'
301) call printErrMsg(option)
302) endif
303) enddo
304)
305) if (compute_derivative) then
306)
307) do idof = 1, option%nflowdof
308) if (option%iflowmode == RICHARDS_MODE .and. idof == ONE_INTEGER) then
309) ! regulate liquid pressure in Richards mode
310) Jacobian(idof,idof) = -1.0d0 * rate * drate_regulator
311) else
312) option%io_buffer = 'srcsink_sandbox_downreg is implemented ' // &
313) 'only for RICHARDS mode.'
314) call printErrMsg(option)
315) endif
316) enddo
317) endif
318)
319) end subroutine DownregSrcSink
320)
321) ! ************************************************************************** !
322)
323) subroutine DownregDestroy(this)
324) !
325) ! Destroys allocatable or pointer objects created in this module
326) !
327) ! Author: Guoping Tang
328) ! Date: 06/04/14
329)
330) use Dataset_module
331) use Dataset_Ascii_class
332)
333) implicit none
334)
335) class(srcsink_sandbox_downreg_type) :: this
336) class(dataset_ascii_type), pointer :: dataset_ascii
337)
338) dataset_ascii => DatasetAsciiCast(this%dataset)
339) call DatasetAsciiDestroy(dataset_ascii)
340)
341) call SSSandboxBaseDestroy(this)
342)
343) end subroutine DownregDestroy
344)
345) end module SrcSink_Sandbox_Downreg_class