immis_aux.F90 coverage: 0.00 %func 0.00 %block
1) module Immis_Aux_module
2)
3) use Mphase_pckr_module
4)
5) use PFLOTRAN_Constants_module
6)
7) implicit none
8)
9) private
10)
11) #include "petsc/finclude/petscsys.h"
12)
13) type, public :: Immis_auxvar_elem_type
14) PetscReal :: pres
15) PetscReal :: temp
16) PetscReal , pointer :: sat(:)
17) PetscReal , pointer :: den(:)
18) PetscReal , pointer :: avgmw(:)
19) PetscReal , pointer :: vis(:)
20) PetscReal , pointer :: h(:)
21) PetscReal , pointer :: u(:)
22) PetscReal , pointer :: pc(:)
23) PetscReal , pointer :: kvr(:)
24) ! PetscReal , pointer :: xmol(:)
25) ! PetscReal , pointer :: diff(:)
26) PetscReal , pointer :: hysdat(:)
27) PetscReal :: zco2
28) ! PetscReal :: vis
29) ! PetscReal :: dvis_dp
30) ! PetscReal :: kr
31) ! PetscReal :: dkr_dp
32) end type Immis_auxvar_elem_type
33)
34) type, public :: Immis_auxvar_type
35)
36) type(Immis_auxvar_elem_type), pointer :: auxvar_elem(:)
37) #if 0
38) PetscReal , pointer :: davgmw_dx(:)
39) PetscReal , pointer :: dden_dp(:)
40) PetscReal , pointer :: dden_dt(:)
41) PetscReal , pointer :: dden_dx(:)
42) PetscReal , pointer :: dkvr_dp(:)
43) PetscReal , pointer :: dkvr_dt(:)
44) PetscReal , pointer :: dkvr_ds(:)
45) PetscReal , pointer :: dkvr_dx(:)
46) PetscReal , pointer :: dh_dp(:)
47) PetscReal , pointer :: dh_dt(:)
48) PetscReal , pointer :: dh_dx(:)
49) PetscReal , pointer :: du_dp(:)
50) PetscReal , pointer :: du_dt(:)
51) PetscReal , pointer :: du_dx(:)
52) #endif
53) end type Immis_auxvar_type
54)
55) type, public :: Immis_parameter_type
56) PetscReal, pointer :: dencpr(:)
57) PetscReal, pointer :: ckwet(:)
58) PetscReal, pointer :: ckdry(:)
59) PetscReal, pointer :: sir(:,:)
60) end type Immis_parameter_type
61)
62) type, public :: Immis_type
63) PetscInt :: n_zero_rows
64) PetscInt, pointer :: zero_rows_local(:), zero_rows_local_ghosted(:)
65)
66) PetscBool :: auxvars_up_to_date
67) PetscBool :: inactive_cells_exist
68) PetscInt :: num_aux, num_aux_bc, num_aux_ss
69) type(Immis_parameter_type), pointer :: immis_parameter
70) type(Immis_auxvar_type), pointer :: auxvars(:)
71) type(Immis_auxvar_type), pointer :: auxvars_bc(:)
72) type(Immis_auxvar_type), pointer :: auxvars_ss(:)
73)
74) PetscReal, pointer :: res_old_AR(:,:), res_old_FL(:,:), delx(:,:)
75) end type Immis_type
76)
77)
78)
79) public :: ImmisAuxCreate, ImmisAuxDestroy, &
80) ImmisAuxVarCompute_NINC, ImmisAuxVarCompute_WINC,&
81) ImmisAuxVarInit, ImmisAuxVarCopy
82)
83) contains
84)
85) ! ************************************************************************** !
86)
87) function ImmisAuxCreate()
88) !
89) ! ImmisAuxVarCreate: Allocate and initialize auxiliary object
90) !
91) ! Author: Chuan Lu
92) ! Date: 02/27/08
93) !
94)
95) use Option_module
96)
97) implicit none
98)
99) type(Immis_type), pointer :: ImmisAuxCreate
100)
101) type(Immis_type), pointer :: aux
102)
103) allocate(aux)
104) aux%auxvars_up_to_date = PETSC_FALSE
105) aux%inactive_cells_exist = PETSC_FALSE
106) aux%num_aux = 0
107) aux%num_aux_bc = 0
108) nullify(aux%auxvars)
109) nullify(aux%auxvars_bc)
110) aux%n_zero_rows = 0
111) allocate(aux%immis_parameter)
112) nullify(aux%immis_parameter%sir)
113) nullify(aux%immis_parameter%ckwet)
114) nullify(aux%immis_parameter%dencpr)
115) nullify(aux%zero_rows_local)
116) nullify(aux%zero_rows_local_ghosted)
117) nullify(aux%res_old_FL)
118) nullify(aux%res_old_AR)
119) nullify(aux%delx)
120)
121) ImmisAuxCreate => aux
122)
123) end function ImmisAuxCreate
124)
125) ! ************************************************************************** !
126)
127) subroutine ImmisAuxVarInit(auxvar,option)
128) !
129) ! Initialize auxiliary object
130) !
131) ! Author: Chuan Lu
132) ! Date: 02/14/08
133) !
134)
135) use Option_module
136)
137) implicit none
138)
139) type(Immis_auxvar_type) :: auxvar
140) type(option_type) :: option
141)
142) PetscInt :: var_elem_size, var_node_size
143) PetscInt :: nvar
144)
145) allocate(auxvar%auxvar_elem(0 : option%nflowdof))
146) allocate(auxvar%auxvar_elem(0)%hysdat(4))
147)
148) do nvar = 0, option%nflowdof
149) allocate ( auxvar%auxvar_elem(nvar)%sat(option%nphase))
150) allocate ( auxvar%auxvar_elem(nvar)%den(option%nphase))
151) allocate ( auxvar%auxvar_elem(nvar)%avgmw(option%nphase))
152) allocate ( auxvar%auxvar_elem(nvar)%vis(option%nphase))
153) allocate ( auxvar%auxvar_elem(nvar)%h(option%nphase))
154) allocate ( auxvar%auxvar_elem(nvar)%u(option%nphase))
155) allocate ( auxvar%auxvar_elem(nvar)%pc(option%nphase))
156) allocate ( auxvar%auxvar_elem(nvar)%kvr(option%nphase))
157) ! allocate ( auxvar%auxvar_elem(nvar)%xmol(option%nphase*option%nflowspec))
158) ! allocate ( auxvar%auxvar_elem(nvar)%diff(option%nphase*option%nflowspec))
159) if (nvar>0) &
160) auxvar%auxvar_elem(nvar)%hysdat => auxvar%auxvar_elem(0)%hysdat
161)
162) auxvar%auxvar_elem(nvar)%pres = 0.d0
163) auxvar%auxvar_elem(nvar)%temp = 0.d0
164) auxvar%auxvar_elem(nvar)%sat = 0.d0
165) auxvar%auxvar_elem(nvar)%den = 0.d0
166) auxvar%auxvar_elem(nvar)%avgmw = 0.d0
167) auxvar%auxvar_elem(nvar)%vis = 0.d0
168) auxvar%auxvar_elem(nvar)%h = 0.d0
169) auxvar%auxvar_elem(nvar)%u = 0.d0
170) auxvar%auxvar_elem(nvar)%pc = 0.d0
171) auxvar%auxvar_elem(nvar)%kvr = 0.d0
172) ! auxvar%auxvar_elem(nvar)%xmol = 0.d0
173) ! auxvar%auxvar_elem(nvar)%diff = 0.d0
174) #if 0
175) auxvar%auxvar_elem(nvar)%dsat_dp = 0.d0
176) auxvar%auxvar_elem(nvar)%dden_dp = 0.d0
177) auxvar%auxvar_elem(nvar)%dkvr_dp = 0.d0
178) #endif
179) enddo
180)
181) end subroutine ImmisAuxVarInit
182)
183) ! ************************************************************************** !
184)
185) subroutine ImmisAuxVarCopy(auxvar,auxvar2,option)
186) !
187) ! Copies an auxiliary variable
188) !
189) ! Author: Chuan Lu
190) ! Date: 10/13/0
191) !
192)
193) use Option_module
194)
195) implicit none
196)
197) type(Immis_auxvar_elem_type) :: auxvar, auxvar2
198) type(option_type) :: option
199)
200) auxvar2%pres = auxvar%pres
201) auxvar2%temp = auxvar%temp
202) auxvar2%sat = auxvar%sat
203) auxvar2%den = auxvar%den
204) auxvar2%avgmw = auxvar%avgmw
205) auxvar2%vis = auxvar%vis
206) auxvar2%h = auxvar%h
207) auxvar2%u = auxvar%u
208) auxvar2%pc = auxvar%pc
209) ! auxvar2%kr = auxvar%kr
210) ! auxvar2%dkr_dp = auxvar%dkr_dp
211) ! auxvar2%vis = auxvar%vis
212) ! auxvar2%dvis_dp = auxvar%dvis_dp
213) auxvar2%kvr = auxvar%kvr
214) #if 0
215) auxvar2%dsat_dp = auxvar%dsat_dp
216) auxvar2%dden_dp = auxvar%dden_dp
217) auxvar2%dden_dt = auxvar%dden_dt
218) auxvar2%dkvr_dp = auxvar%dkvr_dp
219) auxvar2%dkvr_dt = auxvar%dkvr_dt
220) auxvar2%dh_dp = auxvar%dh_dp
221) auxvar2%dh_dt = auxvar%dh_dt
222) auxvar2%du_dp = auxvar%du_dp
223) auxvar2%du_dt = auxvar%du_dt
224) #endif
225) ! auxvar2%xmol = auxvar%xmol
226) ! auxvar2%diff = auxvar%diff
227)
228) end subroutine ImmisAuxVarCopy
229)
230) ! ************************************************************************** !
231)
232) subroutine ImmisAuxVarCompute_NINC(x,auxvar,saturation_function, &
233) fluid_properties,option)
234) !
235) ! ImmisAuxVarCompute_NI: Computes auxiliary variables for each grid cell
236) ! No increments
237) !
238) ! Author: Chuan Lu
239) ! Date: 10/12/08
240) !
241)
242) use Option_module
243) use EOS_Water_module
244) use Gas_EOS_module
245) use co2eos_module
246) use co2_span_wagner_module
247) use co2_span_wagner_spline_module, only: sw_prop
248) use co2_sw_module, only: co2_sw_interp
249) use Saturation_Function_module
250) use Fluid_module
251) use Mphase_pckr_module
252)
253) implicit none
254)
255) type(option_type) :: option
256) type(fluid_property_type) :: fluid_properties
257) type(saturation_function_type) :: saturation_function
258) PetscReal :: x(option%nflowdof)
259) type(Immis_auxvar_elem_type) :: auxvar
260) PetscInt :: iphase
261)
262) PetscErrorCode :: ierr
263) PetscReal :: pw,dw_kg,dw_mol,hw,sat_pressure,visl
264) PetscReal :: p,t,temp,p2,err
265) PetscReal :: henry
266) PetscReal :: dg,dddp,dddt
267) PetscReal :: fg,dfgdp,dfgdt,xphi
268) PetscReal :: eng,hg,dhdp,dhdt
269) PetscReal :: visg, dvdp, dvdt
270) PetscReal :: h(option%nphase),u(option%nphase),kr(option%nphase)
271) PetscReal :: xm_nacl, x_nacl, vphi
272) PetscReal :: aux(1)
273) PetscInt :: iflag
274)
275)
276)
277) auxvar%sat = 0.d0
278) auxvar%h = 0.d0
279) auxvar%u = 0.d0
280) auxvar%den = 0.d0
281) auxvar%avgmw = 0.d0
282) ! auxvar%vis = 0.d0
283) auxvar%pc = 0.d0
284) auxvar%kvr = 0.d0
285) ! auxvar%xmol = 0.d0
286) ! auxvar%diff = 0.d0
287) kr = 0.d0
288)
289) auxvar%pres = x(1)
290) auxvar%temp = x(2)
291) p = auxvar%pres
292) t = auxvar%temp
293)
294) if (x(3)<0.D0)x(3) = 0.D0
295) if (x(3)>1.D0)x(3) = 1.D0
296)
297) auxvar%sat(2) = x(3)
298) if (auxvar%sat(2) < 0.D0) then
299) ! print *,'tran:',iphase, x(1:3)
300) auxvar%sat(2) = 0.D0
301) endif
302) ! if (auxvar%sat(2) > 1.D0) print *,'tran:',iphase, x(1:3)
303) auxvar%sat(1) = 1.D0 - auxvar%sat(2)
304) auxvar%pc(:) = 0.D0
305) temp = 1.D-2
306)
307)
308) ! ********************* Gas phase properties ***********************
309) call EOSWaterSaturationPressure(t, sat_pressure, ierr)
310) err = 1.D0
311) p2 = p
312)
313) if (p2 >= 5.d4) then
314) if (option%co2eos == EOS_SPAN_WAGNER) then
315) ! ************ Span-Wagner EOS ********************
316) select case(option%itable)
317) case(0,1,2,4,5)
318) if (option%itable >=4) then
319) ! print *,' interp', itable
320) call co2_sw_interp(p2*1.D-6, t,dg,dddt,dddp,fg,&
321) dfgdp,dfgdt,eng,hg,dhdt,dhdp,visg,dvdt,dvdp,option%itable)
322) else
323) iflag = 1
324) call co2_span_wagner(p2*1.D-6,t+273.15D0,dg,dddt,dddp,fg,&
325) dfgdp,dfgdt,eng,hg,dhdt,dhdp,visg,dvdt,dvdp,iflag, &
326) option%itable)
327) endif
328) dg = dg / FMWCO2
329) fg = fg*1.D6
330) hg = hg*FMWCO2
331) xphi = fg/p2
332)
333) ! ************* Span-Wagner EOS with Bi-Cubic Spline interpolation ********
334) case(3)
335) call sw_prop(t,p2*1D-6,dg,hg,eng,fg)
336) call visco2(t,dg,visg)
337) dg = dg / FMWCO2
338) fg = fg*1.D6
339) hg = hg*FMWCO2
340) xphi = fg/p2
341) end select
342)
343) elseif (option%co2eos == EOS_MRK) then
344) ! MRK eos [modified version from Kerrick and Jacobs (1981) and Weir et al. (1996).]
345) call CO2(t,p2,dg,fg,xphi,hg)
346) call visco2(t,dg,visg)
347) dg = dg / FMWCO2
348) hg = hg*FMWCO2*option%scale
349) ! print *, 'translator', p2, t, dg,hg,visg
350) else
351) call printErrMsg(option,'pflow Immis ERROR: Need specify CO2 EOS')
352) endif
353) else
354) call ideal_gaseos_noderiv(p2,t,dg,hg,eng)
355) ! J/kmol -> whatever
356) hg = hg * option%scale
357) eng = eng * option%scale
358) call visco2(t,dg*FMWCO2,visg)
359) fg=p2
360) xphi = 1.D0
361) endif
362)
363) ! call Henry_duan_sun(t,p2*1D-5,henry,lngamco2, &
364) ! option%m_nacl,option%m_nacl)
365) ! henry= 1D0 / (FMWH2O*1.D-3) / (henry*1.D-5)/xphi
366)
367) pw = p
368) auxvar%den(2) = dg
369) auxvar%h(2) = hg
370) auxvar%u(2) = hg - p/dg * option%scale
371) auxvar%pc(2) = 0.D0
372)
373) ! auxvar%diff(option%nflowspec+1:option%nflowspec*2) = 2.13D-5
374) ! fluid_properties%diff_base(2)
375) ! Note: not temperature dependent yet.
376) auxvar%zco2 = auxvar%den(2)/(p/IDEAL_GAS_CONSTANT/(t+273.15D0)*1D-3)
377) !*************** Liquid phase properties **************************
378)
379) ! avgmw(1)= xmol(1)*FMWH2O + xmol(2)*FMWCO2
380) call EOSWaterDensity(t,pw,dw_kg,dw_mol,ierr)
381) call EOSWaterEnthalpy(t,pw,hw,ierr)
382) ! J/kmol -> whatever units
383) hw = hw * option%scale
384)
385) auxvar%h(1) = hw
386) auxvar%u(1) = auxvar%h(1) - pw /dw_mol*option%scale
387) ! auxvar%diff(1:option%nflowspec) = 1D-9
388) ! fluid_properties%diff_base(1)
389)
390) xm_nacl = option%m_nacl*FMWNACL
391) xm_nacl = xm_nacl /(1.D3 + xm_nacl)
392) aux(1) = xm_nacl
393) call EOSWaterDensityExt(t,p,aux,dw_kg,dw_mol,ierr)
394) call EOSWaterViscosityNaCl(t,p,xm_nacl,visl)
395)
396) !FEHM mixing ****************************
397) ! den(1) = xmol(2)*dg + xmol(1)*dw_mol
398) ! ideal mixing
399) ! den(1) = 1.D0/(xmol(2)/dg + xmol(1)/dw_mol) !*c+(1-c)*
400)
401) ! Garcia mixing **************************
402) x_nacl = option%m_nacl/(option%m_nacl + 1D3/FMWH2O)
403) ! ** xmol(1) = xh2o + xnacl
404)
405) auxvar%avgmw(1) = (1.D0 - x_nacl)*FMWH2O + x_nacl*FMWNACL
406) auxvar%avgmw(2) = FMWCO2
407) auxvar%den(1) = dw_kg/auxvar%avgmw(1)
408)
409) ! Hebach, J. Chem.Eng.Data 2004 (49),p950 ***********
410) ! den(1) = 949.7109D0 + p * (0.559684D-6 - 0.00097D-12 * p) &
411) ! + (t+273.15)*(0.883148 - 0.00228*(t+273.15))
412) ! den(1) = dw_kg + (den(1)-dw_kg)*xmol(2)/p*henry
413) ! den(1) = den(1)/avgmw(1)
414) !****************************** 2 phase S-Pc-kr relation *********************************
415) auxvar%pc = 0.D0
416)
417) if (saturation_function%hysteresis_id <= 0.1D0) then
418) call pckrNH_noderiv(auxvar%sat,auxvar%pc,kr, &
419) saturation_function, &
420) option)
421) pw=p !-pc(1)
422)
423) else
424) call pckrHY_noderiv(auxvar%sat,auxvar%hysdat,auxvar%pc,kr, &
425) saturation_function, &
426) option)
427) end if
428)
429) ! call SaturationFunctionCompute(auxvar%pres,auxvar%sat,kr, &
430) ! ds_dp,dkr_dp, &
431) ! saturation_function, &
432) ! por,perm, &
433) ! option)
434) auxvar%kvr(2) = kr(2)/visg
435) auxvar%kvr(1) = kr(1)/visl
436) auxvar%vis(2) = visg
437) auxvar%vis(1) = visl
438)
439)
440) ! print *,'immis_aux: ',auxvar%den,auxvar%avgmw,auxvar%vis,auxvar%kvr
441)
442) end subroutine ImmisAuxVarCompute_NINC
443)
444) ! ************************************************************************** !
445)
446) subroutine ImmisAuxVarCompute_WINC(x, delx, auxvar,saturation_function, &
447) fluid_properties,option)
448)
449) use Option_module
450)
451) use Saturation_Function_module
452) use Fluid_module
453)
454) implicit none
455)
456) type(option_type) :: option
457) type(fluid_property_type) :: fluid_properties
458) type(saturation_function_type) :: saturation_function
459) PetscReal :: x(option%nflowdof), xx(option%nflowdof), delx(option%nflowdof)
460) type(Immis_auxvar_elem_type) :: auxvar(1:option%nflowdof)
461) ! PetscInt :: iphase
462)
463) PetscInt :: n
464)
465) do n=1, option%nflowdof
466) xx=x; xx(n)=x(n)+ delx(n)
467) ! *** note: var_node here starts from 1 to option%flowdof ***
468) call ImmisAuxVarCompute_NINC(xx,auxvar(n),saturation_function, &
469) fluid_properties, option)
470) enddo
471)
472) end subroutine ImmisAuxVarCompute_WINC
473)
474) ! ************************************************************************** !
475)
476) subroutine ImmisAuxVarDestroy(auxvar)
477) !
478) ! AuxVarDestroy: Deallocates a richards auxiliary object
479) !
480) ! Author: Glenn Hammond
481) ! Date: 02/14/08
482) !
483)
484) implicit none
485)
486) type(Immis_auxvar_elem_type) :: auxvar
487)
488) ! if (associated(auxvar%xmol)) deallocate(auxvar%xmol)
489) ! nullify(auxvar%xmol)
490) ! if (associated(auxvar%diff))deallocate(auxvar%diff)
491) ! nullify(auxvar%diff)
492) if (associated(auxvar%pc))deallocate(auxvar%pc)
493) nullify(auxvar%pc)
494) if (associated(auxvar%sat))deallocate(auxvar%sat)
495) nullify(auxvar%sat)
496) if (associated(auxvar%u))deallocate(auxvar%u)
497) nullify(auxvar%u)
498) if (associated(auxvar%h))deallocate(auxvar%h)
499) nullify(auxvar%h)
500) if (associated(auxvar%den))deallocate(auxvar%den)
501) nullify(auxvar%den)
502) if (associated(auxvar%avgmw))deallocate(auxvar%avgmw)
503) nullify(auxvar%avgmw)
504) if (associated(auxvar%vis))deallocate(auxvar%vis)
505) nullify(auxvar%vis)
506) end subroutine ImmisAuxVarDestroy
507)
508) ! ************************************************************************** !
509)
510) subroutine ImmisAuxDestroy(aux, option)
511) !
512) ! RichardsAuxDestroy: Deallocates a richards auxiliary object
513) !
514) ! Author: Glenn Hammond
515) ! Date: 02/14/08
516) !
517)
518) use Option_module
519) implicit none
520)
521) type(Immis_type), pointer :: aux
522) type(option_type) :: option
523) PetscInt :: iaux, ielem
524)
525) if (.not.associated(aux)) return
526)
527) do iaux = 1, aux%num_aux
528) do ielem= 0, option%nflowdof
529) call ImmisAuxVarDestroy(aux%auxvars(iaux)%auxvar_elem(ielem))
530) enddo
531) deallocate(aux%auxvars)
532) enddo
533) nullify(aux%auxvars)
534)
535) do iaux = 1, aux%num_aux_bc
536) do ielem= 0, option%nflowdof
537) call ImmisAuxVarDestroy(aux%auxvars_bc(iaux)%auxvar_elem(ielem))
538) enddo
539) deallocate(aux%auxvars_bc)
540) enddo
541) nullify(aux%auxvars_bc)
542)
543) do iaux = 1, aux%num_aux_ss
544) do ielem = 0, option%nflowdof
545) call ImmisAuxVarDestroy(aux%auxvars_ss(iaux)%auxvar_elem(ielem))
546) enddo
547) deallocate(aux%auxvars_ss)
548) enddo
549) nullify(aux%auxvars_ss)
550)
551) if (associated(aux%auxvars)) deallocate(aux%auxvars)
552) nullify(aux%auxvars)
553) if (associated(aux%auxvars_bc)) deallocate(aux%auxvars_bc)
554) nullify(aux%auxvars_bc)
555) if (associated(aux%zero_rows_local)) deallocate(aux%zero_rows_local)
556) nullify(aux%zero_rows_local)
557) if (associated(aux%zero_rows_local_ghosted)) deallocate(aux%zero_rows_local_ghosted)
558) nullify(aux%zero_rows_local_ghosted)
559) if (associated(aux%res_old_AR)) deallocate(aux%res_old_AR)
560) nullify(aux%res_old_AR)
561) if (associated(aux%res_old_FL)) deallocate(aux%res_old_FL)
562) nullify(aux%res_old_FL)
563) if (associated(aux%delx)) deallocate(aux%delx)
564) nullify(aux%delx)
565) if (associated(aux%immis_parameter)) then
566) if (associated(aux%immis_parameter%dencpr)) deallocate(aux%immis_parameter%dencpr)
567) nullify(aux%immis_parameter%dencpr)
568) if (associated(aux%immis_parameter%ckwet)) deallocate(aux%immis_parameter%ckwet)
569) nullify(aux%immis_parameter%ckwet)
570) if (associated(aux%immis_parameter%ckdry)) deallocate(aux%immis_parameter%ckdry)
571) nullify(aux%immis_parameter%ckdry)
572) if (associated(aux%immis_parameter%sir)) deallocate(aux%immis_parameter%sir)
573) nullify(aux%immis_parameter%sir)
574) deallocate(aux%immis_parameter)
575) endif
576) nullify(aux%immis_parameter%dencpr)
577)
578) deallocate(aux)
579) nullify(aux)
580)
581) end subroutine ImmisAuxDestroy
582)
583)
584)
585) end module Immis_Aux_module
586)
587)