reaction_sandbox_simple.F90 coverage: 0.00 %func 0.00 %block
1) module Reaction_Sandbox_Simple_class
2)
3) use Reaction_Sandbox_Base_class
4) use PFLOTRAN_Constants_module
5)
6) implicit none
7)
8) private
9)
10) #include "petsc/finclude/petscsys.h"
11)
12) type, public, &
13) extends(reaction_sandbox_base_type) :: reaction_sandbox_simple_type
14) ! Aqueous species
15) PetscInt :: species_Aaq_id
16) PetscInt :: species_Baq_id
17) PetscInt :: species_Caq_id
18) PetscInt :: species_Daq_id
19) PetscInt :: species_Eaq_id
20) PetscInt :: species_Faq_id
21) ! Immobile species (e.g. biomass)
22) PetscInt :: species_Xim_id
23) PetscInt :: species_Yim_id
24) contains
25) procedure, public :: Setup => SimpleSetup
26) procedure, public :: Evaluate => SimpleReact
27) end type reaction_sandbox_simple_type
28)
29) public :: SimpleCreate
30)
31) contains
32)
33) ! ************************************************************************** !
34)
35) function SimpleCreate()
36) !
37) ! Allocates simple reaction object.
38) !
39) ! Author: Glenn Hammond
40) ! Date: 12/03/15
41)
42) implicit none
43)
44) class(reaction_sandbox_simple_type), pointer :: SimpleCreate
45)
46) allocate(SimpleCreate)
47) nullify(SimpleCreate%next)
48)
49) end function SimpleCreate
50)
51) ! ************************************************************************** !
52)
53) subroutine SimpleSetup(this,reaction,option)
54) !
55) ! Sets up the simple reaction with hardwired parameters
56) !
57) ! Author: Glenn Hammond
58) ! Date: 12/03/15
59)
60) use Reaction_Aux_module, only : reaction_type, GetPrimarySpeciesIDFromName
61) use Reaction_Immobile_Aux_module, only : GetImmobileSpeciesIDFromName
62) use Option_module
63)
64) implicit none
65)
66) class(reaction_sandbox_simple_type) :: this
67) type(reaction_type) :: reaction
68) type(option_type) :: option
69)
70) character(len=MAXWORDLENGTH) :: word
71)
72) ! Aqueous species
73) word = 'Aaq'
74) this%species_Aaq_id = &
75) GetPrimarySpeciesIDFromName(word,reaction,option)
76) word = 'Baq'
77) this%species_Baq_id = &
78) GetPrimarySpeciesIDFromName(word,reaction,option)
79) word = 'Caq'
80) this%species_Caq_id = &
81) GetPrimarySpeciesIDFromName(word,reaction,option)
82) word = 'Daq'
83) this%species_Daq_id = &
84) GetPrimarySpeciesIDFromName(word,reaction,option)
85) word = 'Eaq'
86) this%species_Eaq_id = &
87) GetPrimarySpeciesIDFromName(word,reaction,option)
88) word = 'Faq'
89) this%species_Faq_id = &
90) GetPrimarySpeciesIDFromName(word,reaction,option)
91)
92) ! Immobile species
93) word = 'Xim'
94) this%species_Xim_id = &
95) GetImmobileSpeciesIDFromName(word,reaction%immobile,option)
96) word = 'Yim'
97) this%species_Yim_id = &
98) GetImmobileSpeciesIDFromName(word,reaction%immobile,option)
99)
100)
101) end subroutine SimpleSetup
102)
103) ! ************************************************************************** !
104)
105) subroutine SimpleReact(this,Residual,Jacobian,compute_derivative, &
106) rt_auxvar,global_auxvar,material_auxvar,reaction, &
107) option)
108) !
109) ! Evaluates reaction storing residual and/or Jacobian
110) !
111) ! Author: Glenn Hammond
112) ! Date: 12/03/15
113) !
114)
115) use Option_module
116) use Reaction_Aux_module
117) use Reactive_Transport_Aux_module
118) use Global_Aux_module
119) use Material_Aux_class
120)
121) implicit none
122)
123) class(reaction_sandbox_simple_type) :: this
124) type(option_type) :: option
125) type(reaction_type) :: reaction
126) PetscBool :: compute_derivative
127) ! the following arrays must be declared after reaction
128) PetscReal :: Residual(reaction%ncomp)
129) PetscReal :: Jacobian(reaction%ncomp,reaction%ncomp)
130) type(reactive_transport_auxvar_type) :: rt_auxvar
131) type(global_auxvar_type) :: global_auxvar
132) class(material_auxvar_type) :: material_auxvar
133)
134) PetscInt, parameter :: iphase = 1
135) PetscReal :: volume ! m^3 bulk
136) PetscReal :: porosity ! m^3 pore space / m^3 bulk
137) PetscReal :: liquid_saturation ! m^3 water / m^3 pore space
138) PetscReal :: L_water ! L water
139)
140) PetscReal :: Aaq, Baq, Caq, Daq, Eaq, Faq ! mol/L
141) PetscReal :: Xim, Yim ! mol/m^3
142) PetscReal :: Rate
143) PetscReal :: RateA, RateB, RateC, RateD, RateE, RateF, RateX, RateY ! mol/sec
144) PetscReal :: stoichA, stoichB, stoichC, stoichD, stoichE, stoichF
145) PetscReal :: stoichX, stoichY
146) PetscReal :: k, kr ! units are problem specific
147) PetscReal :: K_Aaq, K_Baq ! [mol/L]
148)
149) porosity = material_auxvar%porosity
150) liquid_saturation = global_auxvar%sat(iphase)
151) volume = material_auxvar%volume
152) L_water = porosity*liquid_saturation*volume*1.d3 ! 1.d3 converts m^3 water -> L water
153)
154) Aaq = rt_auxvar%total(this%species_Aaq_id,iphase)
155) Baq = rt_auxvar%total(this%species_Baq_id,iphase)
156) Caq = rt_auxvar%total(this%species_Caq_id,iphase)
157) Daq = rt_auxvar%total(this%species_Daq_id,iphase)
158) Eaq = rt_auxvar%total(this%species_Eaq_id,iphase)
159) Faq = rt_auxvar%total(this%species_Faq_id,iphase)
160)
161) Xim = rt_auxvar%immobile(this%species_Xim_id)
162) Yim = rt_auxvar%immobile(this%species_Yim_id)
163)
164) ! initialize all rates to zero
165) Rate = 0.d0
166)
167) RateA = 0.d0
168) RateB = 0.d0
169) RateC = 0.d0
170) RateD = 0.d0
171) RateE = 0.d0
172) RateF = 0.d0
173) RateX = 0.d0
174) RateY = 0.d0
175)
176) ! stoichiometries
177) ! reactants have negative stoichiometry
178) ! products have positive stoichiometry
179) stoichA = 0.d0
180) stoichB = 0.d0
181) stoichC = 0.d0
182) stoichD = 0.d0
183) stoichE = 0.d0
184) stoichF = 0.d0
185) stoichX = 0.d0
186) stoichY = 0.d0
187)
188) ! kinetic rate constants
189) k = 0.d0
190) kr = 0.d0
191)
192) ! Monod half-saturation constants
193) K_Aaq = 0.d0
194) K_Baq = 0.d0
195)
196) ! zero-order (A -> C)
197) !k = 0.d0 ! WARNING: Too high a rate can result in negative concentrations
198) ! which are non-physical and the code will not converge.
199) !stoichA = -1.d0
200) !stoichC = 1.d0
201) !Rate = k * L_water
202) !RateA = stoichA * Rate
203) !RateC = stoichC * Rate
204)
205) ! first-order (A -> C)
206) !k = 1.d-10
207) !stoichA = -1.d0
208) !stoichC = 1.d0
209) !Rate = k * Aaq * L_water
210) !RateA = stoichA * Rate
211) !RateC = stoichC * Rate
212)
213) ! second-order (A + B -> C)
214) !k = 1.d-10
215) !stoichA = -1.d0
216) !stoichB = -1.d0
217) !stoichC = 1.d0
218) !Rate = k * Aaq * Baq * L_water
219) !RateA = stoichA * Rate
220) !RateB = stoichB * Rate
221) !RateC = stoichC * Rate
222)
223) ! Monod (A -> C)
224) !k = 1.d-12
225) !K_Aaq = 5.d-4
226) !stoichA = -1.d0
227) !stoichC = 1.d0
228) !Rate = k * Aaq / (K_Aaq + Aaq) * L_water
229) !RateA = stoichA * Rate
230) !RateC = stoichC * Rate
231)
232) ! multiplicative Monod w/biomass
233) ! A + 2B -> C
234) !k = 1.d-7
235) !K_Aaq = 5.d-4
236) !K_Baq = 5.d-4
237) !stoichA = -1.d0
238) !stoichB = -2.d0
239) !stoichC = 1.d0
240) !Rate = k * Xim * Aaq / (K_Aaq + Aaq) * Baq / (K_Baq + Baq) * L_water
241) !RateA = stoichA * Rate
242) !RateB = stoichB * Rate
243) !RateC = stoichC * Rate
244)
245) ! first-order forward - reverse (A <-> C)
246) !k = 1.d-6
247) !kr = 1.d-7
248) !stoichA = -1.d0
249) !stoichC = 1.d0
250) !Rate = (k * Aaq - kr * Caq) * L_water
251) !RateA = stoichA * Rate
252) !RateC = stoichC * Rate
253)
254) ! NOTES
255) ! 1. Always subtract contribution from residual
256) ! 2. Units of residual are moles/second
257) Residual(this%species_Aaq_id) = Residual(this%species_Aaq_id) - RateA
258) Residual(this%species_Baq_id) = Residual(this%species_Baq_id) - RateB
259) Residual(this%species_Caq_id) = Residual(this%species_Caq_id) - RateC
260) Residual(this%species_Daq_id) = Residual(this%species_Daq_id) - RateD
261) Residual(this%species_Eaq_id) = Residual(this%species_Eaq_id) - RateE
262) Residual(this%species_Faq_id) = Residual(this%species_Faq_id) - RateF
263) Residual(this%species_Xim_id) = Residual(this%species_Xim_id) - RateX
264) Residual(this%species_Yim_id) = Residual(this%species_Yim_id) - RateY
265)
266) end subroutine SimpleReact
267)
268) end module Reaction_Sandbox_Simple_class