saturation_function.F90 coverage: 73.33 %func 42.11 %block
1) module Saturation_Function_module
2)
3)
4) use PFLOTRAN_Constants_module
5)
6) implicit none
7)
8) private
9)
10) #include "petsc/finclude/petscsys.h"
11)
12) type, public :: saturation_function_type
13) PetscInt :: id
14) character(len=MAXWORDLENGTH) :: name
15) character(len=MAXWORDLENGTH) :: saturation_function_ctype
16) PetscInt :: saturation_function_itype
17) character(len=MAXWORDLENGTH) :: permeability_function_ctype
18) PetscInt :: permeability_function_itype
19) PetscBool :: print_me
20) PetscReal, pointer :: Sr(:)
21) PetscReal, pointer :: Kr0(:)
22) PetscReal :: m
23) PetscReal :: lambda
24) PetscReal :: alpha
25) PetscReal :: pcwmax
26) PetscReal :: betac
27) PetscReal :: power
28) PetscInt :: hysteresis_id
29) PetscInt :: hysteresis_params(6)
30) PetscReal :: sat_spline_low
31) PetscReal :: sat_spline_high
32) PetscReal :: sat_spline_coefficients(3)
33) PetscReal :: pres_spline_low
34) PetscReal :: pres_spline_high
35) PetscReal :: pres_spline_coefficients(4)
36) PetscReal :: ani_A ! parameters for anisotropic relative permeability model
37) PetscReal :: ani_B ! parameters for anisotropic relative permeability model
38) PetscReal :: ani_C ! parameters for anisotropic relative permeability model
39)
40) type(saturation_function_type), pointer :: next
41) end type saturation_function_type
42)
43) type, public :: saturation_function_ptr_type
44) type(saturation_function_type), pointer :: ptr
45) end type saturation_function_ptr_type
46)
47) interface SaturationFunctionCompute
48) module procedure SaturationFunctionCompute1
49) module procedure SaturationFunctionCompute2
50) module procedure SaturationFunctionCompute3
51) end interface
52)
53) public :: SaturationFunctionCreate, &
54) SaturationFunctionDestroy, &
55) SaturationFunctionAddToList, &
56) SaturationFunctionCompute, &
57) SaturatFuncConvertListToArray, &
58) SatFunctionComputePolynomial, &
59) PermFunctionComputePolynomial, &
60) SaturationFunctionRead, &
61) SatFuncGetLiqRelPermFromSat, &
62) SatFuncGetGasRelPermFromSat, &
63) SatFuncGetCapillaryPressure, &
64) SaturationFunctionGetID, &
65) SatFuncComputeIcePExplicit, &
66) CapillaryPressureThreshold, &
67) SatFuncComputeIcePKImplicit, &
68) SatFuncComputeIcePKExplicit, &
69) SatFuncComputeIceDallAmico, &
70) SatFuncComputeIcePKExplicitNoCryo
71)
72) ! Saturation function
73) PetscInt, parameter, public :: VAN_GENUCHTEN = 1
74) PetscInt, parameter, public :: BROOKS_COREY = 2
75) PetscInt, parameter, public :: THOMEER_COREY = 3
76) PetscInt, parameter, public :: NMT_EXP = 4
77) PetscInt, parameter, public :: PRUESS_1 = 5
78) PetscInt, parameter, public :: LINEAR_MODEL = 6
79) PetscInt, parameter, public :: VAN_GENUCHTEN_PARKER = 7
80) PetscInt, parameter, public :: VAN_GENUCHTEN_DOUGHTY = 8
81) PetscInt, parameter, public :: LEVERETT = 9
82)
83) ! Permeability function
84) PetscInt, parameter :: DEFAULT = 0
85) PetscInt, parameter, public :: BURDINE = 1
86) PetscInt, parameter, public :: MUALEM = 2
87) PetscInt, parameter, public :: FATT_KLIKOFF = 3
88)
89) contains
90)
91) ! ************************************************************************** !
92)
93) function SaturationFunctionCreate(option)
94) !
95) ! Creates a saturation function
96) !
97) ! Author: Glenn Hammond
98) ! Date: 11/02/07
99) !
100)
101) use Option_module
102)
103) implicit none
104)
105) type(saturation_function_type), pointer :: SaturationFunctionCreate
106) type(option_type) :: option
107)
108) type(saturation_function_type), pointer :: saturation_function
109)
110) allocate(saturation_function)
111) saturation_function%id = 0
112) saturation_function%name = ''
113) saturation_function%saturation_function_ctype = 'VAN_GENUCHTEN'
114) saturation_function%saturation_function_itype = VAN_GENUCHTEN
115) saturation_function%permeability_function_ctype = 'MUALEM'
116) saturation_function%permeability_function_itype = MUALEM
117) saturation_function%print_me = PETSC_FALSE
118) allocate(saturation_function%Sr(option%nphase))
119) saturation_function%Sr = 0.d0
120) allocate(saturation_function%Kr0(option%nphase))
121) saturation_function%Kr0 = 1.0d0
122) saturation_function%m = 0.d0
123) saturation_function%lambda = 0.d0
124) saturation_function%alpha = 0.d0
125) saturation_function%pcwmax = 1.d9
126) saturation_function%betac = 0.d0
127) saturation_function%power = 0.d0
128) saturation_function%hysteresis_id = 0
129) saturation_function%hysteresis_params = 0
130) saturation_function%sat_spline_low = 0.d0
131) saturation_function%sat_spline_high = 0.d0
132) saturation_function%sat_spline_coefficients = 0.d0
133) saturation_function%pres_spline_low = 0.d0
134) saturation_function%pres_spline_high = 0.d0
135) saturation_function%pres_spline_coefficients = 0.d0
136) saturation_function%ani_A = 0.d0
137) saturation_function%ani_B = 0.d0
138) saturation_function%ani_C = 0.d0
139) nullify(saturation_function%next)
140) SaturationFunctionCreate => saturation_function
141)
142) end function SaturationFunctionCreate
143)
144) ! ************************************************************************** !
145)
146) subroutine SaturationFunctionRead(saturation_function,input,option)
147) !
148) ! Reads in contents of a saturation_function card
149) !
150) ! Author: Glenn Hammond
151) ! Date: 01/21/09
152) !
153)
154) use Option_module
155) use Input_Aux_module
156) use String_module
157)
158) implicit none
159)
160) type(saturation_function_type) :: saturation_function
161) type(input_type), pointer :: input
162) type(option_type) :: option
163) PetscInt :: iphase
164)
165) character(len=MAXWORDLENGTH) :: keyword, word
166) character(len=MAXSTRINGLENGTH) :: string
167) PetscReal :: tempreal
168)
169) select case(option%iflowmode)
170) case(G_MODE)
171) option%io_buffer = 'SATURATION_FUNCTION card is no longer ' // &
172) 'supported for GENERAL mode. Please use CHARACTERISTIC_' // &
173) 'CURVES card defined on the PFLOTRAN wiki.'
174) call printErrMsg(option)
175) end select
176)
177) input%ierr = 0
178) do
179)
180) call InputReadPflotranString(input,option)
181)
182) if (InputCheckExit(input,option)) exit
183)
184) call InputReadWord(input,option,keyword,PETSC_TRUE)
185) call InputErrorMsg(input,option,'keyword','SATURATION_FUNCTION')
186) call StringToUpper(keyword)
187)
188) select case(trim(keyword))
189)
190) case('PERMEABILITY_FUNCTION_TYPE')
191) call InputReadWord(input,option, &
192) saturation_function%permeability_function_ctype, &
193) PETSC_TRUE)
194) call InputErrorMsg(input,option,'permeability function type', &
195) 'SATURATION_FUNCTION')
196) case('SATURATION_FUNCTION_TYPE')
197) call InputReadWord(input,option, &
198) saturation_function%saturation_function_ctype, &
199) PETSC_TRUE)
200) call InputErrorMsg(input,option,'saturation function type', &
201) 'SATURATION_FUNCTION')
202) case('PERMEABILITY_END_POINT')
203) select case(option%iflowmode)
204) case(FLASH2_MODE)
205) call InputReadWord(input,option,keyword,PETSC_TRUE)
206) call InputErrorMsg(input,option,'keyword','PERMEABILITY_FUNCTION')
207) call StringToUpper(keyword)
208) select case(trim(keyword))
209) case('WATER','WATER_PHASE','LIQUID','LIQUID_PHASE')
210) iphase = 1
211) case('CO2','CO2_PHASE','GAS','GAS_PHASE')
212) iphase = 2
213) end select
214) call InputReadDouble(input,option,saturation_function%Kr0(iphase))
215) word = trim(keyword) // 'permeabiliy end point'
216) call InputErrorMsg(input,option,word,'PERMEABILITY_FUNCTION')
217) case(MPH_MODE)
218) call InputReadWord(input,option,keyword,PETSC_TRUE)
219) call InputErrorMsg(input,option,'keyword','PERMEABILITY_FUNCTION')
220) call StringToUpper(keyword)
221) select case(trim(keyword))
222) case('WATER','WATER_PHASE','LIQUID','LIQUID_PHASE')
223) iphase = 1
224) case('CO2','CO2_PHASE','GAS','GAS_PHASE')
225) iphase = 2
226) end select
227) call InputReadDouble(input,option,saturation_function%Kr0(iphase))
228) word = trim(keyword) // 'permeabiliy end point'
229) call InputErrorMsg(input,option,word,'PERMEABILITY_FUNCTION')
230) case(IMS_MODE)
231) call InputReadWord(input,option,keyword,PETSC_TRUE)
232) call InputErrorMsg(input,option,'keyword','PERMEABILITY_FUNCTION')
233) call StringToUpper(keyword)
234) select case(trim(keyword))
235) case('WATER','WATER_PHASE','LIQUID','LIQUID_PHASE')
236) iphase = 1
237) case('CO2','CO2_PHASE','GAS','GAS_PHASE')
238) iphase = 2
239) end select
240) call InputReadDouble(input,option,saturation_function%Kr0(iphase))
241) word = trim(keyword) // 'permeabiliy end point'
242) call InputErrorMsg(input,option,word,'PERMEABILITY_FUNCTION')
243) end select
244) case('RESIDUAL_SATURATION')
245) select case(option%iflowmode)
246) case(FLASH2_MODE)
247) call InputReadWord(input,option,keyword,PETSC_TRUE)
248) call InputErrorMsg(input,option,'keyword','SATURATION_FUNCTION')
249) call StringToUpper(keyword)
250) select case(trim(keyword))
251) case('WATER','WATER_PHASE','LIQUID','LIQUID_PHASE')
252) iphase = 1
253) case('CO2','CO2_PHASE','GAS','GAS_PHASE')
254) iphase = 2
255) end select
256) call InputReadDouble(input,option,saturation_function%Sr(iphase))
257) word = trim(keyword) // ' residual saturation'
258) call InputErrorMsg(input,option,word,'SATURATION_FUNCTION')
259) case(MPH_MODE)
260) call InputReadWord(input,option,keyword,PETSC_TRUE)
261) call InputErrorMsg(input,option,'keyword','SATURATION_FUNCTION')
262) call StringToUpper(keyword)
263) select case(trim(keyword))
264) case('WATER','WATER_PHASE','LIQUID','LIQUID_PHASE')
265) iphase = 1
266) case('CO2','CO2_PHASE','GAS','GAS_PHASE')
267) iphase = 2
268) end select
269) call InputReadDouble(input,option,saturation_function%Sr(iphase))
270) word = trim(keyword) // ' residual saturation'
271) call InputErrorMsg(input,option,word,'SATURATION_FUNCTION')
272) case(IMS_MODE)
273) call InputReadWord(input,option,keyword,PETSC_TRUE)
274) call InputErrorMsg(input,option,'keyword','SATURATION_FUNCTION')
275) call StringToUpper(keyword)
276) select case(trim(keyword))
277) case('WATER','WATER_PHASE','LIQUID','LIQUID_PHASE')
278) iphase = 1
279) case('CO2','CO2_PHASE','GAS','GAS_PHASE')
280) iphase = 2
281) case('OIL','OIL_PHASE','NAPL','NAPL_PHASE')
282) iphase = 3
283) end select
284) call InputReadDouble(input,option,saturation_function%Sr(iphase))
285) word = trim(keyword) // ' residual saturation'
286) call InputErrorMsg(input,option,word,'SATURATION_FUNCTION')
287) case(G_MODE)
288) iphase = 0
289) string = input%buf
290) call InputReadDouble(input,option,tempreal)
291) ! call InputErrorMsg(input,option,'residual saturation','SATURATION_FUNCTION')
292) if (input%ierr /= 0) then
293) input%ierr = 0
294) input%buf = string
295) call InputReadWord(input,option,keyword,PETSC_TRUE)
296) call InputErrorMsg(input,option,'phase', &
297) 'SATURATION_FUNCTION,RESIDUAL_SATURATION')
298) call StringToUpper(keyword)
299) select case(trim(keyword))
300) case('LIQUID','LIQUID_PHASE')
301) iphase = 1
302) case('GAS','GAS_PHASE')
303) iphase = 2
304) case default
305) call InputKeywordUnrecognized(keyword, &
306) 'SATURATION_FUNCTION,RESIDUAL_SATURATION',option)
307) end select
308) call InputReadDouble(input,option,tempreal)
309) word = trim(keyword) // ' residual saturation'
310) call InputErrorMsg(input,option,word,'SATURATION_FUNCTION')
311) else
312) ! if missing phase keyword, assume for all phases and set
313) ! buffer to value
314) endif
315) if (iphase > 0) then
316) saturation_function%Sr(iphase) = tempreal
317) else
318) saturation_function%Sr(:) = tempreal
319) endif
320) case(RICHARDS_MODE,TH_MODE)
321) call InputReadDouble(input,option,saturation_function%Sr(1))
322) call InputErrorMsg(input,option,'residual saturation','SATURATION_FUNCTION')
323) end select
324) case('LAMBDA')
325) call InputReadDouble(input,option,saturation_function%lambda)
326) call InputErrorMsg(input,option,'lambda','SATURATION_FUNCTION')
327) saturation_function%m = saturation_function%lambda
328) case('ALPHA')
329) call InputReadDouble(input,option,saturation_function%alpha)
330) call InputErrorMsg(input,option,'alpha','SATURATION_FUNCTION')
331) case('MAX_CAPILLARY_PRESSURE')
332) call InputReadDouble(input,option,saturation_function%pcwmax)
333) call InputErrorMsg(input,option,'maximum capillary pressure','SATURATION_FUNCTION')
334) case('BETAC')
335) call InputReadDouble(input,option,saturation_function%betac)
336) call InputErrorMsg(input,option,'betac','SATURATION_FUNCTION')
337) case('POWER')
338) call InputReadDouble(input,option,saturation_function%power)
339) call InputErrorMsg(input,option,'power','SATURATION_FUNCTION')
340) case('ANI_A')
341) call InputReadDouble(input,option,saturation_function%ani_A)
342) call InputErrorMsg(input,option,'ani_A','SATURATION_FUNCTION')
343) case('ANI_B')
344) call InputReadDouble(input,option,saturation_function%ani_B)
345) call InputErrorMsg(input,option,'ani_B','SATURATION_FUNCTION')
346) case('ANI_C')
347) call InputReadDouble(input,option,saturation_function%ani_C)
348) call InputErrorMsg(input,option,'ani_C','SATURATION_FUNCTION')
349) case('VERIFY')
350) saturation_function%print_me = PETSC_TRUE
351) case default
352) call InputKeywordUnrecognized(keyword,'SATURATION_FUNCTION',option)
353) end select
354)
355) enddo
356)
357) if (saturation_function%m < 1.d-40 .and. &
358) .not. &
359) (StringCompare(saturation_function%name,'default',SEVEN_INTEGER) .or. &
360) StringCompareIgnoreCase(saturation_function%saturation_function_ctype, &
361) 'leverett',EIGHT_INTEGER))) then
362) option%io_buffer = 'Saturation function parameter "m" not set ' // &
363) 'properly in saturation function "' // &
364) trim(saturation_function%name) // '".'
365) call printErrMsg(option)
366) endif
367)
368) call SaturationFunctionSetTypes(saturation_function,option)
369)
370) end subroutine SaturationFunctionRead
371)
372) ! ************************************************************************** !
373)
374) subroutine SaturationFunctionSetTypes(saturation_function,option)
375) !
376) ! Sets the integer values that map the saturation and permeability
377) ! functions to a specific type
378) !
379) ! Author: Glenn Hammond
380) ! Date: 04/23/14
381)
382) use Option_module
383) use String_module
384)
385) implicit none
386)
387) type(saturation_function_type) :: saturation_function
388) type(option_type) :: option
389)
390) ! set permeability function integer type
391) call StringToUpper(saturation_function%permeability_function_ctype)
392) select case(trim(saturation_function%permeability_function_ctype))
393) case('DEFAULT')
394) saturation_function%permeability_function_itype = DEFAULT
395) case('BURDINE')
396) saturation_function%permeability_function_itype = BURDINE
397) case('MUALEM')
398) saturation_function%permeability_function_itype = MUALEM
399) case('NMT_EXP')
400) saturation_function%permeability_function_itype = NMT_EXP
401) case('PRUESS_1')
402) saturation_function%permeability_function_itype = PRUESS_1
403) case('VAN_GENUCHTEN_PARKER')
404) saturation_function%permeability_function_itype = VAN_GENUCHTEN_PARKER
405) case('FATT_KLIKOFF')
406) saturation_function%permeability_function_itype = FATT_KLIKOFF
407) case default
408) option%io_buffer = 'Permeability function type "' // &
409) trim(saturation_function%permeability_function_ctype) // &
410) '" not recognized ' // &
411) ' in saturation function ' // &
412) trim(saturation_function%name)
413) call printErrMsg(option)
414) end select
415)
416) ! set saturation function integer type
417) call StringToUpper(saturation_function%saturation_function_ctype)
418) select case(trim(saturation_function%saturation_function_ctype))
419) case('VAN_GENUCHTEN')
420) saturation_function%saturation_function_itype = VAN_GENUCHTEN
421) case('BROOKS_COREY')
422) saturation_function%saturation_function_itype = BROOKS_COREY
423) case('LINEAR_MODEL')
424) saturation_function%saturation_function_itype = LINEAR_MODEL
425) case('THOMEER_COREY')
426) saturation_function%saturation_function_itype = THOMEER_COREY
427) case('NMT_EXP')
428) saturation_function%saturation_function_itype = NMT_EXP
429) case('PRUESS_1')
430) saturation_function%saturation_function_itype = PRUESS_1
431) case('VAN_GENUCHTEN_PARKER')
432) saturation_function%saturation_function_itype = VAN_GENUCHTEN_PARKER
433) case('VAN_GENUCHTEN_DOUGHTY')
434) saturation_function%saturation_function_itype = VAN_GENUCHTEN_DOUGHTY
435) case('LEVERETT')
436) saturation_function%saturation_function_itype = LEVERETT
437) case default
438) option%io_buffer = 'Saturation function type "' // &
439) trim(saturation_function%saturation_function_ctype) // &
440) '" not recognized ' // &
441) ' in saturation function ' // &
442) trim(saturation_function%name)
443) call printErrMsg(option)
444) end select
445)
446) end subroutine SaturationFunctionSetTypes
447)
448) ! ************************************************************************** !
449)
450) subroutine SatFunctionComputePolynomial(option,saturation_function)
451) !
452) ! Computes a spline spanning the
453) ! discontinuity in Brooks Corey
454) !
455) ! Author: Glenn Hammond
456) ! Date: 02/27/08
457) !
458)
459) use Option_module
460) use Utility_module
461)
462) implicit none
463)
464) type(option_type) :: option
465) type(saturation_function_type) :: saturation_function
466)
467) PetscReal :: b(4)
468) PetscReal :: pressure_1, pressure_2
469) PetscReal :: saturation_1, saturation_2
470)
471) PetscReal :: n
472)
473) select case(saturation_function%saturation_function_itype)
474) case(BROOKS_COREY)
475)
476) ! polynomial fitting pc as a function of saturation
477) ! 1.05 is essentially pc*alpha (i.e. pc = 1.05/alpha)
478) saturation_1 = 1.05d0**(-saturation_function%lambda)
479) saturation_2 = 1.d0
480)
481) saturation_function%sat_spline_low = saturation_1
482) saturation_function%sat_spline_high = saturation_2
483)
484) b = 0.d0
485) ! fill right hand side
486) ! capillary pressure at 1
487) b(1) = 1.05d0/saturation_function%alpha
488) ! capillary pressure at 2
489) b(2) = 0.d0
490) ! derivative of pressure at saturation_1
491) ! pc = Se**(-1/lambda)/alpha
492) ! dpc_dSe = -1/lambda*Se**(-1/lambda-1)/alpha
493) b(3) = -1.d0/saturation_function%lambda* &
494) saturation_1**(-1.d0/saturation_function%lambda-1.d0)/ &
495) saturation_function%alpha
496)
497) call QuadraticPolynomialSetup(saturation_1,saturation_2,b(1:3), &
498) ! indicates derivative given at 1
499) PETSC_TRUE)
500)
501) saturation_function%sat_spline_coefficients(1:3) = b(1:3)
502)
503) ! polynomial fitting saturation as a function of pc
504) !geh: cannot invert the pressure/saturation relationship above
505) ! since it can result in saturations > 1 with both
506) ! quadratic and cubic polynomials
507) ! fill matix with values
508) pressure_1 = 0.95/saturation_function%alpha
509) pressure_2 = 1.05/saturation_function%alpha
510)
511) saturation_function%pres_spline_low = pressure_1
512) saturation_function%pres_spline_high = pressure_2
513)
514) b = 0.d0
515) ! Se at 1
516) b(1) = 1.d0
517) ! Se at 2
518) b(2) = (pressure_2*saturation_function%alpha)** &
519) (-saturation_function%lambda)
520) ! derivative of Se at 1
521) b(3) = 0.d0
522) ! derivative of Se at 2
523) b(4) = -saturation_function%lambda/pressure_2* &
524) (pressure_2*saturation_function%alpha)** &
525) (-saturation_function%lambda)
526)
527) call CubicPolynomialSetup(pressure_1,pressure_2,b)
528)
529) saturation_function%pres_spline_coefficients(1:4) = b(1:4)
530)
531) case(VAN_GENUCHTEN)
532)
533) ! return for now
534) return
535)
536) !geh: keep for now
537) #if 0
538) ! geh: needs to be refactored to consider capillary pressure
539) ! fill matix with values
540) ! these are capillary pressures
541) pressure_low = 0 ! saturated
542) pressure_high = 0.01d0*option%reference_pressure ! just below saturated
543)
544) saturation_function%spline_low = pressure_low
545) saturation_function%spline_high = pressure_high
546)
547) n = 1.d0/(1.d0 - saturation_function%m)
548) b(1) = (1.d0+(pressure_high*saturation_function%alpha)**n)** &
549) (-saturation_function%m)
550) b(2) = 1.d0
551) b(3) = -saturation_function%m*n*saturation_function%alpha* &
552) (saturation_function%alpha*pressure_high)**(n-1.d0)* &
553) (1.d0+(saturation_function%alpha*pressure_high)**n)** &
554) (-saturation_function%m-1.d0)
555) b(4) = 0.d0
556)
557) call CubicPolynomialSetup(pressure_high,pressure_low,b)
558)
559) saturation_function%spline_coefficients(1:4) = b(1:4)
560) #endif
561)
562) end select
563)
564) end subroutine SatFunctionComputePolynomial
565)
566) ! ************************************************************************** !
567)
568) subroutine PermFunctionComputePolynomial(option,saturation_function)
569) !
570) ! Computes a spline spanning the
571) ! discontinuity in Brooks Corey
572) !
573) ! Author: Glenn Hammond
574) ! Date: 02/27/12
575) !
576)
577) use Option_module
578) use Utility_module
579)
580) implicit none
581)
582) type(option_type) :: option
583) type(saturation_function_type) :: saturation_function
584)
585) PetscReal :: b(4)
586) PetscReal :: Se_high, Se_low, one_over_m, Se_one_over_m, m
587)
588) select case(saturation_function%saturation_function_itype)
589)
590) case(BROOKS_COREY)
591)
592) case(VAN_GENUCHTEN)
593)
594) #ifdef MUALEM_SPLINE
595) ! fill matix with values
596) Se_low = 0.99d0 ! saturated
597) Se_high = 1.d0 ! just below saturated
598)
599) saturation_function%spline_low = Se_low
600) saturation_function%spline_high = Se_high
601)
602) m = saturation_function%m
603) one_over_m = 1.d0/m
604) Se_one_over_m = Se_low**one_over_m
605) b(1) = 1.d0
606) b(2) = sqrt(Se_low)*(1.d0-(1.d0-Se_one_over_m)**m)**2.d0
607) b(3) = 0.d0
608) b(4) = 0.5d0*b(2)/Se_low+ &
609) 2.d0*Se_low**(one_over_m-0.5d0)* &
610) (1.d0-Se_one_over_m)**(m-1.d0)* &
611) (1.d0-(1.d0-Se_one_over_m)**m)
612)
613) call CubicPolynomialSetup(Se_high,Se_low,b)
614)
615) saturation_function%spline_coefficients(1:4) = b(1:4)
616) #endif
617)
618) end select
619)
620) end subroutine PermFunctionComputePolynomial
621)
622) ! ************************************************************************** !
623)
624) subroutine SaturationFunctionAddToList(saturation_function,list)
625) !
626) ! Adds a saturation function to linked list
627) !
628) ! Author: Glenn Hammond
629) ! Date: 11/02/07
630) !
631)
632) implicit none
633)
634) type(saturation_function_type), pointer :: saturation_function
635) type(saturation_function_type), pointer :: list
636)
637) type(saturation_function_type), pointer :: cur_saturation_function
638)
639) if (associated(list)) then
640) cur_saturation_function => list
641) ! loop to end of list
642) do
643) if (.not.associated(cur_saturation_function%next)) exit
644) cur_saturation_function => cur_saturation_function%next
645) enddo
646) cur_saturation_function%next => saturation_function
647) else
648) list => saturation_function
649) endif
650)
651) end subroutine SaturationFunctionAddToList
652)
653) ! ************************************************************************** !
654)
655) subroutine SaturatFuncConvertListToArray(list,array,option)
656) !
657) ! Creates an array of pointers to the
658) ! saturation functions in the list
659) !
660) ! Author: Glenn Hammond
661) ! Date: 12/11/07
662) !
663)
664) use String_module
665) use Option_module
666)
667) implicit none
668)
669) type(saturation_function_type), pointer :: list
670) type(saturation_function_ptr_type), pointer :: array(:)
671) type(option_type) :: option
672)
673) type(saturation_function_type), pointer :: cur_saturation_function
674) PetscInt :: count
675)
676) count = 0
677) cur_saturation_function => list
678) do
679) if (.not.associated(cur_saturation_function)) exit
680) count = count + 1
681) cur_saturation_function => cur_saturation_function%next
682) enddo
683)
684) if (associated(array)) deallocate(array)
685) allocate(array(count))
686)
687) count = 0
688) cur_saturation_function => list
689) do
690) if (.not.associated(cur_saturation_function)) exit
691) count = count + 1
692) cur_saturation_function%id = count
693) array(count)%ptr => cur_saturation_function
694) if (cur_saturation_function%print_me .and. &
695) option%myrank == option%io_rank) then
696) call SaturationFunctionVerify(cur_saturation_function,option)
697) endif
698) cur_saturation_function => cur_saturation_function%next
699) enddo
700)
701) end subroutine SaturatFuncConvertListToArray
702)
703) ! ************************************************************************** !
704)
705) subroutine SaturationFunctionCompute1(capillary_pressure,saturation, &
706) relative_perm, &
707) dsat_dpres,dkr_dpres, &
708) saturation_function, &
709) auxvar1,auxvar2, &
710) option)
711) !
712) ! Computes the saturation and relative permeability
713) ! (and associated derivatives) as a function of
714) ! capillary pressure
715) !
716) ! Author: Glenn Hammond
717) ! Date: 2/9/12
718) !
719) use Option_module
720)
721) implicit none
722)
723) PetscReal :: capillary_pressure, saturation, relative_perm
724) PetscReal :: dsat_dpres, dkr_dpres
725) type(saturation_function_type) :: saturation_function
726) PetscReal :: auxvar1,auxvar2
727) type(option_type) :: option
728)
729) PetscBool :: switch_to_saturated
730)
731) call SaturationFunctionCompute2(capillary_pressure,saturation, &
732) relative_perm, &
733) dsat_dpres,dkr_dpres, &
734) saturation_function, &
735) auxvar1,auxvar2, &
736) switch_to_saturated,option)
737)
738) end subroutine SaturationFunctionCompute1
739)
740) ! ************************************************************************** !
741)
742) subroutine SaturationFunctionCompute2(capillary_pressure,saturation, &
743) relative_perm, &
744) dsat_dpres,dkr_dpres, &
745) saturation_function, &
746) auxvar1,auxvar2, &
747) switch_to_saturated,option)
748) !
749) ! Computes the saturation and relative permeability
750) ! (and associated derivatives) as a function of
751) ! capillary pressure
752) !
753) ! (1) Chen, J., J.W. Hopmans, M.E. Grismer (1999) "Parameter estimation of
754) ! of two-fluid capillary pressure-saturation and permeability functions",
755) ! Advances in Water Resources, Vol. 22, No. 5, pp 479-493,
756) ! http://dx.doi.org/10.1016/S0309-1708(98)00025-6.
757) !
758) ! Author: Glenn Hammond
759) ! Date: 12/11/07
760) !
761) use Option_module
762) use Utility_module, only:CubicPolynomialEvaluate
763)
764) implicit none
765)
766) PetscReal :: capillary_pressure, saturation, relative_perm
767) PetscReal :: dsat_dpres, dkr_dpres
768) type(saturation_function_type) :: saturation_function
769) PetscReal :: auxvar1,auxvar2
770) PetscBool :: switch_to_saturated
771) type(option_type) :: option
772)
773) PetscInt :: iphase
774) PetscReal :: alpha, lambda, m, n, Sr, one_over_alpha, pcmax
775) PetscReal :: pc, Se, one_over_m, Se_one_over_m, dSe_dpc, dsat_dpc, dkr_dpc
776) PetscReal :: dkr_dSe, power
777) PetscReal :: pc_alpha, pc_alpha_n, one_plus_pc_alpha_n
778) PetscReal :: pc_alpha_neg_lambda
779) PetscReal :: por, perm
780) PetscReal :: Fg, a, Pd, PHg, pct_over_pcmax, pc_over_pcmax, pc_log_ratio
781)
782) PetscReal, parameter :: pc_alpha_n_epsilon = 1.d-15
783)
784) iphase = 1
785) dsat_dpres = 0.d0
786) dkr_dpres = 0.d0
787)
788) ! compute saturation
789) select case(saturation_function%saturation_function_itype)
790) case(VAN_GENUCHTEN)
791) ! reference #1
792) #if 0
793) if (pc < saturation_function%spline_low) then
794) saturation = 1.d0
795) relative_perm = 1.d0
796) switch_to_saturated = PETSC_TRUE
797) return
798) else if (pc < saturation_function%spline_high) then
799) Sr = saturation_function%Sr(iphase)
800) call CubicPolynomialEvaluate(saturation_function%spline_coefficients, &
801) pc,Se,dSe_dpc)
802) saturation = Sr + (1.d0-Sr)*Se
803) dsat_dpc = (1.d0-Sr)*dSe_dpc
804) #else
805) if (capillary_pressure <= 0.d0) then
806) saturation = 1.d0
807) relative_perm = 1.d0
808) switch_to_saturated = PETSC_TRUE
809) return
810) #endif
811) else
812) alpha = saturation_function%alpha
813) pc = capillary_pressure
814) m = saturation_function%m
815) n = 1.d0/(1.d0-m)
816) pc_alpha = pc*alpha
817) pc_alpha_n = pc_alpha**n
818) !geh: This conditional does not catch potential cancelation in
819) ! the dkr_dSe deriviative calculation. Therefore, I am setting
820) ! an epsilon here
821) ! if (1.d0 + pc_alpha_n == 1.d0) then ! check for zero perturbation
822) if (pc_alpha_n < pc_alpha_n_epsilon) then
823) saturation = 1.d0
824) relative_perm = 1.d0
825) switch_to_saturated = PETSC_TRUE
826) return
827) endif
828) one_plus_pc_alpha_n = 1.d0+pc_alpha_n
829) Sr = saturation_function%Sr(iphase)
830) Se = one_plus_pc_alpha_n**(-m)
831) dSe_dpc = -m*n*alpha*pc_alpha_n/(pc_alpha*one_plus_pc_alpha_n**(m+1))
832) saturation = Sr + (1.d0-Sr)*Se
833) dsat_dpc = (1.d0-Sr)*dSe_dpc
834) endif
835) if (saturation > 1.d0) then
836) print *, option%myrank, 'vG Saturation > 1:', saturation
837) else if (saturation > 1.d0 .or. saturation < Sr) then
838) print *, option%myrank, 'vG Saturation < Sr:', saturation, Sr
839) endif
840) ! compute relative permeability
841) select case(saturation_function%permeability_function_itype)
842) case(BURDINE)
843) ! reference #1
844) one_over_m = 1.d0/m
845) Se_one_over_m = Se**one_over_m
846) relative_perm = Se*Se*(1.d0-(1.d0-Se_one_over_m)**m)
847) dkr_dSe = 2.d0*relative_perm/Se + &
848) Se*Se_one_over_m*(1.d0-Se_one_over_m)**(m-1.d0)
849) dkr_dpc = dkr_dSe*dSe_dpc
850) case(MUALEM)
851) ! reference #1
852) #ifdef MUALEM_SPLINE
853) if (Se > saturation_function%spline_low) then
854) call CubicPolynomialEvaluate( &
855) saturation_function%spline_coefficients, &
856) Se,relative_perm,dkr_dSe)
857) else
858) #endif
859) one_over_m = 1.d0/m
860) Se_one_over_m = Se**one_over_m
861) relative_perm = sqrt(Se)*(1.d0-(1.d0-Se_one_over_m)**m)**2.d0
862) dkr_dSe = 0.5d0*relative_perm/Se+ &
863) 2.d0*Se**(one_over_m-0.5d0)* &
864) (1.d0-Se_one_over_m)**(m-1.d0)* &
865) (1.d0-(1.d0-Se_one_over_m)**m)
866) #ifdef MUALEM_SPLINE
867) endif
868) #endif
869) dkr_dpc = dkr_dSe*dSe_dpc
870) case default
871) option%io_buffer = 'Unknown relative permeabilty function'
872) call printErrMsg(option)
873) end select
874) case(BROOKS_COREY)
875) ! reference #1
876) alpha = saturation_function%alpha
877) one_over_alpha = 1.d0/alpha
878) pc = capillary_pressure
879) if (pc < saturation_function%pres_spline_low) then
880) saturation = 1.d0
881) relative_perm = 1.d0
882) switch_to_saturated = PETSC_TRUE
883) return
884) else if (pc < saturation_function%pres_spline_high) then
885) lambda = saturation_function%lambda
886) Sr = saturation_function%Sr(iphase)
887) call CubicPolynomialEvaluate(saturation_function% &
888) pres_spline_coefficients, &
889) pc,Se,dSe_dpc)
890) saturation = Sr + (1.d0-Sr)*Se
891) dsat_dpc = (1.d0-Sr)*dSe_dpc
892) else
893) lambda = saturation_function%lambda
894) Sr = saturation_function%Sr(iphase)
895) pc_alpha_neg_lambda = (pc*alpha)**(-lambda)
896) Se = pc_alpha_neg_lambda
897) dSe_dpc = -lambda/pc*pc_alpha_neg_lambda
898) saturation = Sr + (1.d0-Sr)*Se
899) ! dsat_dpc = -lambda*(1.d0-Sr)/pc*pc_alpha_neg_lambda
900) dsat_dpc = (1.d0-Sr)*dSe_dpc
901) endif
902) if (saturation > 1.d0) then
903) print *, option%myrank, 'BC Saturation > 1:', saturation
904) else if (saturation > 1.d0 .or. saturation < Sr) then
905) print *, option%myrank, 'BC Saturation < Sr:', saturation, Sr
906) endif
907) ! compute relative permeability
908) select case(saturation_function%permeability_function_itype)
909) case(BURDINE)
910) ! reference #1
911) power = 3.d0+2.d0/lambda
912) relative_perm = Se**power
913) dkr_dSe = power*relative_perm/Se
914) dkr_dpc = dkr_dSe*dSe_dpc
915) case(MUALEM)
916) ! reference #1
917) power = 2.5d0+2.d0/lambda
918) relative_perm = Se**power
919) dkr_dSe = power*relative_perm/Se
920) dkr_dpc = dkr_dSe*dSe_dpc
921) case default
922) option%io_buffer = 'Unknown relative permeabilty function'
923) call printErrMsg(option)
924) end select
925) case(LINEAR_MODEL)
926) ! Added by Bwalya Malama 01/30/2014
927) alpha = saturation_function%alpha
928) one_over_alpha = 1.d0/alpha
929) lambda = saturation_function%lambda
930) pcmax = saturation_function%pcwmax
931) pc = capillary_pressure
932)
933) if (capillary_pressure <= 0.d0) then
934) saturation = 1.d0
935) relative_perm = 1.d0
936) switch_to_saturated = PETSC_TRUE
937) return
938) else
939) Sr = saturation_function%Sr(iphase)
940) Se = (pcmax-pc)/(pcmax-one_over_alpha)
941) dSe_dpc = -1.d0/(pcmax-one_over_alpha)
942) saturation = Sr + (1.d0-Sr)*Se
943) dsat_dpc = (1.d0-Sr)*dSe_dpc
944) endif
945) if (saturation > 1.d0) then
946) print *, option%myrank, 'BC Saturation > 1:', saturation
947) else if (saturation < Sr) then
948) print *, option%myrank, 'BC Saturation < Sr:', saturation, Sr
949) endif
950) select case(saturation_function%permeability_function_itype)
951) case(BURDINE)
952) relative_perm = Se
953) dkr_dSe = 1.d0
954) case(MUALEM)
955) power = 5.d-1
956) pct_over_pcmax = one_over_alpha/pcmax
957) pc_over_pcmax = 1.d0-(1.d0-pct_over_pcmax)*Se
958) pc_log_ratio = log(pc_over_pcmax)/log(pct_over_pcmax)
959) relative_perm = (Se**power)*(pc_log_ratio**2.d0)
960) case default
961) option%io_buffer = 'Unknown relative permeabilty function'
962) call printErrMsg(option)
963) end select
964) case(THOMEER_COREY)
965) pc = capillary_pressure
966) por = auxvar1
967) perm = auxvar2*1.013202d15 ! convert from m^2 to mD
968) Fg = saturation_function%alpha
969) a = saturation_function%m
970) Pd = 100.d0*por/sqrt(perm/(3.8068d0*(Fg**(-1.334d0)))) ! psi
971) PHg = 9.63051d-4*pc
972) if (PHg > Pd) then
973) saturation = 1.d0-exp(-Fg/log10(PHg/Pd))
974) #if 0
975) alpha = pc*(1.d0+1.d-8)
976) m = 9.63051d-4*alpha
977) n = 1.d0-exp(-Fg/log10(m/Pd))
978) n = (n-saturation)/(alpha-pc)
979) #endif
980) dsat_dpc = (saturation-1.d0)*Fg/(log10(PHg/Pd)**2.d0)/(pc*2.30258509d0)
981) ! Sr assumed to be zero
982) relative_perm = saturation**a
983) dkr_dpc = a*saturation**(a-1.d0)*dsat_dpc
984) else
985) saturation = 1.d0
986) relative_perm = 1.d0
987) switch_to_saturated = PETSC_TRUE
988) return
989) endif
990) case default
991) option%io_buffer = 'Unknown saturation function'
992) call printErrMsg(option)
993) end select
994)
995) dsat_dpres = -dsat_dpc
996) dkr_dpres = -dkr_dpc
997)
998) end subroutine SaturationFunctionCompute2
999)
1000) ! ************************************************************************** !
1001)
1002) subroutine SaturationFunctionCompute3(capillary_pressure,saturation, &
1003) saturation_function, &
1004) option)
1005) !
1006) ! Computes just saturation as a function of
1007) ! capillary pressure
1008) !
1009) ! Author: Glenn Hammond
1010) ! Date: 2/9/12
1011) !
1012) use Option_module
1013)
1014) implicit none
1015)
1016) PetscReal :: capillary_pressure, saturation
1017) PetscReal :: relative_perm_dummy, dsat_dpres_dummy, dkr_dpres_dummy
1018) type(saturation_function_type) :: saturation_function
1019) PetscReal :: auxvar1_dummy, auxvar2_dummy
1020) type(option_type) :: option
1021)
1022) PetscBool :: switch_to_saturated_dummy
1023)
1024) call SaturationFunctionCompute2(capillary_pressure,saturation, &
1025) relative_perm_dummy, &
1026) dsat_dpres_dummy,dkr_dpres_dummy, &
1027) saturation_function, &
1028) auxvar1_dummy,auxvar2_dummy, &
1029) switch_to_saturated_dummy,option)
1030)
1031) end subroutine SaturationFunctionCompute3
1032)
1033) ! ************************************************************************** !
1034)
1035) subroutine SatFuncComputeIcePExplicit(liquid_pressure, temperature, &
1036) ice_saturation, &
1037) liquid_saturation, gas_saturation, &
1038) liquid_relative_perm, dsl_pl, &
1039) dsl_temp, dsg_pl, dsg_temp, dsi_pl, &
1040) dsi_temp, dkr_pl, dkr_temp, &
1041) saturation_function, pth, option)
1042) !
1043) ! Computes the saturation of ice, water vapor
1044) ! and liquid water given the saturation function
1045) ! temperature, water vapor pressure and liquid
1046) ! water pressure
1047) !
1048) ! Model based on Painter, Comp. Geosci, 2011
1049) !
1050) ! Author: Satish Karra, LANL
1051) ! Date: 11/14/11
1052) !
1053)
1054) use Option_module
1055) use PFLOTRAN_Constants_module
1056)
1057) implicit none
1058)
1059) PetscReal :: liquid_pressure, temperature
1060) PetscReal :: ice_saturation, liquid_saturation, gas_saturation
1061) PetscReal :: liquid_relative_perm
1062) PetscReal :: dsl_pl, dsl_temp
1063) PetscReal :: dsg_pl, dsg_temp
1064) PetscReal :: dsi_pl, dsi_temp
1065) PetscReal :: dkr_pl
1066) type(saturation_function_type) :: saturation_function
1067) type(option_type) :: option
1068)
1069) PetscReal :: alpha, lambda, m, n
1070) PetscReal :: pc, Se, one_over_m, Se_one_over_m, dSe_dpc, dkr_dpc
1071) PetscReal :: dkr_dSe, power
1072) PetscReal :: pc_alpha, pc_alpha_n, one_plus_pc_alpha_n
1073) PetscReal :: pc_alpha_neg_lambda
1074) PetscReal :: function_A, function_B
1075) PetscReal :: pc_il, gamma, pc_il_alpha, pc_il_alpha_n, Se_temp
1076) PetscReal :: one_plus_pc_il_alpha_n
1077) PetscReal :: dfunc_A_temp
1078) PetscReal :: dfunc_B_pl
1079) PetscReal :: liq_sat_one_over_m, dkr_ds_liq, dkr_temp
1080) PetscReal :: pth, dSe_dpc_at_pth
1081) PetscReal, parameter :: den_ice = 9.167d2 !in kg/m3 at 273.15K
1082) PetscReal, parameter :: interfacial_tensions_ratio = 2.33
1083) PetscReal, parameter :: T_0 = 273.15d0 !in K
1084)
1085) PetscReal :: dsi_dpl, dsg_dpl, dsl_dpl
1086) PetscReal :: dsi_dT, dsg_dT, dsl_dT
1087)
1088) dsl_pl = 0.d0
1089) dsl_temp = 0.d0
1090) dsg_pl = 0.d0
1091) dsg_temp = 0.d0
1092) dsi_pl = 0.d0
1093) dsi_temp = 0.d0
1094) dkr_pl = 0.d0
1095) dkr_temp = 0.d0
1096) dkr_ds_liq = 0.d0
1097)
1098) ! compute saturation
1099) select case(saturation_function%saturation_function_itype)
1100) case(VAN_GENUCHTEN)
1101) if (liquid_pressure >= option%reference_pressure) then
1102) function_B = 1.d0
1103) dfunc_B_pl = 0.d0
1104) else
1105) alpha = saturation_function%alpha
1106) pc = option%reference_pressure - liquid_pressure
1107) m = saturation_function%m
1108) n = 1.d0/(1.d0 - m)
1109) pc_alpha = pc*alpha
1110) pc_alpha_n = pc_alpha**n
1111) one_plus_pc_alpha_n = 1.d0 + pc_alpha_n
1112) Se = one_plus_pc_alpha_n**(-m)
1113) dSe_dpc = -m*n*alpha*pc_alpha_n/(pc_alpha*one_plus_pc_alpha_n**(m+1))
1114) if (pc >= pth) then
1115) dSe_dpc_at_pth = -m*n*(1.d0 + (alpha*pth)**n)**(-1.d0-m)*(alpha**n*pth**(n-1.d0))
1116) Se = (pc - 1.d8)*dSe_dpc_at_pth
1117) dSe_dpc = dSe_dpc_at_pth
1118) ! write (*,*) option%myrank, 'pc:', pc, 'Se:', Se, 'dSe_dpc', dSe_dpc
1119) endif
1120) function_B = 1.d0/Se
1121) dfunc_B_pl = 1.d0/(Se**(2.d0))*dSe_dpc
1122) endif
1123) if (temperature >= 0.d0) then
1124) function_A = 1.d0
1125) dfunc_A_temp = 0.d0
1126) else
1127) gamma = den_ice*HEAT_OF_FUSION*interfacial_tensions_ratio
1128) pc_il = gamma*(-(temperature))/T_0
1129) alpha = saturation_function%alpha
1130) m = saturation_function%m
1131) n = 1.d0/(1.d0 - m)
1132) pc_il_alpha = pc_il*alpha
1133) pc_il_alpha_n = pc_il_alpha**n
1134) one_plus_pc_il_alpha_n = 1.d0 + pc_il_alpha_n
1135) Se_temp = one_plus_pc_il_alpha_n**(-m)
1136) function_A = 1.d0/Se_temp
1137) dfunc_A_temp = (gamma/T_0)*1.d0/(Se_temp**(2.d0))*(-m)* &
1138) ((one_plus_pc_il_alpha_n)**(-m - 1.d0))*n* &
1139) (pc_il**(n - 1.d0))*(alpha**n)
1140) endif
1141) case default
1142) option%io_buffer = 'Ice module only supports Van Genuchten'
1143) call printErrMsg(option)
1144) end select
1145)
1146) liquid_saturation = 1.d0/(function_A + function_B - 1.d0)
1147) gas_saturation = liquid_saturation*(function_B - 1.d0)
1148) ice_saturation = liquid_saturation*(function_A - 1.d0)
1149)
1150) dsl_pl = - 1.d0/(function_A + function_B - 1.d0)**(2.d0)*(dfunc_B_pl)
1151) dsl_temp = - 1.d0/(function_A + function_B - 1.d0)**(2.d0)*(dfunc_A_temp)
1152)
1153) dsg_pl = dsl_pl*(function_B - 1.d0) + liquid_saturation*dfunc_B_pl
1154) dsg_temp = dsl_temp*(function_B - 1.d0)
1155)
1156) dsi_pl = dsl_pl*(function_A - 1.d0)
1157) dsi_temp = dsl_temp*(function_A - 1.d0) + liquid_saturation*dfunc_A_temp
1158)
1159) if (liquid_saturation > 1.d0) then
1160) print *, option%myrank, 'vG Liquid Saturation > 1:', liquid_saturation
1161) else if (liquid_saturation < 0.d0) then
1162) print *, option%myrank, 'vG Liquid Saturation < 0:', liquid_saturation
1163) endif
1164)
1165) if (gas_saturation > 1.d0) then
1166) print *, option%myrank, 'vG Gas Saturation > 1:', gas_saturation
1167) else if (gas_saturation < 0.d0) then
1168) print *, option%myrank, 'vG Gas Saturation < 0:', gas_saturation
1169) endif
1170)
1171) if (ice_saturation > 1.d0) then
1172) print *, option%myrank, 'vG Ice Saturation > 1:', ice_saturation
1173) else if (ice_saturation < 0.d0) then
1174) print *, option%myrank, 'vG Ice Saturation < 0:', ice_saturation
1175) endif
1176)
1177) select case(saturation_function%permeability_function_itype)
1178) case(MUALEM)
1179) if (liquid_saturation == 1.d0) then
1180) liquid_relative_perm = 1.d0
1181) dkr_ds_liq = 0.d0
1182) else
1183) m = saturation_function%m
1184) one_over_m = 1.d0/m
1185) liq_sat_one_over_m = liquid_saturation**one_over_m
1186) liquid_relative_perm = sqrt(liquid_saturation)* &
1187) (1.d0 - (1.d0 - liq_sat_one_over_m)**m)**2.d0
1188) dkr_ds_liq = 0.5d0*liquid_relative_perm/liquid_saturation + &
1189) 2.d0*liquid_saturation**(one_over_m - 0.5d0)* &
1190) (1.d0 - liq_sat_one_over_m)**(m - 1.d0)* &
1191) (1.d0 - (1.d0 - liq_sat_one_over_m)**m)
1192) endif
1193) dkr_pl = dkr_ds_liq*dsl_pl
1194) dkr_temp = dkr_ds_liq*dsl_temp
1195) case default
1196) option%io_buffer = 'Ice module only supports Mualem'
1197) call printErrMsg(option)
1198) end select
1199)
1200) end subroutine SatFuncComputeIcePExplicit
1201)
1202) ! ************************************************************************** !
1203)
1204) subroutine ComputeVGAndDerivative(alpha,lambda,Pc,S,dS_dPc)
1205) !
1206) ! Evaluates van Genunchten saturation function and
1207) ! its derivative with respect to capillary pressure
1208) ! for given capillary pressure
1209) !
1210) ! Note: Derivative wrt capillary pressure and not liquid pressure
1211) ! is evaluated
1212) !
1213) ! Author: Satish Karra
1214) ! Date: 10/16/12
1215) !
1216)
1217) implicit none
1218)
1219) PetscReal :: alpha, lambda, gamma
1220) PetscReal :: Pc, S, dS_dPc
1221)
1222) gamma = 1.d0/(1.d0 - lambda)
1223) if (Pc > 0.d0) then
1224) S = (1.d0 + (alpha*Pc)**gamma)**(-lambda)
1225) dS_dPc = (-lambda)*((1.d0 + (alpha*Pc)**gamma)**(-lambda - 1.d0))* &
1226) (gamma*alpha*(alpha*Pc)**(gamma - 1.d0))
1227) else
1228) S = 1.d0
1229) dS_dPc = 0.d0
1230) endif
1231)
1232) end subroutine ComputeVGAndDerivative
1233)
1234) ! ************************************************************************** !
1235)
1236) subroutine ComputeInvVGAndDerivative(alpha,lambda,sat,Sinv,dSinv_dsat)
1237) !
1238) ! Evaluates inverse of van Genunchten saturation function
1239) ! and its derivative with respect to saturation at a given saturation
1240) !
1241) ! Author: Satish Karra
1242) ! Date: 10/16/12
1243) !
1244)
1245) implicit none
1246)
1247) PetscReal :: alpha, lambda, gamma
1248) PetscReal :: sat, Sinv, dSinv_dsat
1249)
1250) gamma = 1.d0/(1.d0 - lambda)
1251) if (sat == 1.d0) then
1252) Sinv = 0.d0
1253) dSinv_dsat = 0.d0
1254) else
1255) Sinv = 1.d0/alpha*((sat)**(-1.d0/lambda) - 1.d0)**(1.d0/gamma)
1256) dSinv_dsat = 1.d0/alpha*1.d0/gamma*((sat**(-1.d0/lambda) - 1.d0)**(1.d0/gamma - &
1257) 1.d0))*(-1.d0/lambda)*(sat**(-1.d0/lambda - 1.d0))
1258) endif
1259)
1260) end subroutine ComputeInvVGAndDerivative
1261)
1262) ! ************************************************************************** !
1263)
1264) subroutine CalculateImplicitIceFunc(alpha,lambda,Pcgl,T,s_i,func_val,dfunc_val)
1265) !
1266) ! Evaluates the value of implicit equation whose
1267) ! solution is ice saturation
1268) !
1269) ! Author: Satish Karra
1270) ! Date: 10/16/12
1271) !
1272)
1273) implicit none
1274)
1275) PetscReal :: alpha, lambda
1276) PetscReal :: Pcgl, T, s_i, func_val
1277) PetscReal :: temp_term, PC
1278) PetscReal :: sat, dsat, sat_term
1279) PetscReal :: sat_inv, dsat_inv
1280) PetscReal :: sat_PC, dsat_dpc, dfunc_val
1281) PetscReal, parameter :: beta = 2.33 ! dimensionless -- ratio of surf. tens
1282) PetscReal, parameter :: rho_i = 9.167d2 ! in kg/m^3
1283) PetscReal, parameter :: T_0 = 273.15 ! in K
1284)
1285) if (T >= 0.d0) then ! T is in C
1286) temp_term = 0.d0
1287) else
1288) temp_term = -beta*rho_i*HEAT_OF_FUSION*T/T_0
1289) endif
1290) call ComputeVGAndDerivative(alpha,lambda,Pcgl,sat,dsat)
1291) sat_term = (s_i + (1.d0 - s_i)*sat)
1292) call ComputeInvVGAndDerivative(alpha,lambda,sat_term,sat_inv,dsat_inv)
1293) PC = temp_term + sat_inv
1294) call ComputeVGAndDerivative(alpha,lambda,PC,sat_PC,dsat_dPC)
1295) func_val = (1.d0 - s_i)*sat - sat_PC
1296) dfunc_val = -sat - dsat_dpc*dsat_inv*(1.d0 - sat)
1297)
1298) end subroutine CalculateImplicitIceFunc
1299)
1300) ! ************************************************************************** !
1301)
1302) subroutine CalcPhasePartitionIceNewt(alpha,lambda,Pcgl,T,s_g,s_i,s_l)
1303) !
1304) ! Solves the implicit constitutive relation
1305) ! to calculate saturations of ice, liquid
1306) ! and vapor phases
1307) !
1308) ! Author: Satish Karra
1309) ! Date: 10/16/12
1310) !
1311)
1312) implicit none
1313) PetscReal :: alpha, lambda
1314) PetscReal :: Pcgl, T, s_g, s_i, s_l
1315) PetscReal :: func_val, dfunc_val
1316) PetscReal :: x, x_new, sat, dsat
1317) PetscInt :: iter
1318) PetscReal, parameter :: eps = 1.d-12
1319) PetscInt, parameter :: maxit = 10
1320)
1321)
1322) if (T >= 0.d0) then
1323) x = 0.d0
1324) else
1325) x = 5.d-1 ! Initial guess
1326) do iter = 1,maxit
1327) call CalculateImplicitIceFunc(alpha,lambda,Pcgl,T,x,func_val,dfunc_val)
1328) ! print *, 'iteration:', iter, 'value:', x, 'inormr:', abs(func_val), &
1329) ! 'dfunc_val:', dfunc_val, 'dfunc_val_num:', dfunc_val_num
1330) if (abs(func_val) < eps) exit
1331) x_new = x - func_val/dfunc_val
1332) if (x_new >= 1.d0) then
1333) x_new = 1.d0 - 1.d-8
1334) endif
1335) if (x_new <= 0.d0) then
1336) x_new = 1.d-8
1337) endif
1338) x = x_new
1339) enddo
1340) endif
1341)
1342) call ComputeVGAndDerivative(alpha,lambda,Pcgl,sat,dsat)
1343) s_i = x
1344) s_l = (1.d0 - x)*sat
1345) s_g = 1.d0 - s_l - s_i
1346)
1347) end subroutine CalcPhasePartitionIceNewt
1348)
1349) ! ************************************************************************** !
1350)
1351) subroutine CalcPhasePartitionIceBis(alpha,lambda,Pcgl,T,s_g,s_i,s_l)
1352) !
1353) ! Solves the implicit constitutive relation
1354) ! to calculate saturations of ice, liquid
1355) ! and vapor phases using bisection method
1356) !
1357) ! Author: Satish Karra
1358) ! Date: 10/16/12
1359) !
1360)
1361) implicit none
1362)
1363) PetscReal :: alpha, lambda
1364) PetscReal :: Pcgl, T, s_g, s_i, s_l
1365) PetscReal :: max_s, min_s
1366) PetscReal :: x, F_max_s, F_min_s, sol, F_x
1367) PetscReal :: sat, dsat
1368) PetscInt :: i
1369) PetscReal :: dF_min_s, dF_max_s, dF_x
1370) PetscInt, parameter :: maxit = 100
1371) PetscReal, parameter :: tol = 1.d-12
1372)
1373) max_s = 1.d0
1374) min_s = 0.d0
1375)
1376) if (T >= 0.d0) then
1377) sol = 0.d0
1378) else
1379) do i = 1, maxit
1380) call CalculateImplicitIceFunc(alpha,lambda,Pcgl,T,min_s,F_min_s,dF_min_s)
1381) call CalculateImplicitIceFunc(alpha,lambda,Pcgl,T,max_s,F_max_s,dF_max_s)
1382) if (abs(F_min_s) < tol) then
1383) sol = min_s
1384) exit
1385) endif
1386) if (abs(F_max_s) < tol) then
1387) sol = max_s
1388) exit
1389) endif
1390) x = min_s + (max_s - min_s)/2.d0
1391) call CalculateImplicitIceFunc(alpha,lambda,Pcgl,T,x,F_x,dF_x)
1392) ! print *, i, x, F_x, max_s, min_s
1393) if (abs(F_x) < tol) then
1394) sol = x
1395) exit
1396) endif
1397) if (F_min_s*F_x > 0.d0) then
1398) min_s = x
1399) else
1400) max_s = x
1401) endif
1402) enddo
1403) endif
1404)
1405) call ComputeVGAndDerivative(alpha,lambda,Pcgl,sat,dsat)
1406) s_i = sol
1407) s_l = (1.d0 - sol)*sat
1408) s_g = 1.d0 - s_l - s_i
1409)
1410) end subroutine CalcPhasePartitionIceBis
1411)
1412) ! ************************************************************************** !
1413)
1414) subroutine CalcPhasePartitionIceDeriv(alpha,lambda,Pcgl,T,s_g,s_i,s_l, &
1415) dsg_dpl,dsg_dT,dsi_dpl,dsi_dT, &
1416) dsl_dpl,dsl_dT)
1417) !
1418) ! Solves the implicit constitutive relation
1419) ! to calculate saturations of ice, liquid
1420) ! and vapor phases
1421) !
1422) ! Author: Satish Karra
1423) ! Date: 10/16/12
1424) !
1425)
1426) implicit none
1427)
1428) PetscReal :: alpha, lambda
1429) PetscReal :: Pcgl, T
1430) PetscReal :: dsg_dpl, dsg_dT
1431) PetscReal :: dsi_dpl, dsi_dT
1432) PetscReal :: dsl_dpl, dsl_dT
1433) PetscReal :: s_g, s_i, s_l
1434) PetscReal :: sat, dsat
1435) PetscReal :: sat_inv, dsat_inv
1436) PetscReal :: PC, sat_PC, dsat_dpc
1437) PetscReal :: G, dS_dpl, temp_term
1438) PetscReal :: L, M, N
1439) PetscReal, parameter :: beta = 2.33 ! dimensionless -- ratio of surf. tens
1440) PetscReal, parameter :: rho_i = 9.167d2 ! in kg/m^3
1441) PetscReal, parameter :: T_0 = 273.15 ! in K
1442)
1443) #if 0
1444) PetscReal :: dsi_dpl_num, dsi_dT_num
1445) PetscReal :: dsg_dpl_num, dsg_dT_num
1446) PetscReal :: dsl_dpl_num, dsl_dT_num
1447) PetscReal :: s_g_pinc, s_i_pinc, s_l_pinc
1448) PetscReal :: s_g_Tinc, s_i_Tinc, s_l_Tinc
1449) PetscReal, parameter :: delta = 1.d-8
1450) #endif
1451)
1452) dsi_dpl = 0.d0
1453) dsi_dT = 0.d0
1454) dsg_dpl = 0.d0
1455) dsg_dT = 0.d0
1456) dsl_dpl = 0.d0
1457) dsl_dT = 0.d0
1458)
1459) #if 0
1460) dsi_dpl_num = 0.d0
1461) dsg_dpl_num = 0.d0
1462) dsl_dpl_num = 0.d0
1463) dsi_dT_num = 0.d0
1464) dsg_dT_num = 0.d0
1465) dsl_dT_num = 0.d0
1466) #endif
1467)
1468) ! Calculate the derivatives of saturation with respect to pl
1469) call CalcPhasePartitionIceNewt(alpha,lambda,Pcgl,T,s_g,s_i,s_l)
1470) if (T >= 0.d0) then
1471) temp_term = 0.d0
1472) else
1473) temp_term = -beta*rho_i*HEAT_OF_FUSION*T/T_0
1474) endif
1475) call ComputeInvVGAndDerivative(alpha,lambda,(s_i + s_l),sat_inv,dsat_inv)
1476) PC = temp_term + sat_inv
1477) call ComputeVGAndDerivative(alpha,lambda,PC,sat_PC,dsat_dPC)
1478) call ComputeVGAndDerivative(alpha,lambda,Pcgl,sat,dsat)
1479) G = dsat_dpc*dsat_inv
1480) dS_dpl = dsat*(-1.d0)
1481) if (G == 1.d0) then
1482) dsi_dpl = 0.d0
1483) dsl_dpl = (1.d0 - s_i)*dS_dpl
1484) else
1485) dsi_dpl = (1.d0 - s_i)/(G/(1.d0 - G) + sat)*dS_dpl
1486) dsl_dpl = dsi_dpl*G/(1.d0 - G)
1487) endif
1488) dsg_dpl = -dsi_dpl - dsl_dpl
1489)
1490)
1491) ! Calculate the derivatives of saturation with respect to temp
1492) L = dsat_dpc
1493) if (T >= 0.d0) then
1494) M = 0.d0
1495) else
1496) M = temp_term/T
1497) endif
1498) N = dsat_inv
1499) dsi_dT = -L*M/(L*N + (1.d0 - L*N)*sat)
1500) dsl_dT = -dsi_dT*sat
1501) dsg_dT = -dsi_dT - dsl_dT
1502)
1503) #if 0
1504) ! Numerical derivatives
1505)
1506) call CalcPhasePartitionIceNewt(alpha,lambda,Pcgl*(1.d0 + delta),T,s_g_pinc, &
1507) s_i_pinc,s_l_pinc)
1508)
1509) if (T >= 0.d0) then
1510) dsi_dT_num = 0.d0
1511) dsg_dT_num = 0.d0
1512) dsl_dT_num = 0.d0
1513) else
1514) call CalcPhasePartitionIceNewt(alpha,lambda,Pcgl,T*(1.d0 + delta),s_g_Tinc, &
1515) s_i_Tinc,s_l_Tinc)
1516) dsi_dT_num = (s_i_Tinc - s_i)/(T*delta)
1517) dsg_dT_num = (s_g_Tinc - s_g)/(T*delta)
1518) dsl_dT_num = (s_l_Tinc - s_l)/(T*delta)
1519) endif
1520)
1521) dsi_dpl_num = (s_i_pinc - s_i)/(delta*Pcgl)*(-1.d0) ! -1.d0 factor for dPcgl/dpl
1522) dsg_dpl_num = (s_g_pinc - s_g)/(delta*Pcgl)*(-1.d0)
1523) dsl_dpl_num = (s_l_pinc - s_l)/(delta*Pcgl)*(-1.d0)
1524)
1525) print *, 'analytical-press:', 'dsg_dpl:', dsg_dpl, &
1526) 'dsi_dpl:', dsi_dpl, 'dsl_dpl:', dsl_dpl
1527) print *, 'numerical-press:' , 'dsg_dpl:', dsg_dpl_num, &
1528) 'dsi_dpl:', dsi_dpl_num, 'dsl_dpl:', dsl_dpl_num
1529)
1530) print *, 'analytical-temp:', 'dsg_dT:', dsg_dT, &
1531) 'dsi_dT:', dsi_dT, 'dsl_dT:', dsl_dT
1532) print *, 'numerical-temp:', 'dsg_dT:', dsg_dT_num, &
1533) 'dsi_dT:', dsi_dT_num, 'dsl_dT:', dsl_dT_num
1534)
1535) dsi_dpl = dsi_dpl_num
1536) dsg_dpl = dsg_dpl_num
1537) dsl_dpl = dsl_dpl_num
1538)
1539) dsi_dT = dsi_dT_num
1540) dsg_dT = dsg_dT_num
1541) dsl_dT = dsl_dT_num
1542) #endif
1543)
1544)
1545) end subroutine CalcPhasePartitionIceDeriv
1546)
1547) ! ************************************************************************** !
1548)
1549) subroutine SatFuncComputeIcePKImplicit(pl,T,s_i,s_l,s_g,kr,dsl_dpl, &
1550) dsl_dT,dsg_dpl,dsg_dT,dsi_dpl, &
1551) dsi_dT,dkr_dpl,dkr_dT, &
1552) saturation_function,pth,option)
1553) !
1554) ! Calculates the saturations of water phases
1555) ! and their derivative with respect to liquid
1556) ! pressure
1557)
1558) ! Model based on Painter & Karra implicit model VZJ (2014)
1559) !
1560) ! Author: Satish Karra
1561) ! Date: 10/16/12
1562) !
1563)
1564) use Option_module
1565)
1566) implicit none
1567)
1568) PetscReal :: pl, T
1569) PetscReal :: s_i, s_g, s_l, kr
1570) PetscReal :: dkr_dpl, dkr_dT
1571) PetscReal :: dkr_dsl
1572) PetscReal :: alpha, m, Pcgl
1573) PetscReal :: one_over_m
1574) PetscReal :: liq_sat_one_over_m
1575) PetscReal :: dsl_dpl, dsl_dT
1576) PetscReal :: dsi_dpl, dsi_dT
1577) PetscReal :: dsg_dpl, dsg_dT
1578) PetscReal :: pth
1579)
1580) type(saturation_function_type) :: saturation_function
1581) type(option_type) :: option
1582)
1583)
1584) select case(saturation_function%saturation_function_itype)
1585) case(VAN_GENUCHTEN)
1586) alpha = saturation_function%alpha
1587) Pcgl = option%reference_pressure - pl
1588) m = saturation_function%m
1589) call CalcPhasePartitionIceDeriv(alpha,m,Pcgl,T,s_g,s_i,s_l,dsg_dpl, &
1590) dsg_dT,dsi_dpl,dsi_dT,dsl_dpl, &
1591) dsl_dT)
1592) case default
1593) option%io_buffer = 'Only van Genuchten supported with ice'
1594) call printErrMsg(option)
1595) end select
1596)
1597) ! Check for bounds on saturations
1598) if (s_l > 1.d0) then
1599) print *, option%myrank, 'vG Liquid Saturation > 1:', s_l
1600) else if (s_l < 0.d0) then
1601) print *, option%myrank, 'vG Liquid Saturation < 0:', s_l
1602) endif
1603)
1604) if (s_g > 1.d0) then
1605) print *, option%myrank, 'vG Gas Saturation > 1:', s_g
1606) else if (s_g < 0.d0) then
1607) print *, option%myrank, 'vG Gas Saturation < 0:', s_g
1608) endif
1609)
1610) if (s_i > 1.d0) then
1611) print *, option%myrank, 'vG Ice Saturation > 1:', s_i
1612) else if (s_i < 0.d0) then
1613) print *, option%myrank, 'vG Ice Saturation < 0:', s_i
1614) endif
1615)
1616) ! Calculate relative permeability
1617) select case(saturation_function%permeability_function_itype)
1618) case(MUALEM)
1619) if (s_l == 1.d0) then
1620) kr = 1.d0
1621) dkr_dsl = 0.d0
1622) else
1623) m = saturation_function%m
1624) one_over_m = 1.d0/m
1625) liq_sat_one_over_m = s_l**one_over_m
1626) kr = sqrt(s_l)*(1.d0 - (1.d0 - liq_sat_one_over_m)**m)**2.d0
1627) dkr_dsl = 0.5d0*kr/s_l + &
1628) 2.d0*s_l**(one_over_m - 0.5d0)* &
1629) (1.d0 - liq_sat_one_over_m)**(m - 1.d0)* &
1630) (1.d0 - (1.d0 - liq_sat_one_over_m)**m)
1631) endif
1632) dkr_dpl = dkr_dsl*dsl_dpl
1633) dkr_dT = dkr_dsl*dsl_dT
1634) case default
1635) option%io_buffer = 'Ice module only supports Mualem'
1636) call printErrMsg(option)
1637) end select
1638)
1639) #if 0
1640) write(*,*) 'rank:', option%myrank, 'sl:', s_l, &
1641) 'sg:', s_g, 'si:', s_i, 'dsl_pl:', dsl_dpl, &
1642) 'dsl_temp:', dsl_dT, 'dsg_pl:', dsg_dpl, 'dsg_temp:', dsg_dT, &
1643) 'dsi_pl:', dsi_dpl, 'dsi_temp:', dsi_dT, 'kr:', kr, &
1644) 'dkr_pl:', dkr_dpl, 'dkr_temp:', dkr_dT, 'pl:', pl, 'T:', T
1645) #endif
1646)
1647) end subroutine SatFuncComputeIcePKImplicit
1648)
1649) ! ************************************************************************** !
1650)
1651) subroutine SatFuncComputeIcePKExplicit(pl,T,s_i,s_l,s_g,kr,dsl_dpl, &
1652) dsl_dT,dsg_dpl,dsg_dT,dsi_dpl, &
1653) dsi_dT,dkr_dpl,dkr_dT, &
1654) saturation_function,pth,option)
1655) !
1656) ! Calculates the saturations of water phases
1657) ! and their derivative with respect to liquid
1658) ! pressure, temperature
1659) !
1660) ! Model used: explicit model from Painter & Karra, VJZ (2014).
1661) !
1662) ! Author: Satish Karra
1663) ! Date: 01/13/13
1664) !
1665)
1666) use Option_module
1667)
1668) implicit none
1669)
1670) PetscReal :: pl, T
1671) PetscReal :: s_i, s_g, s_l, kr
1672) PetscReal :: dkr_dpl, dkr_dT
1673) PetscReal :: dkr_dsl
1674) PetscReal :: alpha, m, Pcgl
1675) PetscReal :: one_over_m
1676) PetscReal :: liq_sat_one_over_m
1677) PetscReal :: dsl_dpl, dsl_dT
1678) PetscReal :: dsi_dpl, dsi_dT
1679) PetscReal :: dsg_dpl, dsg_dT
1680) PetscReal :: pth
1681) PetscReal :: S, dS, Sinv, dSinv
1682) PetscReal, parameter :: beta = 2.2 ! dimensionless -- ratio of surf. tension
1683) PetscReal, parameter :: rho_l = 9.998d2 ! in kg/m^3
1684) PetscReal, parameter :: T_0 = 273.15 ! in K
1685) PetscReal, parameter :: L_f = 3.34d5 ! in J/kg
1686) PetscReal :: T_f, theta, X, Y, dS_dX
1687)
1688)
1689) type(saturation_function_type) :: saturation_function
1690) type(option_type) :: option
1691)
1692)
1693) select case(saturation_function%saturation_function_itype)
1694) case(VAN_GENUCHTEN)
1695) T = T + T_0 ! convert to K
1696) alpha = saturation_function%alpha
1697) Pcgl = option%reference_pressure - pl
1698) m = saturation_function%m
1699) call ComputeVGAndDerivative(alpha,m,Pcgl,S,dS)
1700) call ComputeInvVGAndDerivative(alpha,m,S,Sinv,dSinv)
1701) T_f = T_0 - 1.d0/beta*T_0/L_f/rho_l*Sinv
1702) theta = (T-T_0)/T_0
1703) if (T < T_f) then
1704) X = -beta*theta*L_f*rho_l
1705) call ComputeVGAndDerivative(alpha,m,X,s_l,dS_dX)
1706) s_i = 1.d0 - s_l/S
1707) dsl_dT = 1.d0/T_0*dS_dX*(-beta*rho_l*L_f)
1708) dsl_dpl = 0.d0
1709) dsi_dT = -1.d0/S*dsl_dT
1710) dsi_dpl = -1.d0*s_l/S**2*dS
1711) else
1712) s_l = S
1713) s_i = 0.d0
1714) dsl_dpl = -dS
1715) dsl_dT = 0.d0
1716) dsi_dpl = 0.d0
1717) dsi_dT = 0.d0
1718) endif
1719) s_g = 1.d0 - s_l - s_i
1720) dsg_dpl = -dsl_dpl - dsi_dpl
1721) dsg_dT = -dsl_dT - dsi_dT
1722) T = T - T_0 ! change back to C
1723) case default
1724) option%io_buffer = 'Only van Genuchten supported with ice'
1725) call printErrMsg(option)
1726) end select
1727)
1728) ! Check for bounds on saturations
1729) if (s_l > 1.d0) then
1730) print *, option%myrank, 'vG Liquid Saturation > 1:', s_l
1731) else if (s_l < 0.d0) then
1732) print *, option%myrank, 'vG Liquid Saturation < 0:', s_l
1733) endif
1734)
1735) if (s_g > 1.d0) then
1736) print *, option%myrank, 'vG Gas Saturation > 1:', s_g
1737) else if (s_g < 0.d0) then
1738) print *, option%myrank, 'vG Gas Saturation < 0:', s_g
1739) endif
1740)
1741) if (s_i > 1.d0) then
1742) print *, option%myrank, 'vG Ice Saturation > 1:', s_i
1743) else if (s_i < 0.d0) then
1744) print *, option%myrank, 'vG Ice Saturation < 0:', s_i
1745) endif
1746)
1747) ! Calculate relative permeability
1748) select case(saturation_function%permeability_function_itype)
1749) case(MUALEM)
1750) if (s_l == 1.d0) then
1751) kr = 1.d0
1752) dkr_dsl = 0.d0
1753) else
1754) m = saturation_function%m
1755) one_over_m = 1.d0/m
1756) liq_sat_one_over_m = s_l**one_over_m
1757) kr = sqrt(s_l)*(1.d0 - (1.d0 - liq_sat_one_over_m)**m)**2.d0
1758) dkr_dsl = 0.5d0*kr/s_l + &
1759) 2.d0*s_l**(one_over_m - 0.5d0)* &
1760) (1.d0 - liq_sat_one_over_m)**(m - 1.d0)* &
1761) (1.d0 - (1.d0 - liq_sat_one_over_m)**m)
1762) endif
1763) dkr_dpl = dkr_dsl*dsl_dpl
1764) dkr_dT = dkr_dsl*dsl_dT
1765) case default
1766) option%io_buffer = 'Ice module only supports Mualem'
1767) call printErrMsg(option)
1768) end select
1769)
1770) #if 0
1771) print *, '======================================'
1772) write(*,*) 'rank:', option%myrank, 'sl:', s_l, &
1773) 'sg:', s_g, 'si:', s_i, 'dsl_pl:', dsl_dpl, &
1774) 'dsl_temp:', dsl_dT, 'dsg_pl:', dsg_dpl, 'dsg_temp:', dsg_dT, &
1775) 'dsi_pl:', dsi_dpl, 'dsi_temp:', dsi_dT, 'kr:', kr, &
1776) 'dkr_pl:', dkr_dpl, 'dkr_temp:', dkr_dT, 'pl:', pl, 'T:', T
1777) #endif
1778)
1779) end subroutine SatFuncComputeIcePKExplicit
1780)
1781) ! ************************************************************************** !
1782)
1783) subroutine SatFuncComputeIcePKExplicitNoCryo(pl,T,s_i,s_l,s_g,kr,dsl_dpl, &
1784) dsl_dT,dsg_dpl,dsg_dT,dsi_dpl, &
1785) dsi_dT,dkr_dpl,dkr_dT, &
1786) saturation_function,pth,option)
1787) !
1788) ! Calculates the saturations of water phases
1789) ! and their derivative with respect to liquid
1790) ! pressure, temperature
1791) !
1792) ! Model used: explicit model from Painter & Karra, VJZ (2014).
1793) ! Removed Cryosuction. For this, we assume that s_i + s_l = S*(pcgl), that is
1794) ! the total water content does not change
1795) !
1796) ! Author: Satish Karra
1797) ! Date: 01/13/13
1798) !
1799)
1800) use Option_module
1801)
1802) implicit none
1803)
1804) PetscReal :: pl, T
1805) PetscReal :: s_i, s_g, s_l, kr
1806) PetscReal :: dkr_dpl, dkr_dT
1807) PetscReal :: dkr_dsl
1808) PetscReal :: alpha, m, Pcgl
1809) PetscReal :: one_over_m
1810) PetscReal :: liq_sat_one_over_m
1811) PetscReal :: dsl_dpl, dsl_dT
1812) PetscReal :: dsi_dpl, dsi_dT
1813) PetscReal :: dsg_dpl, dsg_dT
1814) PetscReal :: pth
1815) PetscReal :: S, dS, Sinv, dSinv
1816) PetscReal, parameter :: beta = 2.2 ! dimensionless -- ratio of surf. tension
1817) PetscReal, parameter :: rho_l = 9.998d2 ! in kg/m^3
1818) PetscReal, parameter :: T_0 = 273.15 ! in K
1819) PetscReal, parameter :: L_f = 3.34d5 ! in J/kg
1820) PetscReal :: T_f, theta, X, Y, dS_dX
1821)
1822)
1823) type(saturation_function_type) :: saturation_function
1824) type(option_type) :: option
1825)
1826)
1827) select case(saturation_function%saturation_function_itype)
1828) case(VAN_GENUCHTEN)
1829) T = T + T_0 ! convert to K
1830) alpha = saturation_function%alpha
1831) Pcgl = option%reference_pressure - pl
1832) m = saturation_function%m
1833) call ComputeVGAndDerivative(alpha,m,Pcgl,S,dS)
1834) call ComputeInvVGAndDerivative(alpha,m,S,Sinv,dSinv)
1835) T_f = T_0 - 1.d0/beta*T_0/L_f/rho_l*Sinv
1836) theta = (T-T_0)/T_0
1837) if (T < T_f) then
1838) X = -beta*theta*L_f*rho_l
1839) call ComputeVGAndDerivative(alpha,m,X,s_l,dS_dX)
1840) s_i = S - s_l
1841) dsl_dT = 1.d0/T_0*dS_dX*(-beta*rho_l*L_f)
1842) dsl_dpl = 0.d0
1843) dsi_dT = -dsl_dT
1844) dsi_dpl = -dS - dsl_dpl
1845) else
1846) s_l = S
1847) s_i = 0.d0
1848) dsl_dpl = -dS
1849) dsl_dT = 0.d0
1850) dsi_dpl = 0.d0
1851) dsi_dT = 0.d0
1852) endif
1853) s_g = 1.d0 - s_l - s_i
1854) dsg_dpl = -dsl_dpl - dsi_dpl
1855) dsg_dT = -dsl_dT - dsi_dT
1856) T = T - T_0 ! change back to C
1857) case default
1858) option%io_buffer = 'Only van Genuchten supported with ice'
1859) call printErrMsg(option)
1860) end select
1861)
1862) ! Check for bounds on saturations
1863) if (s_l > 1.d0) then
1864) print *, option%myrank, 'vG Liquid Saturation > 1:', s_l
1865) else if (s_l < 0.d0) then
1866) print *, option%myrank, 'vG Liquid Saturation < 0:', s_l
1867) endif
1868)
1869) if (s_g > 1.d0) then
1870) print *, option%myrank, 'vG Gas Saturation > 1:', s_g
1871) else if (s_g < 0.d0) then
1872) print *, option%myrank, 'vG Gas Saturation < 0:', s_g
1873) endif
1874)
1875) if (s_i > 1.d0) then
1876) print *, option%myrank, 'vG Ice Saturation > 1:', s_i
1877) else if (s_i < 0.d0) then
1878) print *, option%myrank, 'vG Ice Saturation < 0:', s_i
1879) endif
1880)
1881) ! Calculate relative permeability
1882) select case(saturation_function%permeability_function_itype)
1883) case(MUALEM)
1884) if (s_l == 1.d0) then
1885) kr = 1.d0
1886) dkr_dsl = 0.d0
1887) else
1888) m = saturation_function%m
1889) one_over_m = 1.d0/m
1890) liq_sat_one_over_m = s_l**one_over_m
1891) kr = sqrt(s_l)*(1.d0 - (1.d0 - liq_sat_one_over_m)**m)**2.d0
1892) dkr_dsl = 0.5d0*kr/s_l + &
1893) 2.d0*s_l**(one_over_m - 0.5d0)* &
1894) (1.d0 - liq_sat_one_over_m)**(m - 1.d0)* &
1895) (1.d0 - (1.d0 - liq_sat_one_over_m)**m)
1896) endif
1897) dkr_dpl = dkr_dsl*dsl_dpl
1898) dkr_dT = dkr_dsl*dsl_dT
1899) case default
1900) option%io_buffer = 'Ice module only supports Mualem'
1901) call printErrMsg(option)
1902) end select
1903)
1904) #if 0
1905) print *, '======================================'
1906) write(*,*) 'rank:', option%myrank, 'sl:', s_l, &
1907) 'sg:', s_g, 'si:', s_i, 'dsl_pl:', dsl_dpl, &
1908) 'dsl_temp:', dsl_dT, 'dsg_pl:', dsg_dpl, 'dsg_temp:', dsg_dT, &
1909) 'dsi_pl:', dsi_dpl, 'dsi_temp:', dsi_dT, 'kr:', kr, &
1910) 'dkr_pl:', dkr_dpl, 'dkr_temp:', dkr_dT, 'pl:', pl, 'T:', T
1911) #endif
1912)
1913) end subroutine SatFuncComputeIcePKExplicitNoCryo
1914)
1915) ! ************************************************************************** !
1916)
1917) subroutine SatFuncComputeIceDallAmico(pl, T, &
1918) p_fh2o, &
1919) dp_fh2o_dP, &
1920) dp_fh2o_dT, &
1921) s_i, s_l, s_g, &
1922) kr, &
1923) dsl_dpl, dsl_dT, &
1924) dsg_dpl, dsg_dT, &
1925) dsi_dpl, dsi_dT, &
1926) dkr_dpl, dkr_dT, &
1927) saturation_function, &
1928) option)
1929) !
1930) ! Calculates the saturations of water phases and their derivative with
1931) ! respect to liquid pressure and temperature.
1932) !
1933) ! Model used: Dall'Amico (2010) and Dall' Amico et al. (2011)
1934) !
1935) ! Author: Gautam Bisht
1936) ! Date: 02/24/14
1937) !
1938)
1939) use Option_module
1940)
1941) implicit none
1942)
1943) PetscReal :: pl, T
1944) PetscReal :: p_fh2o, dp_fh2o_dP, dp_fh2o_dT
1945) PetscReal :: s_i, s_g, s_l
1946) PetscReal :: kr
1947) PetscReal :: dsl_dpl, dsl_dT
1948) PetscReal :: dsi_dpl, dsi_dT
1949) PetscReal :: dsg_dpl, dsg_dT
1950) PetscReal :: dkr_dpl, dkr_dT
1951) type(saturation_function_type) :: saturation_function
1952) type(option_type) :: option
1953)
1954) PetscReal :: Se,Sr
1955) PetscReal :: dkr_dsl, dkr_dSe
1956) PetscReal :: alpha
1957) PetscReal :: m
1958) PetscReal :: Pc0, Pc1
1959) PetscReal :: one_over_m
1960) PetscReal :: liq_sat_one_over_m
1961) PetscReal :: S0,S1
1962) PetscReal :: dS0,dS1
1963) PetscReal :: H, dH_dT
1964) PetscReal :: T_star
1965) PetscReal :: theta
1966) PetscReal :: x
1967) PetscReal :: dummy
1968) PetscBool :: switch
1969) PetscReal :: numer
1970) PetscReal :: denom
1971) PetscReal :: fct
1972) PetscReal :: T_star_th
1973) PetscReal :: T_star_min
1974) PetscReal :: T_star_max
1975)
1976) !PetscReal, parameter :: beta = 2.2 ! dimensionless -- ratio of surf. tension
1977) PetscReal, parameter :: beta = 1 ! dimensionless [assumed as 1.d0]
1978) PetscReal, parameter :: rho_l = 9.998d2 ! in kg/m^3
1979) PetscReal, parameter :: T_0 = 273.15 ! in K
1980) PetscReal, parameter :: L_f = 3.34d5 ! in J/kg
1981) PetscReal, parameter :: k = 1.d6
1982)
1983) T_star_th = 5.d-1 ! [K]
1984)
1985) s_g = 0.d0
1986) dsg_dpl = 0.d0
1987) dsg_dT = 0.d0
1988)
1989) select case(saturation_function%saturation_function_itype)
1990) case(VAN_GENUCHTEN)
1991)
1992) T = T + T_0 ! convert to K
1993) alpha = saturation_function%alpha
1994) m = saturation_function%m
1995)
1996) Pc0 = option%reference_pressure - pl
1997)
1998) T_star = T_0 - 1.d0/beta*T_0/L_f/rho_l*Pc0
1999)
2000) if (T<T_star) then
2001) H = 1.d0
2002) dH_dT = 0.d0
2003) else
2004) H = 0.d0
2005) dH_dT = 0.d0
2006) endif
2007)
2008) !GB: Add an option to swich between step-function and smooth approximation
2009) ! of step function.
2010) !x = (T - T_star)*k
2011) !H = 0.5d0 - atan(x)/PI
2012) T_star_min = T_star - T_star_th
2013) T_star_max = T_star
2014)
2015) if (T < T_star_min) then
2016) H = 1.d0
2017) dH_dT = 0.d0
2018) else if (T > T_star_max) then
2019) H = 0.d0
2020) dH_dT = 0.d0
2021) else
2022) numer = T - T_star_min
2023) denom = T_star_max - T_star_min
2024) fct = 1.d0 - (numer/denom)**2.d0
2025) H = fct**2.d0
2026) dH_dT = -4.d0*numer*fct/(denom*denom)
2027) endif
2028)
2029) theta = (T - T_star)/T_star
2030) Pc1 = Pc0 - beta*theta*L_f*rho_l*H
2031)
2032) p_fh2o = option%reference_pressure - Pc1
2033) dp_fh2o_dT = -beta*L_f*rho_l*H/T_star
2034) dp_fh2o_dP = 1.d0 - T*T_0/T_star/T_star*H
2035) p_fh2o = pl
2036) dp_fh2o_dT = 0.d0
2037) dp_fh2o_dP = 1.d0
2038)
2039) ! dummy and switch are not used here
2040) call SaturationFunctionCompute2(Pc0,S0,dummy,dS0,dummy,saturation_function,dummy,dummy,switch,option)
2041) call SaturationFunctionCompute2(Pc1,S1,dummy,dS1,dummy,saturation_function,dummy,dummy,switch,option)
2042)
2043) ! convert dS/dpsi to dS/dpl
2044) dS0 = -dS0
2045) dS1 = -dS1
2046)
2047) s_l = S1
2048) s_i = S0 - S1
2049)
2050) dsl_dpl = -dS1*(1.0d0 - T*T_0/T_star/T_star*H)
2051) dsi_dpl = -dS0 - dsl_dpl
2052)
2053) dsl_dT = dS1*(-beta*L_f*rho_l*H/T_star - beta*theta*L_f*rho_l*dH_dT)
2054) dsi_dT = -dsl_dT
2055)
2056) T = T - T_0 ! change back to C
2057)
2058) case default
2059) option%io_buffer = 'Only van Genuchten supported with ice'
2060) call printErrMsg(option)
2061) end select
2062)
2063) ! Calculate relative permeability
2064) select case(saturation_function%permeability_function_itype)
2065) case(MUALEM)
2066) Sr = saturation_function%Sr(1)
2067) Se = (s_l-Sr)/(1.0d0-Sr)
2068) if ( abs(Se-1.d0) < 1.0d-12 ) then
2069) kr = 1.d0
2070) dkr_dsl = 0.d0
2071) else
2072) m = saturation_function%m
2073) one_over_m = 1.d0/m
2074) liq_sat_one_over_m = Se**one_over_m
2075) kr = sqrt(Se)*(1.d0 - (1.d0 - liq_sat_one_over_m)**m)**2.d0
2076) dkr_dSe = 0.5d0*kr/Se + &
2077) 2.d0*Se**(one_over_m - 0.5d0)* &
2078) (1.d0 - liq_sat_one_over_m)**(m - 1.d0)* &
2079) (1.d0 - (1.d0 - liq_sat_one_over_m)**m)
2080) dkr_dsl = dkr_dSe / ( 1.0d0 - Sr )
2081) endif
2082) dkr_dpl = dkr_dsl*dsl_dpl
2083) dkr_dT = dkr_dsl*dsl_dT
2084) case default
2085) option%io_buffer = 'Ice module only supports Mualem'
2086) call printErrMsg(option)
2087) end select
2088)
2089) end subroutine SatFuncComputeIceDallAmico
2090)
2091) ! ************************************************************************** !
2092)
2093) subroutine SatFuncGetLiqRelPermFromSat(saturation,relative_perm,dkr_dSe, &
2094) saturation_function,iphase, &
2095) derivative,option)
2096) !
2097) ! Calculates relative permeability from
2098) ! phase saturation
2099) !
2100) ! (1) Chen, J., J.W. Hopmans, M.E. Grismer (1999) "Parameter estimation of
2101) ! of two-fluid capillary pressure-saturation and permeability functions",
2102) ! Advances in Water Resources, Vol. 22, No. 5, pp 479-493,
2103) ! http://dx.doi.org/10.1016/S0309-1708(98)00025-6.
2104) !
2105) ! Author: Glenn Hammond
2106) ! Date: 03/05/11
2107) !
2108)
2109) use Option_module
2110)
2111) implicit none
2112)
2113) PetscReal :: saturation, relative_perm, dkr_dSe
2114) PetscReal :: power, pct_over_pcmax, pc_over_pcmax, pc_log_ratio
2115) PetscReal :: pcmax, one_over_alpha, alpha, lambda
2116) PetscInt :: iphase
2117) type(saturation_function_type) :: saturation_function
2118) PetscBool :: derivative
2119) type(option_type) :: option
2120)
2121) PetscReal :: m, Sr
2122) PetscReal :: Se, one_over_m, Se_one_over_m
2123)
2124) Sr = saturation_function%Sr(iphase)
2125) Se = (saturation-Sr)/(1.d0-Sr)
2126)
2127) ! if saturation is below residual saturation (this is possible with
2128) ! two-phase with dry-out), need to bail out.
2129) if (Se <= 0.d0) then
2130) relative_perm = 0.d0
2131) dkr_dSe = 0.d0
2132) return
2133) else if (Se > 1.d0) then
2134) Se = 1.d0
2135) endif
2136)
2137) ! compute relative permeability
2138) select case(saturation_function%saturation_function_itype)
2139) case(VAN_GENUCHTEN)
2140) ! compute relative permeability
2141) select case(saturation_function%permeability_function_itype)
2142) case(BURDINE)
2143) ! reference #1
2144) m = saturation_function%m
2145) one_over_m = 1.d0/m
2146) Se_one_over_m = Se**one_over_m
2147) relative_perm = Se*Se*(1.d0-(1.d0-Se_one_over_m)**m)
2148) if (derivative) then
2149) dkr_dSe = 2.d0*relative_perm/Se + &
2150) Se*Se_one_over_m*(1.d0-Se_one_over_m)**(m-1.d0)
2151) endif
2152) case(MUALEM)
2153) ! reference #1
2154) m = saturation_function%m
2155) one_over_m = 1.d0/m
2156) Se_one_over_m = Se**one_over_m
2157) relative_perm = sqrt(Se)*(1.d0-(1.d0-Se_one_over_m)**m)**2.d0
2158) if (derivative) then
2159) dkr_dSe = 0.5d0*relative_perm/Se+ &
2160) 2.d0*Se**(one_over_m-0.5d0)* &
2161) (1.d0-Se_one_over_m)**(m-1.d0)* &
2162) (1.d0-(1.d0-Se_one_over_m)**m)
2163) endif
2164) case default
2165) option%io_buffer = 'Unknown relative permeabilty function'
2166) call printErrMsg(option)
2167) end select
2168) case(BROOKS_COREY)
2169) select case(saturation_function%permeability_function_itype)
2170) case(BURDINE)
2171) ! reference #1
2172) lambda = saturation_function%lambda
2173) power = 3.d0+2.d0/lambda
2174) relative_perm = Se**power
2175) dkr_dSe = power*relative_perm/Se
2176) case(MUALEM)
2177) ! reference #1
2178) lambda = saturation_function%lambda
2179) power = 2.5d0+2.d0/lambda
2180) relative_perm = Se**power
2181) dkr_dSe = power*relative_perm/Se
2182) case default
2183) option%io_buffer = 'Unknown relative permeabilty function'
2184) call printErrMsg(option)
2185) end select
2186) case(LINEAR_MODEL)
2187) select case(saturation_function%permeability_function_itype)
2188) case(BURDINE)
2189) relative_perm = Se
2190) if (derivative) then
2191) dkr_dSe = 1.d0
2192) endif
2193) case(MUALEM)
2194) power = 5.d-1
2195) alpha = saturation_function%alpha
2196) one_over_alpha = 1.d0/alpha
2197) pcmax = saturation_function%pcwmax
2198)
2199) pct_over_pcmax = one_over_alpha/pcmax
2200) pc_over_pcmax = 1.d0-(1.d0-pct_over_pcmax)*Se
2201) pc_log_ratio = log(pc_over_pcmax)/log(pct_over_pcmax)
2202) relative_perm = (Se**power)*(pc_log_ratio**2.d0)
2203) case default
2204) option%io_buffer = 'Unknown relative permeabilty function'
2205) call printErrMsg(option)
2206) end select
2207) case(LEVERETT)
2208) select case(saturation_function%permeability_function_itype)
2209) case(FATT_KLIKOFF)
2210) relative_perm = Se**3
2211) if (derivative) then
2212) dkr_dSe = 3.d0 * Se**2
2213) endif
2214) case default
2215) option%io_buffer = 'Unknown relative permeabilty function'
2216) call printErrMsg(option)
2217) end select
2218) end select
2219)
2220) end subroutine SatFuncGetLiqRelPermFromSat
2221)
2222) ! ************************************************************************** !
2223)
2224) subroutine SatFuncGetGasRelPermFromSat(liquid_saturation, &
2225) gas_relative_perm, &
2226) saturation_function,option)
2227) !
2228) ! Calculates gas phase relative permeability from liquid saturation
2229) !
2230) ! (1) Chen, J., J.W. Hopmans, M.E. Grismer (1999) "Parameter estimation of
2231) ! of two-fluid capillary pressure-saturation and permeability functions",
2232) ! Advances in Water Resources, Vol. 22, No. 5, pp 479-493,
2233) ! http://dx.doi.org/10.1016/S0309-1708(98)00025-6.
2234) !
2235) ! Author: Glenn Hammond
2236) ! Date: 03/05/11
2237) !
2238)
2239) use Option_module
2240)
2241) implicit none
2242)
2243) PetscReal :: liquid_saturation
2244) PetscReal :: gas_relative_perm
2245) PetscReal :: power, pct_over_pcmax, pc_over_pcmax, pc_log_ratio
2246) PetscReal :: pcmax, one_over_alpha, alpha, liq_relative_perm
2247) type(saturation_function_type) :: saturation_function
2248) PetscBool :: derivative
2249) type(option_type) :: option
2250)
2251) PetscReal :: Srl, Srg
2252) PetscReal :: S_star, S_hat
2253) PetscReal :: Sg
2254) PetscReal :: lambda
2255) PetscReal :: m
2256)
2257) Srl = saturation_function%Sr(LIQUID_PHASE)
2258) Srg = saturation_function%Sr(GAS_PHASE)
2259) S_star = (liquid_saturation-Srl)/(1.d0-Srl)
2260) S_hat = (liquid_saturation-Srl)/(1.d0-Srl-Srg)
2261)
2262) gas_relative_perm = 0.d0
2263)
2264) if (S_hat >= 1.d0) then
2265) return
2266) else if (S_hat <= 0.d0) then
2267) gas_relative_perm = 1.d0
2268) return
2269) endif
2270) Sg = 1.d0-S_hat
2271)
2272) ! compute relative permeability
2273) select case(saturation_function%saturation_function_itype)
2274) case(VAN_GENUCHTEN)
2275) ! compute relative permeability
2276) m = saturation_function%m
2277) select case(saturation_function%permeability_function_itype)
2278) case(BURDINE)
2279) ! reference #1
2280) gas_relative_perm = Sg*Sg*(1.d0-S_hat**(1.d0/m))**m
2281) case(MUALEM)
2282) ! reference #1
2283) gas_relative_perm = sqrt(Sg)*(1.d0-S_hat**(1.d0/m))**(2.d0*m)
2284) case default
2285) option%io_buffer = 'Unknown relative permeabilty function'
2286) call printErrMsg(option)
2287) end select
2288) case(BROOKS_COREY)
2289) lambda = saturation_function%lambda
2290) select case(saturation_function%permeability_function_itype)
2291) case(BURDINE)
2292) ! reference #1
2293) gas_relative_perm = Sg*Sg* &
2294) (1.d0-S_hat**(1.d0+2.d0/lambda))
2295) case(MUALEM)
2296) ! reference #1
2297) gas_relative_perm = sqrt(Sg)* &
2298) (1.d0-S_hat**(1.d0+1.d0/lambda))**2.d0
2299) case default
2300) option%io_buffer = 'Unknown relative permeabilty function'
2301) call printErrMsg(option)
2302) end select
2303) case(LINEAR_MODEL)
2304) option%io_buffer = &
2305) 'Linear model not yet supported in SatFuncGetGasRelPermFromSat.'
2306) call printErrMsg(option)
2307) select case(saturation_function%permeability_function_itype)
2308) case(BURDINE)
2309) gas_relative_perm = Sg
2310) case(MUALEM)
2311) power = 5.d-1
2312) alpha = saturation_function%alpha
2313) one_over_alpha = 1.d0/alpha
2314) pcmax = saturation_function%pcwmax
2315)
2316) pct_over_pcmax = one_over_alpha/pcmax
2317) pc_over_pcmax = 1.d0-(1.d0-pct_over_pcmax)*S_star
2318) pc_log_ratio = log(pc_over_pcmax)/log(pct_over_pcmax)
2319) liq_relative_perm = (S_star**power)*(pc_log_ratio**2.d0)
2320)
2321) gas_relative_perm = Sg**power * liq_relative_perm * S_hat**(-power)
2322) case default
2323) option%io_buffer = 'Unknown relative permeabilty function'
2324) call printErrMsg(option)
2325) end select
2326) case(LEVERETT)
2327) select case(saturation_function%permeability_function_itype)
2328) case(FATT_KLIKOFF)
2329) gas_relative_perm = Sg**3
2330) case default
2331) option%io_buffer = 'Unknown relative permeabilty function'
2332) call printErrMsg(option)
2333) end select
2334) end select
2335)
2336) end subroutine SatFuncGetGasRelPermFromSat
2337)
2338) ! ************************************************************************** !
2339)
2340) subroutine CapillaryPressureThreshold(saturation_function,cap_threshold,option)
2341) !
2342) ! Computes the capillary pressure threshold
2343) ! after which instead of van Genuchten a linear function is used upto 100 Mpa
2344) ! capillary pressure. The saturation at 100 Mpa is set to zero
2345) ! This threshold value depends only on van Genuchten parameters alpha and lambda
2346) ! This is used mainly for ice problem, so that the pressure doesnt go to large
2347) ! negative values
2348) !
2349) ! Author: Satish Karra, LANL
2350) ! Date: 09/12/12
2351) !
2352)
2353) use Option_module
2354)
2355) implicit none
2356)
2357) PetscReal :: alpha, lambda, cap_threshold
2358) type(option_type) :: option
2359) type(saturation_function_type) :: saturation_function
2360)
2361)
2362) PetscReal :: gamma, p_new, res_value, jac_value, p_old
2363) PetscReal, parameter :: eps = 1.d-8
2364) PetscInt, parameter :: maxit = 100
2365) PetscInt :: i
2366)
2367) alpha = saturation_function%alpha
2368) lambda = saturation_function%m
2369) alpha = alpha*1.d6
2370) gamma = 1.d0/(1.d0 - lambda)
2371)
2372) p_old = 99.d0
2373)
2374)
2375) do i = 1, maxit
2376) call ResidualCapPressThre(p_old,alpha,lambda,gamma,res_value)
2377) call JacobianCapPressThre(p_old,alpha,lambda,gamma,jac_value)
2378) p_new = p_old - res_value/jac_value
2379) !write (*,*) 'rank:', option%myrank, 'iter:', i, 'p_new:', p_new, 'p_old:', p_old, &
2380) ! 'residual:', res_value, 'jacobian:', jac_value
2381) if ((abs(p_new - p_old) < eps)) exit
2382) p_old = p_new
2383) enddo
2384)
2385) cap_threshold = p_old*1.d6 ! convert to Pa
2386)
2387) end subroutine CapillaryPressureThreshold
2388)
2389) ! ************************************************************************** !
2390)
2391) subroutine ResidualCapPressThre(p,alpha,lambda,gamma,res)
2392) !
2393) ! Computes the residual to calculate capillary pressure
2394) ! thresold in the subroutine CapillaryPressureThreshold
2395) !
2396) ! Author: Satish Karra, LANL
2397) ! Date: 09/12/12
2398) !
2399)
2400) implicit none
2401)
2402) PetscReal :: p, alpha, lambda, gamma, res
2403)
2404) res = lambda*gamma*alpha**gamma*p**(gamma-1.0)*1.d2 - &
2405) (alpha*p)**gamma*(1.d0 + gamma*lambda) - 1.d0
2406)
2407) end subroutine ResidualCapPressThre
2408)
2409) ! ************************************************************************** !
2410)
2411) subroutine JacobianCapPressThre(p,alpha,lambda,gamma,jac)
2412) !
2413) ! Computes the jacobian to calculate capillary pressure
2414) ! thresold in the subroutine CapillaryPressureThreshold
2415) !
2416) ! Author: Satish Karra, LANL
2417) ! Date: 09/12/12
2418) !
2419)
2420) implicit none
2421)
2422) PetscReal :: p, alpha, lambda, gamma, jac
2423)
2424) jac = lambda*gamma*alpha**gamma*(gamma - 1.d0)*p**(gamma - 2.d0)*1.d2 - &
2425) alpha**gamma*gamma*p**(gamma - 1.d0)*(1.d0 + gamma*lambda)
2426)
2427)
2428) end subroutine JacobianCapPressThre
2429)
2430) ! ************************************************************************** !
2431)
2432) subroutine SatFuncGetCapillaryPressure(capillary_pressure,saturation, &
2433) temp,saturation_function,option)
2434) !
2435) ! Computes the capillary pressure as a function of
2436) ! pressure
2437) !
2438) ! Author: Glenn Hammond
2439) ! Date: 06/03/09
2440) !
2441)
2442) use Option_module
2443) use Utility_module, only : QuadraticPolynomialEvaluate
2444)
2445) implicit none
2446)
2447) PetscReal :: capillary_pressure, saturation, temp
2448) type(saturation_function_type) :: saturation_function
2449) type(option_type) :: option
2450)
2451) PetscInt :: iphase
2452) PetscReal :: alpha, lambda, m, n, Sr, one_over_alpha
2453) PetscReal :: Se, derivative
2454) PetscReal :: pc_alpha, pc_alpha_n, one_plus_pc_alpha_n
2455) PetscReal :: pc_alpha_neg_lambda, pcmax
2456) PetscReal :: f, sigma, tk, os
2457)
2458) iphase = 1
2459)
2460) Sr = saturation_function%Sr(iphase)
2461) if (saturation <= Sr) then
2462) capillary_pressure = saturation_function%pcwmax
2463) return
2464) else if (saturation >= 1.d0) then
2465) capillary_pressure = 0.d0
2466) return
2467) endif
2468)
2469) select case(saturation_function%saturation_function_itype)
2470) case(VAN_GENUCHTEN)
2471) alpha = saturation_function%alpha
2472) m = saturation_function%m
2473) n = 1.d0/(1.d0-m)
2474) Se = (saturation-Sr)/(1.d0-Sr)
2475) one_plus_pc_alpha_n = Se**(-1.d0/m)
2476) pc_alpha_n = one_plus_pc_alpha_n - 1.d0
2477) pc_alpha = pc_alpha_n**(1.d0/n)
2478) capillary_pressure = pc_alpha/alpha
2479) case(BROOKS_COREY)
2480) alpha = saturation_function%alpha
2481) lambda = saturation_function%lambda
2482) Sr = saturation_function%Sr(iphase)
2483) Se = (saturation-Sr)/(1.d0-Sr)
2484) if (Se > saturation_function%sat_spline_low) then
2485) call QuadraticPolynomialEvaluate(saturation_function% &
2486) sat_spline_coefficients(1:3), &
2487) Se,capillary_pressure,derivative)
2488) else
2489) pc_alpha_neg_lambda = Se
2490) capillary_pressure = (pc_alpha_neg_lambda**(-1.d0/lambda))/alpha
2491) endif
2492) case(LEVERETT)
2493) tk = 273.15d0 + temp
2494) Se = (saturation-Sr)/(1.d0-Sr)
2495) os = 1.d0-Se
2496) f = os*(1.417d0 + os*(-2.120d0 + 1.263d0*os))
2497) sigma = 1.d0 - 0.625d0 * (374.15d0 - tk)/H2O_CRITICAL_TEMPERATURE
2498) sigma = sigma * 0.2358d0 * &
2499) ((374.15d0 - tk)/H2O_CRITICAL_TEMPERATURE)**1.256d0
2500) capillary_pressure = 632455.53d0 * sigma * f
2501) case(LINEAR_MODEL)
2502) alpha = saturation_function%alpha
2503) one_over_alpha = 1.d0/alpha
2504) pcmax = saturation_function%pcwmax
2505) Sr = saturation_function%Sr(iphase)
2506) Se = (saturation-Sr)/(1.d0-Sr)
2507) capillary_pressure = (one_over_alpha-pcmax)*Se + pcmax
2508) #if 0
2509) case(THOMEER_COREY)
2510) pc = option%reference_pressure-pressure
2511) por = auxvar1
2512) perm = auxvar2*1.013202d15 ! convert from m^2 to mD
2513) Fg = saturation_function%alpha
2514) a = saturation_function%m
2515) Pd = 100.d0*por/sqrt(perm/(3.8068d0*(Fg**(-1.334d0)))) ! psi
2516) PHg = 9.63051d-4*pc
2517) if (PHg > Pd) then
2518) saturation = 1.d0-exp(-Fg/log10(PHg/Pd))
2519) #if 0
2520) alpha = pc*(1.d0+1.d-8)
2521) m = 9.63051d-4*alpha
2522) n = 1.d0-exp(-Fg/log10(m/Pd))
2523) n = (n-saturation)/(alpha-pc)
2524) #endif
2525) dsat_dpc = (saturation-1.d0)*Fg/(log10(PHg/Pd)**2.d0)/(pc*2.30258509d0)
2526) ! Sr assumed to be zero
2527) relative_perm = saturation**a
2528) dkr_dpc = a*saturation**(a-1.d0)*dsat_dpc
2529) else
2530) saturation = 1.d0
2531) relative_perm = 1.d0
2532) return
2533) endif
2534) #endif
2535) case default
2536) option%io_buffer = 'Unknown saturation function'
2537) call printErrMsg(option)
2538) end select
2539)
2540) capillary_pressure = min(capillary_pressure,saturation_function%pcwmax)
2541)
2542) end subroutine SatFuncGetCapillaryPressure
2543)
2544) ! ************************************************************************** !
2545)
2546) function SaturationFunctionGetID(saturation_function_list, &
2547) saturation_function_name, &
2548) material_property_name, option)
2549) !
2550) ! Returns the ID of the saturation function named
2551) ! "saturation_function_name"
2552) !
2553) ! Author: Glenn Hammond
2554) ! Date: 01/12/11
2555) !
2556)
2557) use Option_module
2558) use String_module
2559)
2560) type(saturation_function_type), pointer :: saturation_function_list
2561) character(len=MAXWORDLENGTH) :: saturation_function_name
2562) character(len=MAXWORDLENGTH) :: material_property_name
2563) type(option_type) :: option
2564)
2565) PetscInt :: SaturationFunctionGetID
2566) PetscBool :: found
2567) type(saturation_function_type), pointer :: cur_saturation_function
2568)
2569) found = PETSC_FALSE
2570) cur_saturation_function => saturation_function_list
2571) do
2572) if (.not.associated(cur_saturation_function)) exit
2573) if (StringCompare(saturation_function_name, &
2574) cur_saturation_function%name,MAXWORDLENGTH)) then
2575) found = PETSC_TRUE
2576) SaturationFunctionGetID = cur_saturation_function%id
2577) return
2578) endif
2579) cur_saturation_function => cur_saturation_function%next
2580) enddo
2581) if (.not.found) then
2582) option%io_buffer = 'Saturation function "' // &
2583) trim(saturation_function_name) // &
2584) '" in material property "' // &
2585) trim(material_property_name) // &
2586) '" not found among available saturation functions.'
2587) call printErrMsg(option)
2588) endif
2589)
2590) end function SaturationFunctionGetID
2591)
2592) ! ************************************************************************** !
2593)
2594) subroutine SaturationFunctionVerify(saturation_function,option)
2595) !
2596) ! Evaluates saturation function curves for plotting external to PFLOTRAN
2597) !
2598) ! Author: Glenn Hammond
2599) ! Date: 04/28/14
2600) !
2601)
2602) use Option_module
2603)
2604) implicit none
2605)
2606) type(saturation_function_type) :: saturation_function
2607) type(option_type) :: option
2608)
2609) character(len=MAXSTRINGLENGTH) :: string
2610) PetscReal :: pc, pc_increment, pc_max
2611) PetscReal :: sat, dummy_real
2612) PetscInt :: count, i
2613) PetscReal :: x(101), y(101), krl(101), krg(101)
2614)
2615) if (.not.(saturation_function%saturation_function_itype == VAN_GENUCHTEN .or.&
2616) saturation_function%saturation_function_itype == BROOKS_COREY)) then
2617) return
2618) endif
2619)
2620) ! calculate saturation as a function of capillary pressure
2621) ! start at 1 Pa up to maximum capillary pressure
2622) pc_max = saturation_function%pcwmax
2623) pc = 1.d0
2624) pc_increment = 1.d0
2625) count = 0
2626) do
2627) if (pc > pc_max) exit
2628) count = count + 1
2629) call SaturationFunctionCompute(pc,sat,saturation_function,option)
2630) x(count) = pc
2631) y(count) = sat
2632) if (pc > 0.99*pc_increment*10.d0) pc_increment = pc_increment*10.d0
2633) pc = pc + pc_increment
2634) enddo
2635)
2636) write(string,*) saturation_function%name
2637) string = trim(saturation_function%name) // '_pc_sat.dat'
2638) open(unit=86,file=string)
2639) write(86,*) '"capillary pressure", "saturation"'
2640) do i = 1, count
2641) write(86,'(2es14.6)') x(i), y(i)
2642) enddo
2643) close(86)
2644)
2645) ! calculate capillary pressure as a function of saturation
2646) do i = 1, 101
2647) sat = dble(i-1)*0.01d0
2648) call SatFuncGetCapillaryPressure(pc,sat,option%reference_temperature, &
2649) saturation_function,option)
2650) x(i) = sat
2651) y(i) = pc
2652) call SatFuncGetLiqRelPermFromSat(sat,krl(i),dummy_real, &
2653) saturation_function,ONE_INTEGER, &
2654) PETSC_FALSE,option)
2655) call SatFuncGetGasRelPermFromSat(sat,krg(i),saturation_function,option)
2656) enddo
2657) count = 101
2658)
2659) write(string,*) saturation_function%name
2660) string = trim(saturation_function%name) // '_sat_pc.dat'
2661) open(unit=86,file=string)
2662) write(86,*) '"saturation", "capillary pressure"'
2663) do i = 1, count
2664) write(86,'(2es14.6)') x(i), y(i)
2665) enddo
2666) close(86)
2667)
2668) write(string,*) saturation_function%name
2669) string = trim(saturation_function%name) // '_krl.dat'
2670) open(unit=86,file=string)
2671) write(86,*) '"saturation", "liquid relative permeability"'
2672) do i = 1, count
2673) write(86,'(2es14.6)') x(i), krl(i)
2674) enddo
2675) close(86)
2676)
2677) write(string,*) saturation_function%name
2678) string = trim(saturation_function%name) // '_krg.dat'
2679) open(unit=86,file=string)
2680) write(86,*) '"saturation", "gas relative permeability"'
2681) do i = 1, count
2682) write(86,'(2es14.6)') x(i), krg(i)
2683) enddo
2684) close(86)
2685)
2686) end subroutine SaturationFunctionVerify
2687)
2688) ! ************************************************************************** !
2689)
2690) recursive subroutine SaturationFunctionDestroy(saturation_function)
2691) !
2692) ! Destroys a saturation function
2693) !
2694) ! Author: Glenn Hammond
2695) ! Date: 11/02/07
2696) !
2697)
2698) implicit none
2699)
2700) type(saturation_function_type), pointer :: saturation_function
2701)
2702) if (.not.associated(saturation_function)) return
2703)
2704) call SaturationFunctionDestroy(saturation_function%next)
2705)
2706) if (associated(saturation_function%Sr)) deallocate(saturation_function%Sr)
2707) nullify(saturation_function%Sr)
2708)
2709) if (associated(saturation_function%Kr0)) deallocate(saturation_function%Kr0)
2710) nullify(saturation_function%Kr0)
2711)
2712) deallocate(saturation_function)
2713) nullify(saturation_function)
2714)
2715) end subroutine SaturationFunctionDestroy
2716)
2717) end module Saturation_Function_module