mphase_pckr_mod.F90 coverage: 50.00 %func 19.06 %block
1) module Mphase_pckr_module
2)
3) use Material_module
4) use Option_module
5) use PFLOTRAN_Constants_module
6)
7) implicit none
8)
9) private
10) #include "petsc/finclude/petscsys.h"
11) PetscReal, private, parameter :: pckr_sat_water_cut = 1.D0 - 5.D-7
12)
13) public :: pckrNH_noderiv, pckrHY_noderiv
14) contains
15)
16) ! ************************************************************************** !
17)
18) subroutine pflow_pckr(ipckrtype,pckr_swir,pckr_lambda,pckr_alpha,&
19) pckr_m ,pckr_pcmax,sg,pc,pc_s,kr,kr_s,pckr_beta,pckr_pwr)
20)
21) implicit none
22)
23) PetscInt :: ipckrtype
24) PetscReal :: sg
25) PetscReal :: pckr_swir,pckr_lambda,pckr_alpha,pckr_m,pckr_pcmax,pckr_pwr
26) PetscReal :: pc(1:2),pc_s(1:2),kr(1:2),kr_s(1:2)
27) PetscReal :: pckr_beta
28)
29) PetscReal :: sw,se,swir,sw0,lam,ala,um,un,upc,upc_s
30) PetscReal :: temp,pcmax,ser
31) PetscReal :: uum,pckr_betac,betac,st
32)
33) ! if (present(pckr_beta))
34) pckr_betac=pckr_beta
35) sw=1.D0-sg
36) swir=pckr_swir
37) sw0=1.d0
38) pcmax=pckr_pcmax
39)
40) ! print *,'pflow_pckr: ',ipckrtype,sg,pckr_pcmax
41)
42) select case(ipckrtype)
43)
44) case(1) ! van Genuchten
45)
46) ala=pckr_alpha
47) um=pckr_m
48) un=1.D0/(1.D0-um)
49)
50) if (sw>pckr_sat_water_cut) then
51) upc=0.D0; upc_s=0.D0
52) kr(1)=1.d0; kr(2)=0.d0;
53) kr_s(1)=0.d0; kr_s(2)=0.d0;
54) elseif (sw > (1.05D0*swir)) then
55) se=(sw-swir)/(1.D0-swir)
56) temp=se**(-1.D0/um)
57) upc=(temp-1.D0)**(1.d0/un)/ala
58) upc_s=-1.D0/um/un*upc*(se**(-1.D0-1.D0/um))/(se**(-1.D0/um)-1.d0)
59) temp=1.D0/temp
60) kr(1)=sqrt(se)*(1.D0-(1.D0-temp)**um)**2.D0
61) kr_s(1)=0.5d0*kr(1)/se+2.D0*sqrt(se)*(1.d0-(1.d0-temp)**um)* &
62) (1.d0-temp)**(um-1.d0)*temp/se
63) kr(2)=1.D0-kr(1)
64) kr_s(2)= -kr_s(1)
65) ! print *,'in pckr ',um, sw , ala, se,upc,kr
66) else ! use linear extropolation
67) se=(0.05d0*swir)/(1.D0-swir)
68) temp=se**(-1.D0/um)
69) upc=(temp-1.D0)**(1.d0/un)/ala
70) upc_s=-1.D0/um/un*upc*(se**(-1.D0-1.D0/um))/(se**(-1.D0/um)-1.d0)
71) if (sw>swir) then
72) temp=1.D0/temp
73) kr(1)=sqrt(se)*(1.D0-(1.D0-temp)**um)**2.D0
74) kr_s(1)=0.5d0*kr(1)/se+2.d0*sqrt(se)*(1.D0-(1.D0-temp)**um)* &
75) (1.D0-temp)**(um-1.d0)*temp/se
76) ser=(sw-swir)/(1.D0-swir)
77) upc=upc+(ser-se)*upc_s
78) ! print *,se,ser,kr(1),kr_s(1)
79) kr(1)=kr(1)+(ser-se)*kr_s(1)
80) kr(2)=1.D0-kr(1)
81) kr_s(2)= -kr_s(1)
82) else
83) upc=upc-upc_s*se
84) upc_s=0.D0
85) kr(1)=0.d0
86) kr_s(1)=0.d0
87) kr(2)=1.d0
88) kr_s(2)=0.d0
89) end if
90) end if
91)
92) pc(1)=upc; pc(2)=0.d0
93) pc_s(1)=upc_s
94) pc_s(2)=0.d0
95)
96) ! since our primary variable for saturation is sg, and all of the
97) ! derivative here is to se based on sw, so need transfer
98) temp=sw0-swir
99) kr_s(:)=-kr_s(:) / temp
100) pc_s(1)=-pc_s(1) / temp
101)
102) case(2) !Brooks-Corey
103)
104) lam=pckr_lambda
105) ala=pckr_alpha
106)
107) if (sw > (1.05D0*swir)) then
108) se=(sw-swir)/(sw0-swir)
109) upc=se**(-1.D0/lam)/ala
110) upc_s=-upc/se/lam
111) kr(1)=se**(2.d0/lam+3.d0)
112) kr_s(1)=(2.d0/lam+3.d0)*kr(1)/se
113) kr(2)=1.D0-kr(1)
114) kr_s(1)= -kr_s(1)
115) else ! use linear extropolation
116) se=(0.05d0*swir)/(1.D0-swir)
117) upc=se**(-1.D0/lam)/ala
118) upc_s=-upc/se/lam
119) if (sw > swir) then
120) kr(1)=se**(2.d0/lam+3.d0)
121) kr_s(1)=(2.d0/lam+3.d0)*kr(1)/se
122) ser=(sw-swir)/(1.D0-swir)
123) upc=upc+(ser-se)*upc_s
124) kr(1)=kr(1)+(ser-se)*kr_s(1)
125) kr(2)=1.D0-kr(1)
126) kr_s(2)=-kr_s(1)
127) else
128) upc=upc-upc_s*se
129) kr(1)=0.D0
130) kr_s(1)=0.D0
131) kr(2)=1.D0
132) kr_s(2)=0.D0
133) end if
134) end if
135) pc(1)=upc; pc(2)=0.d0
136) pc_s(1)=upc_s
137) pc_s(2)=0.d0
138)
139) ! since our primary variable for saturation is sg, and all of the
140) ! derivative here is to se based on sw, so need transfer
141) temp=sw0-swir
142) kr_s(:)=-kr_s(:) / temp
143) pc_s(1)=-pc_s(1) / temp
144)
145) case(3) !linear interpolation, need pcmax, assign krmax=1.
146)
147) if (sw>swir)then
148) se=(sw-swir)/(sw0-swir)
149) upc=pcmax*(1.D0-se)
150) upc_s= - pcmax
151) kr(1)=se
152) kr(2)=1.D0 - kr(1)
153) kr_s(1)=1.D0
154) kr_s(2)=-1.D0
155) else
156) upc=pcmax
157) upc_s=0.D0
158) kr(1)=0.d0
159) kr(2)=1.d0
160) kr_s(1)=0.d0
161) kr_s(2)=0.d0
162) end if
163)
164) pc(1)=upc; pc(2)=0.d0
165) pc_s(1)=upc_s
166) pc_s(2)=0.d0
167)
168) ! since our primary variable for saturation is sg, and all of the
169) ! derivative here is to se based on sw, so need transfer
170) temp=sw0-swir
171) kr_s(:)=-kr_s(:) / temp
172) pc_s(1)=-pc_s(1) / temp
173)
174)
175) case(4) ! po model with gas phase residual
176)
177) if (sw>1.D0) sw=1.D0
178) if (sw<0.D-0) sw= 0.D-0
179)
180) if (sw> pckr_sat_water_cut)then
181) upc=0.D0; kr(1)=1.0; kr(2)=0.D0;
182) upc_s=0.D0; kr_s(1)=0.D0; kr(2)=0.D0;
183) else
184) ala=pckr_alpha
185) betac=pckr_betac
186) um = pckr_m
187) uum = 1.D0/um
188) un=1.D0/(1.D0-um)
189) ! print *,'pckr:', ala, betac, pckr_m, pckr_pwr, pcmax,swir
190) se=sw
191) temp=se**uum
192) upc=(1.D0/temp - 1.D0)**(1.D0 - um) / ala / betac
193) upc_s= -((1.D0/temp -1.D0)**(-um))*(se**(-uum-1.D0))*uum*(1.D0-um)&
194) /ala/betac
195)
196) se=(sw-swir)/(1.D0-swir)
197) st= 1.D0
198) if (sw>=swir)then
199) kr(1)= sqrt(se)*(1.D0-(1.D0-se**uum)**um)**2.D0
200) ! kr(2)= sqrt(st-se)* ((1.D0-se**uum)**um )**7.D0
201) kr(2)= sqrt(st-se)* ((1.D0-se**uum)**um )**pckr_pwr
202)
203) kr_s(1) = 0.5D0 / sqrt(se)*(1.D0 - (1.D0 - se**uum)**um)**2.D0 &
204) + 2.D0 * sqrt(se) * (1.D0- (1.D0-se**uum)**um) &
205) * ((1.D0 - se**uum)**(um-1.D0)) * (Se**(uum-1.d0))
206) ! kr_s(2) = -0.5D0 / sqrt(1.D0-se) * (1.D0 - se**uum)**(7.D0*um) &
207) ! - 7.D0 * sqrt(1.D0- se) * (1.D0 - se**uum)**(7.D0*um-1.D0) &
208) ! * se**(uum -1.D0)
209) kr_s(2) = -0.5D0 / sqrt(1.D0-se) * (1.D0-se**uum)**(pckr_pwr*um) &
210) - pckr_pwr * sqrt(1.D0-se) * (1.D0-se**uum)**(pckr_pwr*um-1.D0) &
211) * se**(uum-1.D0)
212) else
213) kr(1)=0.D0
214) kr_s(1)=0.D0
215)
216) ! kr(2)= sqrt(st-se)* ((1.D0-se**uum)**um)**7.D0
217) kr(2)= sqrt(st-se)* ((1.D0-se**uum)**um)**pckr_pwr
218) ! kr_s(2) = -0.5D0 / sqrt(1.D0-se) * (1.D0 - se**uum)**(7.D0*um) &
219) ! - 7.D0 * sqrt(1.D0- se) * (1.D0 - se**uum)**(7.D0*um-1.D0) &
220) ! * se**(uum -1.D0)
221) kr_s(2) = -0.5D0 / sqrt(1.D0-se) * (1.D0 - se**uum)**(pckr_pwr*um) &
222) - pckr_pwr * sqrt(1.D0- se) * (1.D0 - se**uum)**(pckr_pwr*um-1.D0) &
223) * se**(uum -1.D0)
224) endif
225) endif
226)
227) if (kr(1)<0.D0) kr(1)=0.D0; kr_s(1)=0.D0!; kr(2)=1.D0; kr_s(2)=0.D0
228) if (kr(1)>1.D0) kr(1)=1.D0; kr_s(1)=0.D0!;kr(2)=0.D0; kr_s(2)=0.D0
229) if (kr(2)<0.D0) kr(2)=0.D0; kr_s(2)=0.D0!;kr(1)=1.D0; kr_s(1)=0.D0
230) if (kr(2)>1.D0) kr(2)=1.D0; kr_s(2)=0.D0!;kr(1)=0.D0; kr_s(1)=0.D0
231)
232)
233) pc(1)=upc; pc(2)=0.d0
234) pc_s(1)=upc_s
235) pc_s(2)=0.d0
236)
237) ! since our primary variable for saturation is sg, and all of the
238) ! derivative here is to se based on sw, so need transfer
239) temp=sw0-swir
240) kr_s(:)=-kr_s(:) / temp
241) pc_s(1)=-pc_s(1)
242)
243) ! if (pc(1)>pcmax) print *, 'pckr4: ',sg,pc,kr,pc_s,kr_s
244) ! if (sw<pckr_sat_water_cut) print *, sg,pc,kr,pc_s,kr_s
245)
246) end select
247)
248) return
249)
250) end subroutine pflow_pckr
251)
252) ! ************************************************************************** !
253) !subroutine pflow_pckr_noderiv_exec(ipckrtype,pckr_sir,pckr_lambda, &
254) ! pckr_alpha,pckr_m,pckr_pcmax,sg,pc,kr,pckr_beta,pckr_pwr)
255)
256) subroutine pflow_pckr_noderiv_exec(ipckrtype,ikrtype,pckr_sir,kr0, &
257) pckr_lambda, pckr_alpha,pckr_m,pckr_pcmax,sg,pc,kr,pckr_beta,pckr_pwr)
258) !
259) ! pckrNH_noderiv: Non-hysteric S-Pc-kr relation excuting routine
260) ! Copied from pflotran_orig
261) !
262) ! Author: Chuan Lu
263) ! Date: 05/12/08
264) !
265) use Saturation_Function_module
266) implicit none
267)
268)
269) PetscInt, intent(in) :: ipckrtype, ikrtype
270) PetscReal, intent(in) :: pckr_sir(:)
271) PetscReal, intent(in) :: kr0(:)
272) PetscReal, intent(in) :: pckr_lambda,pckr_alpha,pckr_m,pckr_pcmax
273) PetscReal, intent(in) :: pckr_beta,pckr_pwr
274) PetscReal, intent(in) :: sg
275) PetscReal, intent(out) :: pc(1:2)
276) PetscReal, intent(out) :: kr(1:2)
277)
278) PetscReal :: se,swir,sgir,sw0,lam,ala,um,un,upc,upc_s,kr_s,krg_s
279) PetscReal :: temp,ser,pcmax,sw
280) PetscReal :: uum,pckr_betac,betac,st
281) PetscReal :: se0,upc0,upc_s0
282)
283) pckr_betac = pckr_beta
284) sw = 1.D0 - sg
285) if (sw > 1.D0) sw = 1.D0
286) if (sw < 0.D0) sw = 0.D0
287)
288) sw0 = 1.d0
289) pcmax = pckr_pcmax
290) swir = pckr_sir(1)
291) sgir = pckr_sir(2)
292) upc = 0.d0
293)
294) select case(ipckrtype)
295)
296) case(0) ! kr = 1
297)
298) if (sw >= 1.01D0*swir) then
299) kr(1) = 1.D0
300) elseif (sw <= swir) then
301) kr(1) = 0.D0
302) else
303) kr(1) = (sw - swir)/swir*1.D2
304) endif
305)
306) if (sg >= 1.01D0*sgir) then
307) kr(2) = 1.D0
308) elseif (sg <= sgir) then
309) kr(2) = 0.D0
310) else
311) kr(2) = (sg - sgir)/sgir*1.D2
312) endif
313)
314) kr = 1.D0
315)
316) upc = 0.D0
317)
318) case(VAN_GENUCHTEN) ! van Genuchten
319)
320) ala = pckr_alpha
321) um = pckr_m
322) un = 1.D0/(1.D0 - um)
323) if (sw > pckr_sat_water_cut) then
324) upc = 0.D0; kr(1) = 1.d0; kr(2) = 0.d0;
325) elseif (sw > (1.05D0*swir)) then
326) if (sw <= 0.99D0) then
327) se = (sw - swir)/(1.D0 - swir)
328) temp = se**(-1.D0/um)
329) upc = (temp - 1.D0)**(1.d0/un)/ala
330) kr(1) = sqrt(se)*(1.D0 - (1.D0 - 1.D0/temp)**um)**2.d0
331) ! kr(2) = 1.D0-kr(1)
332) kr(2) = sqrt(1.D0 - se)*((1.D0 - se**(1.D0/um))**um)**2.D0
333) ! print *,'in pckr nond ',sw,se,upc,kr
334) else
335) se = (sw - swir)/(1.D0 - swir)
336) temp = se**(-1.D0/um)
337) kr(1) = sqrt(se)*(1.D0 - (1.D0 - 1.D0/temp)**um)**2.d0
338) kr(2) = sqrt(1.D0 - se)*((1.D0 - se**(1.D0/um))**um)**2.D0
339) se = (0.99D0 - swir)/(1.D0 - swir)
340) temp = se**(-1.D0/um)
341) upc = (temp - 1.D0)**(1.d0/un)/ala
342) ! kr(1)=sqrt(se)*(1.D0-(1.D0-1.D0/temp)**um)**2.d0
343) ! kr(2)=sqrt(1.D0 - se)*((1.D0-se**(1.D0/um))**um)**2.D0
344) upc = upc*(1.D0 - sw)/1.D-2
345) ! kr(1) = kr(1) * (1.D0-sw)/5D-3
346) ! kr(2) = kr(2) * (1.D0-sw)/5D-3
347) endif
348)
349) else ! use linear extropolation
350) se0 = (0.05D0*swir)/(1.D0 - swir)
351) temp = se0**(-1.D0/um)
352) upc0 = (temp - 1.D0)**(1.d0/un)/ala
353) upc_s0 = -1.D0/um/un*upc0*(se0**(-1.D0 - 1.D0/um))/(se0**(-1.D0/um) - 1.d0)
354) upc_s0 = upc_s0 /(1.D0 - swir)
355) if (sw > swir) then
356) se = (sw - swir)/(1.D0 - swir)
357) temp = se**(-1.D0/um)
358) kr(1) = sqrt(se)*(1.D0 - (1.D0 - 1.D0/temp)**um)**2.D0
359) kr(2) = sqrt(1.D0 - se)*((1.D0 - se**(1.D0/um))**um)**2.D0
360) ! temp=1.D0/temp
361) ! kr_s=0.5d0*kr(1)/se+2.d0*sqrt(se)*(1.d0-(1.d0-temp)**um)* &
362) ! (1.d0-temp)**(um-1.d0)*temp/se
363) ! krg_s=-0.5D0/(1.D0-se)*kr(2) -2.D0*sqrt(1.D0-se)*((1.D0-se**(1.D0/um))**um)&
364) ! *((1.D0-se**(1.D0/um))**(um-1.D0)) * (se**(1.D0/um-1.D0))
365) ! ser=(sw-swir)/(1.D0-swir)
366) upc = upc0 + (sw - 1.05D0 * swir) * upc_s0
367) ! kr(1)=kr(1)+(ser-se)*kr_s
368) ! kr(2)=kr(2)+ (ser-se)*krg_s
369)
370) ! kr(2)=1.D0-kr(1)
371) else
372) upc = upc0 + (sw - 1.05D0 * swir) * upc_s0
373) kr(1) = 0.D0
374) kr(2) = 1.D0
375) end if
376) end if
377)
378) case(BROOKS_COREY) !Brooks-Corey
379)
380) lam = pckr_lambda
381) ala = pckr_alpha
382) ! swir=pckr_swir
383)
384) if (sw > (1.05D0*swir)) then
385) se = (sw - swir)/(sw0 - swir)
386) upc = se**(-1.D0/lam)/ala
387) kr(1) = se**(2.d0/lam+3.d0)
388) !kr(2) = 1.D0- kr(1)
389) kr(2) = (1.D0 - se)**2.D0 * (1.D0 - se**(2.D0/lam + 1.D0))
390) else ! use linear extropolation
391) se0 = (0.05d0*swir)/(1.D0-swir)
392) upc0 = se0**(-1.D0/lam)/ala
393) upc_s0 = -upc0/se0/lam
394) upc_s0 = upc_s0 /(1.D0 - swir)
395) if (sw > swir) then
396) se = (sw - swir)/(1.D0 - swir)
397) kr(1) = se**(2.D0/lam+3.d0)
398) kr(2) = (1.D0 - se)**2.D0 * (1.D0 - se**(2.D0/lam + 1.D0))
399) upc = upc0 + (sw - 1.05D0 * swir) * upc_s0
400)
401) ! kr_s=(2.d0/lam+3.d0)*kr(1)/se
402) ! krg_s = -2.D0*kr(2)/(1.D0-se) -(2.D0+lam)/lam*(1.D0-se)**2.D0*(se**(2.D0/lam))
403) ! ser=(sw-swir)/(1.D0-swir)
404)
405) ! kr(1)=kr(1)+(ser-se)*kr_s
406) ! kr(2)=kr(2)+(ser-se)*krg_s
407) ! kr(2)=1.D0-kr(1)
408) else
409) upc = upc0 + (sw - 1.05D0 * swir) * upc_s0
410) kr(1) = 0.D0
411) kr(2) = 1.D0
412) end if
413) end if
414)
415)
416) case(THOMEER_COREY) !linear interpolation, need pcmax, assign krmax=1.
417)
418) if (sw > swir) then
419) se = (sw - swir)/(sw0 - swir)
420) upc = pcmax * (1.D0 - se)
421) kr(1) = se*se
422) ! kr(2)=1.D0 - kr(1)
423) else
424) upc = pcmax
425) kr(1) = 0.d0
426) ! kr(2)=1.d0
427) end if
428)
429) if (sg > sgir) then
430) se = (sg - sgir)/(1.D0 - sgir)
431) ! upc=pcmax*(1.D0-se)
432) ! kr(1)=se
433) kr(2) = se*se
434) else
435) ! upc=pcmax
436) ! kr(1)=0.d0
437) kr(2) = 0.d0
438) end if
439)
440) case(NMT_EXP) ! po model with gas phase residual (NMT)
441)
442) if (sw > 1.D0) sw = 1.D0
443) if (sw < 0.D-0) sw = 1.D-5
444) if (sw > pckr_sat_water_cut) then
445) upc = 0.D0; kr(1) = 1.0d0; kr(2) = 0.D0;
446) else
447) ala = pckr_alpha
448) betac = pckr_betac
449) um = pckr_m
450) uum = 1.D0/um
451) un = 1.D0/(1.D0-um)
452) !print *,'pckr:', ala, betac, pckr_m, pckr_pwr, pcmax,swir
453) se = sw
454) temp = se**uum
455) if (sw < 0.95D0) then
456) upc = (1.D0/temp - 1.D0)**(1.D0 - um) / ala / betac
457) else
458) temp = 0.95D0**uum
459) upc = (1.D0/temp - 1.D0)**(1.D0 - um) / ala / betac
460) upc = upc/0.05D0 * (1.D0 - se)
461) endif
462)
463) if (upc > pcmax) upc = pcmax
464) se = (sw - swir)/(1.D0 - swir)
465) st = 1.D0
466) if (sw >= swir) then
467) kr(1) = sqrt(se)*(1.D0 - (1.D0 - se**uum)**um)**2.D0
468) ! kr(2)= sqrt(st-se)*((1.D0-se**uum)**um)**7.D0
469) kr(2) = sqrt(st-se)*((1.D0-se**uum)**um)**pckr_pwr
470) else
471) ! if (se <= 0.D0) se = 1.D-7
472) kr(1) = 0.D0
473) ! kr(2) = sqrt(st-se)*((1.D0-se**uum)**um)**7.D0
474) kr(2) = sqrt(st-se)!*((1.D0-se**uum)**um)**pckr_pwr
475) endif
476) endif
477)
478) if (kr(1) < 0.D0) kr(1)=0.D0!; kr(2)=1.D0
479) if (kr(1) > 1.D0) kr(1)=1.D0!; kr(2)=0.D0
480) if (kr(2) < 0.D0) kr(2)=0.D0!; kr(1)=1.D0
481) if (kr(2) > 1.D0) kr(2)=1.D0!; kr(1)=0.D0
482)
483) case(PRUESS_1) !linear interpolation, need pcmax, assign krmax=1 (Pruess_1).
484)
485) if (sw > swir) then
486) se = (sw - swir)/(sw0 - swir)
487) upc = pcmax * (1.D0 - se)
488) kr(1) = se
489) else
490) upc = pcmax
491) kr(1) = 0.d0
492) end if
493)
494) if (sg > sgir) then
495) se = (sg - sgir)/(1.D0 - sgir)
496) kr(2) = se
497) if (kr(2) > 1.D0) kr(2) = 1.D0
498) else
499) kr(2) = 0.d0
500) end if
501)
502) case(VAN_GENUCHTEN_PARKER) ! van Genuchten-Parker
503)
504) lam = pckr_lambda
505) ala = pckr_alpha
506)
507) ! Water phase using van Genuchten
508) um = pckr_m
509) un = 1.D0/(1.D0 - um)
510)
511) se = (sw - swir)/(sw0 - swir)
512)
513) if (sw > swir) then
514) temp = se**(-1.D0/um)
515) ! if (temp < 1.D0+1e-6) temp = 1.D0+1e-6
516) upc = (temp - 1.D0)**(1.d0/un)/ala
517) if (upc > pcmax) upc = pcmax
518)
519) ! Mualem rel. perm.
520) kr(1) = sqrt(se)*(1.D0 - (1.D0-1.D0/temp)**um)**2.d0
521) else
522) upc = pcmax
523) se = 0.D0
524) kr(1) = 0.D0
525) endif
526)
527) ! Gas phase using BC
528)
529) se = (sw - swir)/(1.D0 - swir - sgir)
530)
531) if (sw < swir) then
532) se = 0.D0
533) kr(2) = 1.0
534) elseif (se >= 1.D0-1d-6) then
535) kr(2) = 0.D0
536) else
537) kr(2) = (1.D0 - se)**0.5D0 * (1.D0 - se**(1.D0/um))**(2.D0*um)
538) ! kr(2) = (1.D0 - se)**0.33333333D0 * (1.D0 - se**(1.D0/um))**(2.D0*um)
539) endif
540)
541) case(VAN_GENUCHTEN_DOUGHTY) ! Doughty (2007) - Van Genutchen Mualem model adjusted for sgir/=0
542) ! Modeling geologic storage of carbon dioxide: Comparison of non-hysteretic and
543) ! hysteretic characteristic curves, Doughty (2007)
544) ala = pckr_alpha
545) um = pckr_m
546) un = 1.D0/(1.D0 - um)
547) if (sw > pckr_sat_water_cut) then
548) upc = 0.D0; kr(1) = 1.d0; kr(2) = 0.d0;
549) elseif (sw > (1.05D0*swir)) then
550) select case(ikrtype)
551) case(2)
552) ! Krl Van Genuchten-Mualem
553) se = (sw - swir)/(1.D0 - swir)
554) temp = se**(-1.D0/um)
555) kr(1) = sqrt(se)*(1.D0 - (1.D0 - 1.D0/temp)**um)**2.d0
556) ! Krg Corey - TOUGH2
557) se = (sw - swir)/(1.D0 - sgir - swir)
558) kr(2) = (1.D0 - se)**2.0D0 * (1.D0 - se**(2.D0))
559) end select
560) !! capillary pressure
561) if (sw <= ( 0.99*(1 - sgir) )) then
562) se = (sw - swir)/(1.D0 - sgir - swir)
563) temp = se**(-1.D0/um)
564) upc = (temp - 1.D0)**(1.d0/un)/ala
565) else ! for sw > (1 - sgir) linear interpolation to pc = 0
566) se = ( 0.99*(1.0D0 - sgir) - swir)/(1.D0 - sgir - swir)
567) temp = se**(-1.D0/um)
568) upc = (temp - 1.D0)**(1.d0/un)/ala
569) upc = upc*(1.D0 - sw)/(1.0D0 - 0.99*(1.0D0 - sgir))
570) endif
571)
572) else ! use linear extropolation
573) se0 = (0.05D0*swir)/(1.D0 - sgir - swir)
574) temp = se0**(-1.D0/um)
575) upc0 = (temp - 1.D0)**(1.d0/un)/ala
576) upc_s0 = -1.D0/um/un*upc0*(se0**(-1.D0 - 1.D0/um))/(se0**(-1.D0/um) - 1.d0)
577) upc_s0 = upc_s0 /(1.D0 - sgir - swir)
578) if (sw > swir) then
579) select case(ikrtype)
580) case(2)
581) ! Krl Van Genutchen-Mualem
582) se = (sw - swir)/(1.D0 - swir)
583) temp = se**(-1.D0/um)
584) kr(1) = sqrt(se)*(1.D0 - (1.D0 - 1.D0/temp)**um)**2.d0
585) ! Krg Corey - TOUGH2
586) se = (sw - swir)/(1.D0 - sgir - swir)
587) kr(2) = (1.D0 - se)**2.0D0 * (1.D0 - se**(2.D0))
588) end select
589) !! pressure
590) upc = upc0 + (sw - 1.05D0 * swir) * upc_s0
591) else
592) upc = upc0 + (sw - 1.05D0 * swir) * upc_s0
593) kr(1) = 0.D0
594) kr(2) = 1.D0
595) end if
596) end if
597) !if (upc > pcmax) upc = pcmax
598) end select
599)
600) ! scaling kr with end points
601) kr(1) = kr(1) * Kr0(1)
602) kr(2) = kr(2) * Kr0(2)
603)
604) pc(1) = upc; pc(2) = 0.d0;
605)
606) return
607)
608) end subroutine pflow_pckr_noderiv_exec
609)
610) ! ************************************************************************** !
611)
612) subroutine pckrNH_noderiv(sat, pc, kr, saturation_function, option)
613) !
614) ! pckrHY_noderiv: Hysteric S-Pc-kr relation driver
615) !
616) ! Author: Chuan Lu
617) ! Date: 05/12/08
618) !
619)
620) use Saturation_Function_module
621)
622) implicit none
623)
624) type(saturation_function_type) :: saturation_function
625) type(option_type) :: option
626) PetscReal :: sat(option%nphase),pc(option%nphase),kr(option%nphase)
627)
628) PetscReal :: pckr_sir(option%nphase)
629) PetscReal :: Kr0(option%nphase)
630) PetscReal :: pckr_lambda, &
631) pckr_alpha,pckr_m,pckr_pcmax,sg ,pckr_beta,pckr_pwr
632)
633)
634) pckr_sir(:) = saturation_function%Sr(:)
635) Kr0(:) = saturation_function%Kr0(:)
636) pckr_m = saturation_function%m
637) pckr_lambda = saturation_function%lambda
638) pckr_alpha = saturation_function%alpha
639) pckr_pcmax = saturation_function%pcwmax
640) pckr_beta = saturation_function%betac
641) pckr_pwr = saturation_function%power
642)
643) sg = sat(2)
644)
645) ! call pflow_pckr_noderiv_exec(saturation_function%saturation_function_itype,&
646) ! pckr_sir,pckr_lambda, pckr_alpha,pckr_m,pckr_pcmax,sg,pc,kr,pckr_beta,pckr_pwr)
647)
648) call pflow_pckr_noderiv_exec(saturation_function%saturation_function_itype,&
649) saturation_function%permeability_function_itype, pckr_sir, Kr0, &
650) pckr_lambda, pckr_alpha,pckr_m,pckr_pcmax, sg,pc,kr,pckr_beta,pckr_pwr)
651)
652) end subroutine pckrNH_noderiv
653)
654) ! ************************************************************************** !
655)
656) subroutine pckrHY_noderiv(sat, hysdat, pc, kr, saturation_function, option)
657) !
658) ! Hysteric S-Pc-kr relation driver
659) !
660) ! Author: Chuan Lu
661) ! Date: 05/12/08
662) !
663)
664) use Saturation_Function_module
665)
666) implicit none
667)
668) type(saturation_function_type) :: saturation_function
669) type(option_type) :: option
670) PetscReal :: sat(option%nphase),pc(option%nphase),kr(option%nphase)
671) PetscReal :: hysdat(:)
672)
673) pc=0.D0
674) kr=1.D0
675)
676) end subroutine pckrHY_noderiv
677)
678)
679)
680) end module Mphase_pckr_module