transport.F90 coverage: 37.50 %func 21.83 %block
1) module Transport_module
2)
3) use Reactive_Transport_Aux_module
4) use Global_Aux_module
5) use Material_Aux_class
6) use Matrix_Block_Aux_module
7)
8) use PFLOTRAN_Constants_module
9)
10) implicit none
11)
12) private
13)
14) #include "petsc/finclude/petscsys.h"
15)
16) #include "petsc/finclude/petscvec.h"
17) #include "petsc/finclude/petscvec.h90"
18) #include "petsc/finclude/petscmat.h"
19) #include "petsc/finclude/petscmat.h90"
20) #include "petsc/finclude/petscsnes.h"
21) #include "petsc/finclude/petscviewer.h"
22) #include "petsc/finclude/petsclog.h"
23)
24) ! Cutoff parameters
25) PetscReal, parameter :: eps = 1.D-8
26)
27) public :: TDispersion, &
28) TDispersionBC, &
29) TFlux, &
30) TFluxDerivative, &
31) TFluxCoef, &
32) TSrcSinkCoef, &
33) TFlux_CD, &
34) TFluxDerivative_CD, &
35) TFluxCoef_CD, &
36) TFluxTVD
37)
38) ! this interface is required for the pointer to procedure employed
39) ! for flux limiters below
40) interface
41) function TFluxLimiterDummy(d)
42) PetscReal :: d
43) PetscReal :: TFluxLimiterDummy
44) end function TFluxLimiterDummy
45) end interface
46)
47) public :: TFluxLimiterDummy, &
48) TFluxLimiter, &
49) TFluxLimitUpwind, &
50) TFluxLimitMinmod, &
51) TFluxLimitMC, &
52) TFluxLimitSuperBee, &
53) TFluxLimitVanLeer
54)
55) PetscInt, parameter, public :: TVD_LIMITER_UPWIND = 1
56) PetscInt, parameter, public :: TVD_LIMITER_MC = 2
57) PetscInt, parameter, public :: TVD_LIMITER_MINMOD = 3
58) PetscInt, parameter, public :: TVD_LIMITER_SUPERBEE = 4
59) PetscInt, parameter, public :: TVD_LIMITER_VAN_LEER = 5
60)
61) contains
62)
63) ! ************************************************************************** !
64)
65) subroutine TDispersion(global_auxvar_up,material_auxvar_up, &
66) cell_centered_velocity_up,dispersivity_up, &
67) global_auxvar_dn,material_auxvar_dn, &
68) cell_centered_velocity_dn,dispersivity_dn,dist, &
69) rt_parameter,option,qdarcy,dispersion)
70) !
71) ! Computes dispersion term at cell interface
72) !
73) ! Author: Glenn Hammond
74) ! Date: 02/24/10
75) !
76)
77) use Option_module
78) use Connection_module
79)
80) implicit none
81)
82) type(global_auxvar_type) :: global_auxvar_up, global_auxvar_dn
83) class(material_auxvar_type) :: material_auxvar_up, material_auxvar_dn
84) PetscReal :: dispersivity_up(3), dispersivity_dn(3)
85) PetscReal :: cell_centered_velocity_up(3), cell_centered_velocity_dn(3)
86) PetscReal :: dist(-1:3)
87) PetscReal :: qdarcy(*)
88) type(option_type) :: option
89) type(reactive_transport_param_type) :: rt_parameter
90) PetscReal :: dispersion(option%nphase)
91)
92) PetscInt :: iphase
93) PetscReal :: stp_ave_over_dist, disp_ave_over_dist
94) PetscReal :: dist_up, dist_dn
95) PetscReal :: sat_up, sat_dn
96) PetscReal :: stp_up, stp_dn
97) PetscReal :: velocity_dn(3), velocity_up(3)
98) PetscReal :: distance_gravity, upwind_weight ! both are dummy variables
99) PetscReal :: q
100) PetscReal :: Dxx_up, Dyy_up, Dzz_up, D_up
101) PetscReal :: Dxx_dn, Dyy_dn, Dzz_dn, D_dn
102)
103) #if defined(TEMP_DEPENDENT_LOGK) || defined (CHUAN_HPT)
104) PetscReal :: temp_up, temp_dn
105) PetscReal :: Ddiff_up, Ddiff_dn
106) PetscReal :: Ddiff_avg
107) PetscReal :: T_ref_inv
108) PetscReal :: weight_new
109) #endif
110)
111) dispersion(:) = 0.d0
112) iphase = 1
113) q = qdarcy(iphase)
114)
115) sat_up = global_auxvar_up%sat(iphase)
116) sat_dn = global_auxvar_dn%sat(iphase)
117)
118) call ConnectionCalculateDistances(dist,option%gravity,dist_up, &
119) dist_dn,distance_gravity, &
120) upwind_weight)
121)
122) if (dispersivity_up(2) > 0.d0 .and. dispersivity_dn(2) > 0.d0) then
123) velocity_dn = q*dist(1:3) + (1.d0-dist(1:3))*cell_centered_velocity_dn
124) velocity_up = q*dist(1:3) + (1.d0-dist(1:3))*cell_centered_velocity_up
125) Dxx_up = dispersivity_up(1)*dabs(velocity_up(X_DIRECTION)) + &
126) dispersivity_up(2)*dabs(velocity_up(Y_DIRECTION)) + &
127) dispersivity_up(3)*dabs(velocity_up(Z_DIRECTION))
128) Dxx_dn = dispersivity_dn(1)*dabs(velocity_dn(X_DIRECTION)) + &
129) dispersivity_dn(2)*dabs(velocity_dn(Y_DIRECTION)) + &
130) dispersivity_dn(3)*dabs(velocity_dn(Z_DIRECTION))
131) Dyy_up = dispersivity_up(2)*dabs(velocity_up(X_DIRECTION)) + &
132) dispersivity_up(1)*dabs(velocity_up(Y_DIRECTION)) + &
133) dispersivity_up(3)*dabs(velocity_up(Z_DIRECTION))
134) Dyy_dn = dispersivity_dn(2)*dabs(velocity_dn(X_DIRECTION)) + &
135) dispersivity_dn(1)*dabs(velocity_dn(Y_DIRECTION)) + &
136) dispersivity_dn(3)*dabs(velocity_dn(Z_DIRECTION))
137) Dzz_up = dispersivity_up(3)*dabs(velocity_up(X_DIRECTION)) + &
138) dispersivity_up(3)*dabs(velocity_up(Y_DIRECTION)) + &
139) dispersivity_up(1)*dabs(velocity_up(Z_DIRECTION))
140) Dzz_dn = dispersivity_dn(3)*dabs(velocity_dn(X_DIRECTION)) + &
141) dispersivity_dn(3)*dabs(velocity_dn(Y_DIRECTION)) + &
142) dispersivity_dn(1)*dabs(velocity_dn(Z_DIRECTION))
143) D_up = max(dist(1)*Dxx_up+dist(2)*Dyy_up+dist(3)*Dzz_up,1.d-40)
144) D_dn = max(dist(1)*Dxx_dn+dist(2)*Dyy_dn+dist(3)*Dzz_dn,1.d-40)
145) dispersion(iphase) = D_up*D_dn/(D_up*dist_dn+D_dn*dist_up)
146) else
147)
148) ! Weighted harmonic mean of dispersivity divided by distance
149) if (dispersivity_up(1) > eps .and. dispersivity_dn(1) > eps) then
150) disp_ave_over_dist = (dispersivity_up(1)*dispersivity_dn(1)) / &
151) (dispersivity_up(1)*dist_dn+dispersivity_dn(1)*dist_up)
152) else
153) ! still need to set this as it is used later in CO2 below.
154) disp_ave_over_dist = 0.d0
155) endif
156) dispersion(iphase) = disp_ave_over_dist*dabs(q)
157) endif
158)
159) if (sat_up > eps .and. sat_dn > eps) then
160) stp_up = sat_up*material_auxvar_up%tortuosity*material_auxvar_up%porosity
161) stp_dn = sat_dn*material_auxvar_dn%tortuosity*material_auxvar_dn%porosity
162) ! units = (m^3 water/m^3 por)*(m^3 por/m^3 bulk)/(m bulk) = m^3 water/m^4 bulk
163) stp_ave_over_dist = (stp_up*stp_dn)/(stp_up*dist_dn+stp_dn*dist_up)
164) ! need to account for multiple phases
165) ! units = (m^3 water/m^4 bulk)*(m^2 bulk/sec) = m^3 water/m^2 bulk/sec
166) dispersion(iphase) = dispersion(iphase) + &
167) stp_ave_over_dist*rt_parameter%diffusion_coefficient(iphase)
168)
169) ! Add the effect of temperature on diffusivity, Satish Karra, LANL, 10/29/2011
170)
171) #if defined(TEMP_DEPENDENT_LOGK) || defined (CHUAN_HPT)
172) T_ref_inv = 1.d0/(25.d0 + 273.15d0)
173) temp_up = global_auxvar_up%temp ! getting data from global to local variables
174) temp_dn = global_auxvar_dn%temp
175) Ddiff_up = rt_parameter%diffusion_coefficient(iphase)* &
176) exp(rt_parameter%diffusion_activation_energy(iphase) &
177) /IDEAL_GAS_CONSTANT*(T_ref_inv - 1.d0/(temp_up + 273.15d0)))
178) Ddiff_dn = rt_parameter%diffusion_coefficient(iphase)* &
179) exp(rt_parameter%diffusion_activation_energy(iphase) &
180) /IDEAL_GAS_CONSTANT*(T_ref_inv - 1.d0/(temp_dn + 273.15d0)))
181) weight_new = (stp_up*Ddiff_up*stp_dn*Ddiff_dn)/ &
182) (stp_up*Ddiff_up*dist_dn + stp_dn*Ddiff_dn*dist_up)
183) dispersion(iphase) = dispersion(iphase) + weight_new - &
184) stp_ave_over_dist*rt_parameter%diffusion_coefficient(iphase)
185) #endif
186) endif
187)
188) ! CO2-specific
189) ! Add in multiphase, clu 12/29/08
190) if (option%iflowmode == MPH_MODE .or. option%iflowmode == IMS_MODE &
191) .or. option%iflowmode == FLASH2_MODE) then
192) do
193) iphase = iphase +1
194) if (iphase > option%nphase) exit
195) ! super critical CO2 phase have the index 2: need implementation
196) q = qdarcy(iphase)
197) sat_up = global_auxvar_up%sat(iphase)
198) sat_dn = global_auxvar_dn%sat(iphase)
199) if (sat_up > eps .and. sat_dn > eps) then
200) stp_up = sat_up*material_auxvar_up%tortuosity*material_auxvar_up%porosity
201) stp_dn = sat_dn*material_auxvar_dn%tortuosity*material_auxvar_dn%porosity
202) ! units = (m^3 water/m^3 por)*(m^3 por/m^3 bulk)/(m bulk) = m^3 water/m^4 bulk
203) stp_ave_over_dist = (stp_up*stp_dn)/(stp_up*dist_dn+stp_dn*dist_up)
204) ! need to account for multiple phases
205) ! units = (m^3 water/m^4 bulk)*(m^2 bulk/sec) = m^3 water/m^2 bulk/sec
206) if (iphase == 2) then
207) dispersion(iphase) = &
208) disp_ave_over_dist*dabs(q) + &
209) stp_ave_over_dist*rt_parameter%diffusion_coefficient(iphase)
210)
211) ! Add the effect of temperature on diffusivity, Satish Karra, LANL, 11/1/2011
212) #if defined(TEMP_DEPENDENT_LOGK) || defined (CHUAN_HPT)
213) T_ref_inv = 1.d0/(25.d0 + 273.15d0)
214) temp_up = global_auxvar_up%temp
215) temp_dn = global_auxvar_dn%temp
216) Ddiff_up = rt_parameter%diffusion_coefficient(iphase)* &
217) exp(rt_parameter%diffusion_activation_energy(iphase) &
218) /IDEAL_GAS_CONSTANT*(T_ref_inv - 1.d0/(temp_up + 273.15d0)))
219) Ddiff_dn = rt_parameter%diffusion_coefficient(iphase)* &
220) exp(rt_parameter%diffusion_activation_energy(iphase) &
221) /IDEAL_GAS_CONSTANT*(T_ref_inv - 1.d0/(temp_dn + 273.15d0)))
222) weight_new = (stp_up*Ddiff_up*stp_dn*Ddiff_dn)/ &
223) (stp_up*Ddiff_up*dist_dn + stp_dn*Ddiff_dn*dist_up)
224) dispersion(iphase) = dispersion(iphase) + weight_new - &
225) stp_ave_over_dist* &
226) rt_parameter%diffusion_coefficient(iphase)
227) #endif
228) endif
229) endif
230) enddo
231) endif
232)
233) end subroutine TDispersion
234)
235) ! ************************************************************************** !
236)
237) subroutine TDispersionBC(ibndtype, &
238) global_auxvar_up, &
239) global_auxvar_dn, &
240) material_auxvar_dn, &
241) cell_centered_velocity_dn, &
242) dispersivity_dn,dist_dn, &
243) rt_parameter,option,qdarcy,dispersion)
244) !
245) ! Computes dispersion term at cell boundary interface
246) !
247) ! Author: Glenn Hammond
248) ! Date: 02/15/08
249) !
250)
251) use Option_module
252)
253) implicit none
254)
255) PetscInt :: ibndtype
256) type(global_auxvar_type) :: global_auxvar_up, global_auxvar_dn
257) class(material_auxvar_type) :: material_auxvar_dn
258) PetscReal :: cell_centered_velocity_dn(3)
259) PetscReal :: dispersivity_dn(3), dist_dn(-1:3)
260) PetscReal :: qdarcy(1)
261) type(reactive_transport_param_type) :: rt_parameter
262) type(option_type) :: option
263) PetscReal :: dispersion(option%nphase)
264)
265) PetscInt :: icomp
266) PetscInt :: iphase
267) PetscReal :: stp_ave_over_dist
268) PetscReal :: q
269) PetscReal :: sat_up
270) PetscReal :: temp_dispersion(option%nphase)
271) PetscReal :: Dxx_dn, Dyy_dn, Dzz_dn, D_dn
272) PetscReal :: velocity_dn(3)
273)
274) #if defined(TEMP_DEPENDENT_LOGK) || defined (CHUAN_HPT)
275) PetscReal :: temp_up ! variable to store temperature at the boundary
276) PetscReal :: T_ref_inv
277) #endif
278)
279) temp_dispersion(:) = 0.d0
280) dispersion(:) = 0.d0
281) iphase = 1
282) q = qdarcy(iphase)
283)
284) ! we use upwind saturation as that is the saturation at the boundary face
285) sat_up = global_auxvar_up%sat(iphase)
286)
287) if (dispersivity_dn(2) > 0.d0) then
288) velocity_dn = q*dist_dn(1:3) + (1.d0-dist_dn(1:3))*cell_centered_velocity_dn
289) Dxx_dn = dispersivity_dn(1)*dabs(velocity_dn(X_DIRECTION)) + &
290) dispersivity_dn(2)*dabs(velocity_dn(Y_DIRECTION)) + &
291) dispersivity_dn(3)*dabs(velocity_dn(Z_DIRECTION))
292) Dyy_dn = dispersivity_dn(2)*dabs(velocity_dn(X_DIRECTION)) + &
293) dispersivity_dn(1)*dabs(velocity_dn(Y_DIRECTION)) + &
294) dispersivity_dn(3)*dabs(velocity_dn(Z_DIRECTION))
295) Dzz_dn = dispersivity_dn(3)*dabs(velocity_dn(X_DIRECTION)) + &
296) dispersivity_dn(3)*dabs(velocity_dn(Y_DIRECTION)) + &
297) dispersivity_dn(1)*dabs(velocity_dn(Z_DIRECTION))
298) D_dn = max(Dxx_dn+Dyy_dn+Dzz_dn,1.d-40)
299) temp_dispersion(iphase) = D_dn
300) else
301) temp_dispersion(iphase) = dispersivity_dn(1)*dabs(q)/dist_dn(0)
302) endif
303)
304) select case(ibndtype)
305) case(DIRICHLET_BC)
306) ! units = (m^3 water/m^3 por)*(m^3 por/m^3 bulk)/(m bulk) =
307) ! m^3 water/m^4 bulk
308)
309) stp_ave_over_dist = (material_auxvar_dn%tortuosity* &
310) material_auxvar_dn%porosity*sat_up) / dist_dn(0)
311)
312) ! need to account for multiple phases
313) ! units = (m^3 water/m^4 bulk)*(m^2 bulk/sec) = m^3 water/m^2 bulk/sec
314) dispersion(iphase) = temp_dispersion(iphase) + &
315) stp_ave_over_dist* &
316) rt_parameter%diffusion_coefficient(iphase)
317)
318) #if defined(TEMP_DEPENDENT_LOGK) || defined (CHUAN_HPT)
319) T_ref_inv = 1.d0/(25.d0 + 273.15d0)
320) temp_up = global_auxvar_up%temp
321) dispersion(iphase) = dispersion(iphase) + &
322) stp_ave_over_dist*rt_parameter%diffusion_coefficient(iphase)* &
323) (exp(rt_parameter%diffusion_activation_energy(iphase)/ &
324) IDEAL_GAS_CONSTANT*(T_ref_inv-1.d0/(temp_up + 273.15d0))) - 1.d0)
325) #endif
326)
327) case(DIRICHLET_ZERO_GRADIENT_BC)
328) if (q >= 0.d0) then
329) ! same as dirichlet above
330) ! units = (m^3 water/m^3 por)*(m^3 por/m^3 bulk)/(m bulk) =
331) ! m^3 water/m^4 bulk
332)
333) stp_ave_over_dist = (material_auxvar_dn%tortuosity* &
334) material_auxvar_dn%porosity*sat_up) / dist_dn(0)
335)
336) ! need to account for multiple phases
337) ! units = (m^3 water/m^4 bulk)*(m^2 bulk/sec) = m^3 water/m^2 bulk/sec
338) dispersion(iphase) = temp_dispersion(iphase) + &
339) stp_ave_over_dist* &
340) rt_parameter%diffusion_coefficient(iphase)
341)
342) #if defined(TEMP_DEPENDENT_LOGK) || defined (CHUAN_HPT)
343) T_ref_inv = 1.d0/(25.d0 + 273.15d0)
344) temp_up = global_auxvar_up%temp
345) dispersion(iphase) = dispersion(iphase) + &
346) stp_ave_over_dist*rt_parameter%diffusion_coefficient(iphase)* &
347) (exp(rt_parameter%diffusion_activation_energy(iphase)/ &
348) IDEAL_GAS_CONSTANT*(T_ref_inv-1.d0/(temp_up + 273.15d0))) - 1.d0)
349) #endif
350) endif
351) case(CONCENTRATION_SS,NEUMANN_BC,ZERO_GRADIENT_BC)
352) end select
353)
354)
355) ! CO2-specific
356) ! Add in multiphase, clu 12/29/08
357) if (option%iflowmode == MPH_MODE .or. option%iflowmode == IMS_MODE &
358) .or. option%iflowmode == FLASH2_MODE) then
359) do
360) iphase = iphase + 1
361) if (iphase > option%nphase) exit
362) q = qdarcy(iphase)
363) sat_up = global_auxvar_up%sat(iphase)
364)
365) select case(ibndtype)
366) case(DIRICHLET_BC)
367) ! units = (m^3 water/m^3 por)*(m^3 por/m^3 bulk)/(m bulk) =
368) ! m^3 water/m^4 bulk
369)
370) stp_ave_over_dist = (material_auxvar_dn%tortuosity* &
371) material_auxvar_dn%porosity*sat_up) / dist_dn(0)
372)
373) ! need to account for multiple phases
374) ! units = (m^3 water/m^4 bulk)*(m^2 bulk/sec) =
375) ! m^3 water/m^2 bulk/sec
376) if ( iphase == 2) then
377) dispersion(iphase) = dispersivity_dn(1)*dabs(q)/dist_dn(0) + &
378) stp_ave_over_dist * &
379) rt_parameter%diffusion_coefficient(iphase)
380)
381) #if defined(TEMP_DEPENDENT_LOGK) || defined (CHUAN_HPT)
382) T_ref_inv = 1.d0/(25.d0 + 273.15d0)
383) temp_up = global_auxvar_up%temp
384) dispersion(iphase) = dispersion(iphase) + &
385) stp_ave_over_dist*rt_parameter%diffusion_coefficient(iphase)* &
386) (exp(rt_parameter%diffusion_activation_energy(iphase)/ &
387) IDEAL_GAS_CONSTANT*(T_ref_inv-1.d0/(temp_up + 273.15d0))) - 1.d0)
388) #endif
389) endif
390)
391) case(DIRICHLET_ZERO_GRADIENT_BC)
392) if (q >= 0.d0) then
393) ! same as dirichlet above
394) stp_ave_over_dist = (material_auxvar_dn%tortuosity* &
395) material_auxvar_dn%porosity*sat_up) / dist_dn(0)
396) if (iphase == 2) then
397) dispersion(iphase) = dispersivity_dn(1)*dabs(q)/dist_dn(0) + &
398) stp_ave_over_dist * &
399) rt_parameter%diffusion_coefficient(iphase)
400) #if defined(TEMP_DEPENDENT_LOGK) || defined (CHUAN_HPT)
401) T_ref_inv = 1.d0/(25.d0 + 273.15d0)
402) temp_up = global_auxvar_up%temp
403) dispersion(iphase) = dispersion(iphase) + &
404) stp_ave_over_dist*rt_parameter%diffusion_coefficient(iphase)* &
405) (exp(rt_parameter%diffusion_activation_energy(iphase)/ &
406) IDEAL_GAS_CONSTANT*(T_ref_inv-1.d0/(temp_up + 273.15d0))) - 1.d0)
407) #endif
408) endif
409) endif
410) case(CONCENTRATION_SS,NEUMANN_BC,ZERO_GRADIENT_BC)
411) end select
412) enddo
413) endif
414)
415) end subroutine TDispersionBC
416)
417) ! ************************************************************************** !
418)
419) subroutine TFlux(rt_parameter, &
420) rt_auxvar_up,global_auxvar_up, &
421) rt_auxvar_dn,global_auxvar_dn, &
422) coef_up,coef_dn,option,Res)
423) !
424) ! Computes flux term in residual function
425) !
426) ! Author: Glenn Hammond
427) ! Date: 02/15/08
428) !
429)
430) use Option_module
431)
432) implicit none
433)
434) type(reactive_transport_param_type) :: rt_parameter
435) type(reactive_transport_auxvar_type) :: rt_auxvar_up, rt_auxvar_dn
436) type(global_auxvar_type) :: global_auxvar_up, global_auxvar_dn
437) PetscReal :: coef_up(*), coef_dn(*)
438) type(option_type) :: option
439) PetscReal :: Res(rt_parameter%ncomp)
440)
441) PetscInt :: iphase
442) PetscInt :: idof
443) PetscInt :: ndof
444) PetscInt :: istart
445) PetscInt :: iend
446) PetscInt :: icollcomp
447) PetscInt :: icoll
448) PetscInt :: iaqcomp
449)
450) iphase = 1
451) ndof = rt_parameter%naqcomp
452)
453) Res = 0.d0
454)
455) ! units = (L water/sec)*(mol/L) = mol/s
456) ! total = mol/L water
457) Res(1:ndof) = coef_up(iphase)*rt_auxvar_up%total(1:ndof,iphase) + &
458) coef_dn(iphase)*rt_auxvar_dn%total(1:ndof,iphase)
459)
460) if (rt_parameter%ncoll > 0) then
461) do icoll = 1, rt_parameter%ncoll
462) idof = rt_parameter%offset_colloid + icoll
463) Res(idof) = &
464) ! conc_mob = mol/L water
465) coef_up(iphase)*rt_auxvar_up%colloid%conc_mob(icoll)+ &
466) coef_dn(iphase)*rt_auxvar_dn%colloid%conc_mob(icoll)
467) enddo
468) endif
469) if (rt_parameter%ncollcomp > 0) then
470) do icollcomp = 1, rt_parameter%ncollcomp
471) iaqcomp = rt_parameter%coll_spec_to_pri_spec(icollcomp)
472) ! total_eq_mob = mol/L water
473) Res(iaqcomp) = Res(iaqcomp) + &
474) coef_up(iphase)*rt_auxvar_up%colloid%total_eq_mob(icollcomp) + &
475) coef_dn(iphase)*rt_auxvar_dn%colloid%total_eq_mob(icollcomp)
476) enddo
477) endif
478)
479) ! CO2-specific
480) ! Add in multiphase, clu 12/29/08
481) if (option%iflowmode == MPH_MODE .or. option%iflowmode == IMS_MODE &
482) .or. option%iflowmode == FLASH2_MODE) then
483) do
484) iphase = iphase +1
485) if (iphase > option%nphase) exit
486) ! super critical CO2 phase have the index 2: need implementation
487)
488) ! units = (L water/sec)*(mol/L) = mol/s
489) Res(1:ndof) = Res (1:ndof) + &
490) coef_up(iphase)*rt_auxvar_up%total(1:ndof,iphase) + &
491) coef_dn(iphase)*rt_auxvar_dn%total(1:ndof,iphase)
492) enddo
493) endif
494)
495) end subroutine TFlux
496)
497) ! ************************************************************************** !
498)
499) subroutine TFlux_CD(rt_parameter, &
500) rt_auxvar_up,global_auxvar_up, &
501) rt_auxvar_dn,global_auxvar_dn, &
502) coef_11,coef_12,coef_21,coef_22,option,Res_1,Res_2)
503) !
504) ! TFlux: Computes flux term in residual function
505) !
506) ! Author: Glenn Hammond
507) ! Date: 02/15/08
508) !
509)
510) use Option_module
511)
512) implicit none
513)
514) type(reactive_transport_param_type) :: rt_parameter
515) type(reactive_transport_auxvar_type) :: rt_auxvar_up, rt_auxvar_dn
516) type(global_auxvar_type) :: global_auxvar_up, global_auxvar_dn
517) PetscReal :: coef_11(*), coef_12(*), coef_21(*), coef_22(*)
518) type(option_type) :: option
519) PetscReal :: Res_1(rt_parameter%ncomp)
520) PetscReal :: Res_2(rt_parameter%ncomp)
521)
522) PetscInt :: iphase
523) PetscInt :: idof
524) PetscInt :: ndof
525) PetscInt :: istart
526) PetscInt :: iend
527) PetscInt :: icollcomp
528) PetscInt :: icoll
529) PetscInt :: iaqcomp
530)
531) iphase = 1
532) ndof = rt_parameter%naqcomp
533)
534) Res_1 = 0.d0
535) Res_2 = 0.d0
536)
537) ! units = (L water/sec)*(mol/L) = mol/s
538) ! total = mol/L water
539) Res_1(1:ndof) = coef_11(iphase)*rt_auxvar_up%total(1:ndof,iphase) + &
540) coef_12(iphase)*rt_auxvar_dn%total(1:ndof,iphase)
541) Res_2(1:ndof) = coef_21(iphase)*rt_auxvar_up%total(1:ndof,iphase) + &
542) coef_22(iphase)*rt_auxvar_dn%total(1:ndof,iphase)
543)
544) if (rt_parameter%ncoll > 0) then
545) do icoll = 1, rt_parameter%ncoll
546) idof = rt_parameter%offset_colloid + icoll
547) ! conc_mob = mol/L water
548) Res_1(idof) = coef_11(iphase)*rt_auxvar_up%colloid%conc_mob(icoll)+ &
549) coef_12(iphase)*rt_auxvar_dn%colloid%conc_mob(icoll)
550) Res_2(idof) = coef_21(iphase)*rt_auxvar_up%colloid%conc_mob(icoll)+ &
551) coef_22(iphase)*rt_auxvar_dn%colloid%conc_mob(icoll)
552) enddo
553) endif
554) if (rt_parameter%ncollcomp > 0) then
555) do icollcomp = 1, rt_parameter%ncollcomp
556) iaqcomp = rt_parameter%coll_spec_to_pri_spec(icollcomp)
557) ! total_eq_mob = mol/L water
558) Res_1(iaqcomp) = Res_1(iaqcomp) + &
559) coef_11(iphase)*rt_auxvar_up%colloid%total_eq_mob(icollcomp) + &
560) coef_12(iphase)*rt_auxvar_dn%colloid%total_eq_mob(icollcomp)
561) Res_2(iaqcomp) = Res_2(iaqcomp) + &
562) coef_21(iphase)*rt_auxvar_up%colloid%total_eq_mob(icollcomp) + &
563) coef_22(iphase)*rt_auxvar_dn%colloid%total_eq_mob(icollcomp)
564) enddo
565) endif
566)
567) ! CO2-specific
568) ! Add in multiphase, clu 12/29/08
569) if (option%iflowmode == MPH_MODE .or. option%iflowmode == IMS_MODE &
570) .or. option%iflowmode == FLASH2_MODE) then
571) do
572) iphase = iphase +1
573) if (iphase > option%nphase) exit
574) ! super critical CO2 phase have the index 2: need implementation
575)
576) ! units = (L water/sec)*(mol/L) = mol/s
577) Res_1(1:ndof) = Res_1(1:ndof) + &
578) coef_11(iphase)*rt_auxvar_up%total(1:ndof,iphase) + &
579) coef_12(iphase)*rt_auxvar_dn%total(1:ndof,iphase)
580) Res_2(1:ndof) = Res_2(1:ndof) + &
581) coef_21(iphase)*rt_auxvar_up%total(1:ndof,iphase) + &
582) coef_22(iphase)*rt_auxvar_dn%total(1:ndof,iphase)
583) enddo
584) endif
585)
586) end subroutine TFlux_CD
587)
588) ! ************************************************************************** !
589)
590) subroutine TFluxDerivative(rt_parameter, &
591) rt_auxvar_up,global_auxvar_up, &
592) rt_auxvar_dn,global_auxvar_dn, &
593) coef_up,coef_dn,option,J_up,J_dn)
594) !
595) ! Computes derivatives of flux term in residual function
596) !
597) ! Author: Glenn Hammond
598) ! Date: 02/15/08
599) !
600)
601) use Option_module
602)
603) implicit none
604)
605) type(reactive_transport_param_type) :: rt_parameter
606) type(reactive_transport_auxvar_type) :: rt_auxvar_up, rt_auxvar_dn
607) type(global_auxvar_type) :: global_auxvar_up, global_auxvar_dn
608) PetscReal :: coef_up(*), coef_dn(*)
609) type(option_type) :: option
610) PetscReal :: J_up(rt_parameter%ncomp,rt_parameter%ncomp), &
611) J_dn(rt_parameter%ncomp,rt_parameter%ncomp)
612)
613) PetscInt :: iphase
614) PetscInt :: icomp
615) PetscInt :: icoll
616) PetscInt :: idof
617) PetscInt :: istart
618) PetscInt :: iendaq
619)
620) iphase = 1
621)
622) ! units = (m^3 water/sec)*(kg water/L water)*(1000L water/m^3 water)
623) ! = kg water/sec
624) istart = 1
625) iendaq = rt_parameter%naqcomp
626) J_up = 0.d0
627) J_dn = 0.d0
628) if (associated(rt_auxvar_dn%aqueous%dtotal)) then
629) J_up(istart:iendaq,istart:iendaq) = &
630) rt_auxvar_up%aqueous%dtotal(:,:,iphase)*coef_up(iphase)
631) J_dn(istart:iendaq,istart:iendaq) = &
632) rt_auxvar_dn%aqueous%dtotal(:,:,iphase)*coef_dn(iphase)
633) else
634) do icomp = istart, iendaq
635) J_up(icomp,icomp) = coef_up(iphase)* &
636) global_auxvar_up%den_kg(iphase)*1.d-3
637) J_dn(icomp,icomp) = coef_dn(iphase)* &
638) global_auxvar_dn%den_kg(iphase)*1.d-3
639) enddo
640) endif
641)
642) if (rt_parameter%ncoll > 0) then
643) do icoll = 1, rt_parameter%ncoll
644) idof = rt_parameter%offset_colloid + icoll
645) J_up(idof,idof) = coef_up(iphase)*global_auxvar_up%den_kg(iphase)*1.d-3
646) J_dn(idof,idof) = coef_dn(iphase)*global_auxvar_dn%den_kg(iphase)*1.d-3
647) enddo
648) endif
649) if (rt_parameter%ncollcomp > 0) then
650) ! dRj_dCj - mobile
651) ! istart & iend same as above
652) J_up(istart:iendaq,istart:iendaq) = J_up(istart:iendaq,istart:iendaq) + &
653) rt_auxvar_up%colloid%dRj_dCj%dtotal(:,:,iphase)*coef_up(iphase)
654) J_dn(istart:iendaq,istart:iendaq) = J_dn(istart:iendaq,istart:iendaq) + &
655) rt_auxvar_dn%colloid%dRj_dCj%dtotal(:,:,iphase)*coef_dn(iphase)
656) ! need the below
657) ! dRj_dSic
658) ! dRic_dSic
659) ! dRic_dCj
660) endif
661)
662) ! CO2-specific
663) ! Add in multiphase, clu 12/29/08
664) if (option%iflowmode == MPH_MODE .or. option%iflowmode == IMS_MODE &
665) .or. option%iflowmode == FLASH2_MODE) then
666) do
667) iphase = iphase + 1
668) if (iphase > option%nphase) exit
669) ! super critical CO2 phase
670)
671) ! units = (m^3 water/sec)*(kg water/L water)*(1000L water/m^3 water) = kg water/sec
672) if (associated(rt_auxvar_dn%aqueous%dtotal)) then
673) J_up(istart:iendaq,istart:iendaq) = J_up(istart:iendaq,istart:iendaq) + &
674) rt_auxvar_up%aqueous%dtotal(:,:,iphase)*coef_up(iphase)
675) J_dn(istart:iendaq,istart:iendaq) = J_dn(istart:iendaq,istart:iendaq) + &
676) rt_auxvar_dn%aqueous%dtotal(:,:,iphase)*coef_dn(iphase)
677) else
678) print *,'Dtotal needed for SC problem. STOP'
679) stop
680) ! J_up = 0.d0
681) ! J_dn = 0.d0
682) ! do icomp = 1, ndof
683) ! J_up(icomp,icomp) = J_up(icomp,icomp) + coef_up*global_auxvar_up%den_kg(iphase)
684) ! J_dn(icomp,icomp) = J_dn(icomp,icomp) + coef_dn*global_auxvar_dn%den_kg(iphase)
685) ! enddo
686) endif
687) enddo
688) endif
689)
690) end subroutine TFluxDerivative
691)
692) ! ************************************************************************** !
693)
694) subroutine TFluxDerivative_CD(rt_parameter, &
695) rt_auxvar_up,global_auxvar_up, &
696) rt_auxvar_dn,global_auxvar_dn, &
697) coef_11,coef_12,coef_21,coef_22,option, &
698) J_11,J_12,J_21,J_22)
699) !
700) ! TFluxDerivative: Computes derivatives of flux term in residual function
701) !
702) ! Author: Glenn Hammond
703) ! Date: 02/15/08
704) !
705)
706) use Option_module
707)
708) implicit none
709)
710) type(reactive_transport_param_type) :: rt_parameter
711) type(reactive_transport_auxvar_type) :: rt_auxvar_up, rt_auxvar_dn
712) type(global_auxvar_type) :: global_auxvar_up, global_auxvar_dn
713) PetscReal :: coef_11(*), coef_12(*), coef_21(*), coef_22(*)
714) type(option_type) :: option
715) PetscReal :: J_11(rt_parameter%ncomp,rt_parameter%ncomp), &
716) J_12(rt_parameter%ncomp,rt_parameter%ncomp), &
717) J_21(rt_parameter%ncomp,rt_parameter%ncomp), &
718) J_22(rt_parameter%ncomp,rt_parameter%ncomp)
719)
720) PetscInt :: iphase
721) PetscInt :: icomp
722) PetscInt :: icoll
723) PetscInt :: idof
724) PetscInt :: istart
725) PetscInt :: iendaq
726)
727) iphase = 1
728)
729) ! units = (m^3 water/sec)*(kg water/L water)*(1000L water/m^3 water) = kg water/sec
730) istart = 1
731) iendaq = rt_parameter%naqcomp
732) J_11 = 0.d0
733) J_12 = 0.d0
734) J_21 = 0.d0
735) J_22 = 0.d0
736) if (associated(rt_auxvar_dn%aqueous%dtotal)) then
737) J_11(istart:iendaq,istart:iendaq) = rt_auxvar_up%aqueous%dtotal(:,:,iphase)*coef_11(iphase)
738) J_12(istart:iendaq,istart:iendaq) = rt_auxvar_dn%aqueous%dtotal(:,:,iphase)*coef_12(iphase)
739) J_21(istart:iendaq,istart:iendaq) = rt_auxvar_up%aqueous%dtotal(:,:,iphase)*coef_21(iphase)
740) J_22(istart:iendaq,istart:iendaq) = rt_auxvar_dn%aqueous%dtotal(:,:,iphase)*coef_22(iphase)
741) else
742) do icomp = istart, iendaq
743) J_11(icomp,icomp) = coef_11(iphase)*global_auxvar_up%den_kg(iphase)*1.d-3
744) J_12(icomp,icomp) = coef_12(iphase)*global_auxvar_dn%den_kg(iphase)*1.d-3
745) J_21(icomp,icomp) = coef_21(iphase)*global_auxvar_up%den_kg(iphase)*1.d-3
746) J_22(icomp,icomp) = coef_22(iphase)*global_auxvar_dn%den_kg(iphase)*1.d-3
747) enddo
748) endif
749)
750) if (rt_parameter%ncoll > 0) then
751) do icoll = 1, rt_parameter%ncoll
752) idof = rt_parameter%offset_colloid + icoll
753) J_11(idof,idof) = coef_11(iphase)*global_auxvar_up%den_kg(iphase)*1.d-3
754) J_12(idof,idof) = coef_12(iphase)*global_auxvar_dn%den_kg(iphase)*1.d-3
755) J_21(idof,idof) = coef_21(iphase)*global_auxvar_up%den_kg(iphase)*1.d-3
756) J_22(idof,idof) = coef_22(iphase)*global_auxvar_dn%den_kg(iphase)*1.d-3
757) enddo
758) endif
759) if (rt_parameter%ncollcomp > 0) then
760) ! dRj_dCj - mobile
761) ! istart & iend same as above
762) J_11(istart:iendaq,istart:iendaq) = J_11(istart:iendaq,istart:iendaq) + &
763) rt_auxvar_up%colloid%dRj_dCj%dtotal(:,:,iphase)*coef_11(iphase)
764) J_12(istart:iendaq,istart:iendaq) = J_12(istart:iendaq,istart:iendaq) + &
765) rt_auxvar_dn%colloid%dRj_dCj%dtotal(:,:,iphase)*coef_12(iphase)
766) J_21(istart:iendaq,istart:iendaq) = J_21(istart:iendaq,istart:iendaq) + &
767) rt_auxvar_up%colloid%dRj_dCj%dtotal(:,:,iphase)*coef_21(iphase)
768) J_22(istart:iendaq,istart:iendaq) = J_22(istart:iendaq,istart:iendaq) + &
769) rt_auxvar_dn%colloid%dRj_dCj%dtotal(:,:,iphase)*coef_22(iphase)
770) ! need the below
771) ! dRj_dSic
772) ! dRic_dSic
773) ! dRic_dCj
774) endif
775)
776) ! CO2-specific
777) ! Add in multiphase, clu 12/29/08
778) if (option%iflowmode == MPH_MODE .or. option%iflowmode == IMS_MODE &
779) .or. option%iflowmode == FLASH2_MODE) then
780) do
781) iphase = iphase + 1
782) if (iphase > option%nphase) exit
783) ! super critical CO2 phase
784)
785) ! units = (m^3 water/sec)*(kg water/L water)*(1000L water/m^3 water) = kg water/sec
786) if (associated(rt_auxvar_dn%aqueous%dtotal)) then
787) J_11(istart:iendaq,istart:iendaq) = J_11(istart:iendaq,istart:iendaq) + &
788) rt_auxvar_up%aqueous%dtotal(:,:,iphase)*coef_11(iphase)
789) J_12(istart:iendaq,istart:iendaq) = J_12(istart:iendaq,istart:iendaq) + &
790) rt_auxvar_dn%aqueous%dtotal(:,:,iphase)*coef_12(iphase)
791) J_21(istart:iendaq,istart:iendaq) = J_21(istart:iendaq,istart:iendaq) + &
792) rt_auxvar_up%aqueous%dtotal(:,:,iphase)*coef_21(iphase)
793) J_22(istart:iendaq,istart:iendaq) = J_22(istart:iendaq,istart:iendaq) + &
794) rt_auxvar_dn%aqueous%dtotal(:,:,iphase)*coef_22(iphase)
795) else
796) print *,'Dtotal needed for SC problem. STOP'
797) stop
798) ! J_up = 0.d0
799) ! J_dn = 0.d0
800) ! do icomp = 1, ndof
801) ! J_up(icomp,icomp) = J_up(icomp,icomp) + coef_up*global_auxvar_up%den_kg(iphase)
802) ! J_dn(icomp,icomp) = J_dn(icomp,icomp) + coef_dn*global_auxvar_dn%den_kg(iphase)
803) ! enddo
804) endif
805) enddo
806) endif
807)
808) end subroutine TFluxDerivative_CD
809)
810) ! ************************************************************************** !
811)
812) subroutine TFluxCoef(option,area,velocity,diffusion,fraction_upwind,T_up,T_dn)
813) !
814) ! Computes flux coefficients for transport matrix
815) !
816) ! Author: Glenn Hammond
817) ! Date: 02/22/10
818) !
819)
820) use Option_module
821)
822) implicit none
823)
824) type(option_type) :: option
825) PetscReal :: area
826) PetscReal :: velocity(*)
827) PetscReal :: diffusion(*)
828) PetscReal :: fraction_upwind
829) PetscReal :: T_up(*), T_dn(*)
830)
831) PetscInt :: iphase
832) PetscReal :: coef_up, coef_dn
833) PetscReal :: q
834)
835) iphase = 1
836)
837) q = velocity(iphase)
838)
839) if (option%use_upwinding) then
840) ! upstream weighting
841) ! units = (m^3 water/m^2 bulk/sec)
842) if (q > 0.d0) then
843) coef_up = diffusion(iphase)+q
844) coef_dn = -diffusion(iphase)
845) else
846) coef_up = diffusion(iphase)
847) coef_dn = -diffusion(iphase)+q
848) endif
849) else
850) ! central difference, currently assuming uniform grid spacing
851) ! units = (m^3 water/m^2 bulk/sec)
852) !
853) coef_up = diffusion(iphase)+ (1.d0-fraction_upwind)*q
854) coef_dn = -diffusion(iphase)+ fraction_upwind*q
855) endif
856)
857) ! units = (m^3 water/m^2 bulk/sec)*(m^2 bulk)*(1000 L water/m^3 water)
858) ! = L water/sec
859) T_up(iphase) = coef_up*area*1000.d0 ! 1000 converts m^3 -> L
860) T_dn(iphase) = coef_dn*area*1000.d0
861)
862) ! CO2-specific
863) ! Add in multiphase, clu 12/29/08
864) if (option%iflowmode == MPH_MODE .or. option%iflowmode == IMS_MODE &
865) .or. option%iflowmode == FLASH2_MODE) then
866) do
867) iphase = iphase +1
868) if (iphase > option%nphase) exit
869) ! super critical CO2 phase have the index 2: need implementation
870) q = velocity(iphase)
871)
872) if (option%use_upwinding) then
873) !upstream weighting
874) ! units = (m^3 water/m^2 bulk/sec)
875) if (q > 0.d0) then
876) coef_up = diffusion(iphase)+q
877) coef_dn = -diffusion(iphase)
878) else
879) coef_up = diffusion(iphase)
880) coef_dn = -diffusion(iphase)+q
881) endif
882) else
883) coef_up = diffusion(iphase)+ (1.d0-fraction_upwind)*q
884) coef_dn = -diffusion(iphase)+ fraction_upwind*q
885) endif
886)
887) ! units = (m^3 water/m^2 bulk/sec)*(m^2 bulk)*(1000 L water/m^3 water)
888) ! = L water/sec
889) T_up(iphase) = coef_up*area*1000.d0 ! 1000 converts m^3 -> L
890) T_dn(iphase) = coef_dn*area*1000.d0
891)
892) enddo
893) endif
894)
895) end subroutine TFluxCoef
896)
897) ! ************************************************************************** !
898)
899) subroutine TFluxCoef_CD(option,area,velocity,diffusion,fraction_upwind, &
900) T_11,T_12,T_21,T_22)
901) !
902) ! Computes flux coefficients for transport matrix
903) !
904) ! Author: Glenn Hammond
905) ! Date: 02/22/10
906) !
907)
908) use Option_module
909)
910) implicit none
911)
912) type(option_type) :: option
913) PetscReal :: area
914) PetscReal :: velocity(*)
915) PetscReal :: diffusion(*)
916) PetscReal :: fraction_upwind
917) PetscReal :: T_11(*), T_12(*), T_21(*), T_22(*)
918)
919) PetscInt :: iphase
920) PetscReal :: coef_up, coef_dn
921) PetscReal :: tempreal
922) PetscReal :: advection_upwind(option%nphase)
923) PetscReal :: advection_downwind(option%nphase)
924) PetscReal :: q
925)
926) ! T_11 = diagonal term for upwind cell (row)
927) ! T_12 = off diagonal term for upwind cell (row)
928) ! T_21 = off diagonal term for downwind cell (row)
929) ! T_22 = diagonal term for downwind cell (row)
930)
931) ! Advection
932) if (option%use_upwinding) then
933) ! upstream weighting
934) ! units = (m^3 water/m^2 bulk/sec)
935) do iphase = 1, option%nphase
936) if (velocity(iphase) > 0.d0) then
937) advection_upwind(iphase) = velocity(iphase)
938) advection_downwind(iphase) = 0.d0
939) else
940) advection_upwind(iphase) = 0.d0
941) advection_downwind(iphase) = velocity(iphase)
942) endif
943) enddo
944) else
945) ! central difference
946) do iphase = 1, option%nphase
947) advection_upwind(iphase) = (1.d0-fraction_upwind)*velocity(iphase)
948) advection_downwind(iphase) = fraction_upwind*velocity(iphase)
949) enddo
950) endif
951)
952) tempreal = area*1000.d0
953) do iphase = 1, option%nphase
954) T_11(iphase) = (diffusion(iphase) + advection_upwind(iphase))*tempreal
955) T_12(iphase) = (-diffusion(iphase) + advection_downwind(iphase))*tempreal
956) ! T_21(iphase) = -(diffusion(iphase) + advection_upwind(iphase))*tempreal
957) ! T_22(iphase) = (diffusion(iphase) - advection_downwind(iphase))*tempreal
958) T_21(iphase) = -T_11(iphase)
959) T_22(iphase) = -T_12(iphase)
960) enddo
961)
962) end subroutine TFluxCoef_CD
963)
964) ! ************************************************************************** !
965)
966) subroutine TSrcSinkCoef(option,qsrc,tran_src_sink_type,T_in,T_out)
967) !
968) ! Computes src/sink coefficients for transport matrix
969) ! Here qsrc [m^3/sec] provided by flow.
970) !
971) ! Author: Glenn Hammond
972) ! Date: 01/12/11
973) !
974)
975) use Option_module
976)
977) implicit none
978)
979) type(option_type) :: option
980) PetscReal :: qsrc
981) PetscInt :: tran_src_sink_type
982) PetscReal :: T_in ! coefficient that scales concentration at cell
983) PetscReal :: T_out ! concentration that scales external concentration
984)
985) T_in = 0.d0
986) T_out = 0.d0
987)
988) select case(tran_src_sink_type)
989) case(EQUILIBRIUM_SS)
990) ! units should be mol/sec
991) ! 1.d-3 is a relatively large rate designed to equilibrate
992) ! the aqueous concentration with the concentrations specified at
993) ! the src/sink
994) T_in = 1.d-3 ! units L water/sec
995) T_out = -1.d0*T_in
996) case(MASS_RATE_SS)
997) ! in this case, rt_auxvar_bc%total actually holds the mass rate
998) T_in = 0.d0
999) T_out = -1.d0
1000) case default
1001) ! qsrc always in m^3/sec
1002) if (qsrc > 0.d0) then ! injection
1003) T_in = 0.d0
1004) T_out = -1.d0*qsrc*1000.d0 ! m^3/sec * 1000 L/m^3 -> L/s
1005) else
1006) T_out = 0.d0
1007) T_in = -1.d0*qsrc*1000.d0 ! m^3/sec * 1000 L/m^3 -> L/s
1008) endif
1009) end select
1010)
1011) ! Units of Tin & Tout should be L/s. When multiplied by Total (M) you get
1012) ! moles/sec, the units of the residual. To get the units of the Jacobian
1013) ! kg/sec, one must either scale by dtotal or den/1000. (kg/L).
1014)
1015) end subroutine TSrcSinkCoef
1016)
1017) ! ************************************************************************** !
1018)
1019) subroutine TFluxTVD(rt_parameter,velocity,area,dist, &
1020) total_up2,rt_auxvar_up, &
1021) rt_auxvar_dn,total_dn2, &
1022) TFluxLimitPtr, &
1023) option,flux)
1024) !
1025) ! Computes TVD flux term
1026) !
1027) ! Author: Glenn Hammond
1028) ! Date: 02/03/12
1029) !
1030)
1031) use Option_module
1032)
1033) implicit none
1034)
1035) type(reactive_transport_param_type) :: rt_parameter
1036) PetscReal :: velocity(:), area
1037) type(reactive_transport_auxvar_type) :: rt_auxvar_up, rt_auxvar_dn
1038) PetscReal, pointer :: total_up2(:,:), total_dn2(:,:)
1039) type(option_type) :: option
1040) PetscReal :: flux(rt_parameter%ncomp)
1041) procedure (TFluxLimiterDummy), pointer :: TFluxLimitPtr
1042) PetscReal :: dist(-1:3) ! list of distance vectors, size(-1:3,num_connections) where
1043) ! -1 = fraction upwind
1044) ! 0 = magnitude of distance
1045) ! 1-3 = components of unit vector
1046)
1047) PetscInt :: iphase
1048) PetscInt :: idof, ndof
1049) PetscReal :: dc, theta, correction, nu, velocity_area
1050)
1051) ndof = rt_parameter%naqcomp
1052)
1053) flux = 0.d0
1054)
1055) ! flux should be in mol/sec
1056)
1057) do iphase = 1, option%nphase
1058) nu = velocity(iphase)*option%tran_dt/dist(0)
1059) ! L/sec = m/sec * m^2 * 1000 [L/m^3]
1060) velocity_area = velocity(iphase)*area*1000.d0
1061) if (velocity_area >= 0.d0) then
1062) ! mol/sec = L/sec * mol/L
1063) flux = velocity_area*rt_auxvar_up%total(1:rt_parameter%naqcomp,iphase)
1064) if (associated(total_up2)) then
1065) do idof = 1, ndof
1066) dc = rt_auxvar_dn%total(idof,iphase) - &
1067) rt_auxvar_up%total(idof,iphase)
1068) if (dabs(dc) < 1.d-20) then
1069) theta = 1.d0
1070) else
1071) theta = (rt_auxvar_up%total(idof,iphase) - &
1072) total_up2(idof,iphase)) / &
1073) dc
1074) endif
1075) ! mol/sec = L/sec * mol/L
1076) correction = 0.5d0*velocity_area*(1.d0-nu)* &
1077) TFluxLimitPtr(theta)* &
1078) dc
1079) flux(idof) = flux(idof) + correction
1080) enddo
1081) endif
1082) else
1083) flux = velocity_area*rt_auxvar_dn%total(1:rt_parameter%naqcomp,iphase)
1084) if (associated(total_dn2)) then
1085) do idof = 1, ndof
1086) dc = rt_auxvar_dn%total(idof,iphase) - &
1087) rt_auxvar_up%total(idof,iphase)
1088) if (dabs(dc) < 1.d-20) then
1089) theta = 1.d0
1090) else
1091) theta = (total_dn2(idof,iphase) - &
1092) rt_auxvar_dn%total(idof,iphase)) / &
1093) dc
1094) endif
1095) correction = 0.5d0*velocity_area*(1.d0+nu)* &
1096) TFluxLimitPtr(theta)* &
1097) dc
1098) flux(idof) = flux(idof) + correction
1099) enddo
1100) endif
1101) endif
1102) enddo
1103)
1104) end subroutine TFluxTVD
1105)
1106) ! ************************************************************************** !
1107)
1108) function TFluxLimiter(theta)
1109) !
1110) ! Applies flux limiter
1111) !
1112) ! Author: Glenn Hammond
1113) ! Date: 02/03/12
1114) !
1115)
1116) implicit none
1117)
1118) PetscReal :: theta
1119)
1120) PetscReal :: TFluxLimiter
1121)
1122) ! Linear
1123) !---------
1124) ! upwind
1125) TFluxLimiter = 0.d0
1126) ! Lax-Wendroff
1127) !TFluxLimiter = 1.d0
1128) ! Beam-Warming
1129) !TFluxLimiter = theta
1130) ! Fromm
1131) !TFluxLimiter = 0.5d0*(1.d0+theta)
1132)
1133) ! Higher-order
1134) !---------
1135) ! minmod
1136) !TFluxLimiter = max(0.d0,min(1.d0,theta))
1137) ! superbee
1138) !TFluxLimiter = max(0.d0,min(1.d0,2.d0*theta),min(2.d0,theta))
1139) ! MC
1140) !TFluxLimiter = max(0.d0,min((1.d0+theta)/2.d0,2.d0,2.d0*theta))
1141) ! van Leer
1142) !TFluxLimiter = (theta+dabs(theta))/(1.d0+dabs(theta)
1143)
1144) end function TFluxLimiter
1145)
1146) ! ************************************************************************** !
1147)
1148) function TFluxLimitUpwind(theta)
1149) !
1150) ! Applies an upwind flux limiter
1151) !
1152) ! Author: Glenn Hammond
1153) ! Date: 02/03/12
1154) !
1155)
1156) implicit none
1157)
1158) PetscReal :: theta
1159)
1160) PetscReal :: TFluxLimitUpwind
1161)
1162) ! upwind
1163) TFluxLimitUpwind = 0.d0
1164)
1165) end function TFluxLimitUpwind
1166)
1167) ! ************************************************************************** !
1168)
1169) function TFluxLimitMinmod(theta)
1170) !
1171) ! Applies a minmod flux limiter
1172) !
1173) ! Author: Glenn Hammond
1174) ! Date: 02/03/12
1175) !
1176)
1177) implicit none
1178)
1179) PetscReal :: theta
1180)
1181) PetscReal :: TFluxLimitMinmod
1182)
1183) ! minmod
1184) TFluxLimitMinmod = max(0.d0,min(1.d0,theta))
1185)
1186) end function TFluxLimitMinmod
1187)
1188) ! ************************************************************************** !
1189)
1190) function TFluxLimitMC(theta)
1191) !
1192) ! Applies an MC flux limiter
1193) !
1194) ! Author: Glenn Hammond
1195) ! Date: 02/03/12
1196) !
1197)
1198) implicit none
1199)
1200) PetscReal :: theta
1201)
1202) PetscReal :: TFluxLimitMC
1203)
1204) ! MC
1205) TFluxLimitMC = max(0.d0,min((1.d0+theta)/2.d0,2.d0,2.d0*theta))
1206)
1207) end function TFluxLimitMC
1208)
1209) ! ************************************************************************** !
1210)
1211) function TFluxLimitSuperBee(theta)
1212) !
1213) ! Applies an superbee flux limiter
1214) !
1215) ! Author: Glenn Hammond
1216) ! Date: 02/03/12
1217) !
1218)
1219) implicit none
1220)
1221) PetscReal :: theta
1222)
1223) PetscReal :: TFluxLimitSuperBee
1224)
1225) ! superbee
1226) TFluxLimitSuperBee = max(0.d0,min(1.d0,2.d0*theta),min(2.d0,theta))
1227)
1228) end function TFluxLimitSuperBee
1229)
1230) ! ************************************************************************** !
1231)
1232) function TFluxLimitVanLeer(theta)
1233) !
1234) ! Applies an van Leer flux limiter
1235) !
1236) ! Author: Glenn Hammond
1237) ! Date: 02/03/12
1238) !
1239)
1240) implicit none
1241)
1242) PetscReal :: theta
1243)
1244) PetscReal :: TFluxLimitVanLeer
1245)
1246) ! superbee
1247) TFluxLimitVanLeer = (theta+dabs(theta))/(1.d0+dabs(theta))
1248)
1249) end function TFluxLimitVanLeer
1250)
1251) ! ************************************************************************** !
1252)
1253) end module Transport_module