reaction_sandbox_example.F90 coverage: 0.00 %func 0.00 %block
1) module Reaction_Sandbox_Example_class
2)
3) ! 1. Change all references to "Example" as desired to rename the module and
4) ! and subroutines within the module.
5)
6) use Reaction_Sandbox_Base_class
7)
8) use Global_Aux_module
9) use Reactive_Transport_Aux_module
10)
11) use PFLOTRAN_Constants_module
12)
13) implicit none
14)
15) private
16)
17) #include "petsc/finclude/petscsys.h"
18)
19) ! 2. Add module variables here. Note that one must use the PETSc data types
20) ! PetscInt, PetscReal, PetscBool to declare variables of type integer
21) ! float/real*8, and logical respectively. E.g.
22) !
23) ! PetscReal, parameter :: formula_weight_of_water = 18.01534d0
24)
25) type, public, &
26) extends(reaction_sandbox_base_type) :: reaction_sandbox_example_type
27) ! 3. Add variables/arrays associated with new reaction.
28) character(len=MAXWORDLENGTH) :: species_name
29) PetscInt :: species_id
30) PetscReal :: rate_constant
31) contains
32) procedure, public :: ReadInput => ExampleRead
33) procedure, public :: Setup => ExampleSetup
34) procedure, public :: Evaluate => ExampleReact
35) procedure, public :: Destroy => ExampleDestroy
36) end type reaction_sandbox_example_type
37)
38) public :: ExampleCreate
39)
40) contains
41)
42) ! ************************************************************************** !
43)
44) function ExampleCreate()
45) !
46) ! Allocates example reaction object.
47) !
48) ! Author: John Doe (replace in all subroutine headers with name of developer)
49) ! Date: 00/00/00 (replace in all subroutine headers with current date)
50) !
51)
52) implicit none
53)
54) class(reaction_sandbox_example_type), pointer :: ExampleCreate
55)
56) ! 4. Add code to allocate object and initialized all variables to zero and
57) ! nullify all pointers. E.g.
58) allocate(ExampleCreate)
59) ExampleCreate%species_name = ''
60) ExampleCreate%species_id = 0
61) ExampleCreate%rate_constant = 0.d0
62) nullify(ExampleCreate%next)
63)
64) end function ExampleCreate
65)
66) ! ************************************************************************** !
67)
68) subroutine ExampleRead(this,input,option)
69) !
70) ! Reads input deck for example reaction parameters (if any)
71) !
72) ! Author: John Doe
73) ! Date: 00/00/00
74) !
75)
76) use Option_module
77) use String_module
78) use Input_Aux_module
79) use Units_module, only : UnitsConvertToInternal
80)
81) implicit none
82)
83) class(reaction_sandbox_example_type) :: this
84) type(input_type), pointer :: input
85) type(option_type) :: option
86)
87) PetscInt :: i
88) character(len=MAXWORDLENGTH) :: word, internal_units
89)
90) do
91) call InputReadPflotranString(input,option)
92) if (InputError(input)) exit
93) if (InputCheckExit(input,option)) exit
94)
95) call InputReadWord(input,option,word,PETSC_TRUE)
96) call InputErrorMsg(input,option,'keyword', &
97) 'CHEMISTRY,REACTION_SANDBOX,TEMPLATE')
98) call StringToUpper(word)
99)
100) select case(trim(word))
101)
102) ! Example Input:
103)
104) ! CHEMISTRY
105) ! ...
106) ! REACTION_SANDBOX
107) ! : begin user-defined input
108) ! TEMPLATE
109) ! EXAMPLE_INTEGER 1
110) ! EXAMPLE_INTEGER_ARRAY 2 3 4
111) ! END
112) ! : end user defined input
113) ! END
114) ! ...
115) ! END
116)
117) ! 5. Add case statement for reading variables.
118) case('SPECIES_NAME')
119) ! 6. Read the variable
120) ! Read the character string indicating which of the primary species
121) ! is being decayed.
122) call InputReadWord(input,option,this%species_name,PETSC_TRUE)
123) ! 7. Inform the user of any errors if not read correctly.
124) call InputErrorMsg(input,option,'species_name', &
125) 'CHEMISTRY,REACTION_SANDBOX,EXAMPLE')
126) ! 8. Repeat for other variables
127) case('RATE_CONSTANT')
128) ! Read the double precision rate constant
129) call InputReadDouble(input,option,this%rate_constant)
130) call InputErrorMsg(input,option,'rate_constant', &
131) 'CHEMISTRY,REACTION_SANDBOX,EXAMPLE')
132) ! Read the units
133) call InputReadWord(input,option,word,PETSC_TRUE)
134) if (InputError(input)) then
135) ! If units do not exist, assume default units of 1/s which are the
136) ! standard internal PFLOTRAN units for this rate constant.
137) input%err_buf = 'REACTION_SANDBOX,EXAMPLE,RATE CONSTANT UNITS'
138) call InputDefaultMsg(input,option)
139) else
140) ! If units exist, convert to internal units of 1/s
141) internal_units = 'unitless/sec'
142) this%rate_constant = this%rate_constant * &
143) UnitsConvertToInternal(word,internal_units,option)
144) endif
145) case default
146) call InputKeywordUnrecognized(word, &
147) 'CHEMISTRY,REACTION_SANDBOX,TEMPLATE',option)
148) end select
149) enddo
150)
151) end subroutine ExampleRead
152)
153) ! ************************************************************************** !
154)
155) subroutine ExampleSetup(this,reaction,option)
156) !
157) ! Sets up the example reaction either with parameters either
158) ! read from the input deck or hardwired.
159) !
160) ! Author: John Doe
161) ! Date: 00/00/00
162) !
163)
164) use Reaction_Aux_module, only : reaction_type, GetPrimarySpeciesIDFromName
165) use Option_module
166)
167) implicit none
168)
169) class(reaction_sandbox_example_type) :: this
170) type(reaction_type) :: reaction
171) type(option_type) :: option
172)
173) ! 9. Add code to initialize
174) this%species_id = &
175) GetPrimarySpeciesIDFromName(this%species_name,reaction,option)
176)
177) end subroutine ExampleSetup
178)
179) ! ************************************************************************** !
180)
181) subroutine ExampleReact(this,Residual,Jacobian,compute_derivative, &
182) rt_auxvar,global_auxvar,material_auxvar,reaction, &
183) option)
184) !
185) ! Evaluates reaction storing residual and/or Jacobian
186) !
187) ! Author: John Doe
188) ! Date: 00/00/00
189) !
190)
191) use Option_module
192) use Reaction_Aux_module
193) use Material_Aux_class
194)
195) implicit none
196)
197) class(reaction_sandbox_example_type) :: this
198) type(option_type) :: option
199) type(reaction_type) :: reaction
200) PetscBool :: compute_derivative
201) ! the following arrays must be declared after reaction
202) PetscReal :: Residual(reaction%ncomp)
203) PetscReal :: Jacobian(reaction%ncomp,reaction%ncomp)
204) type(reactive_transport_auxvar_type) :: rt_auxvar
205) type(global_auxvar_type) :: global_auxvar
206) class(material_auxvar_type) :: material_auxvar
207)
208) PetscInt, parameter :: iphase = 1
209) PetscReal :: L_water
210)
211) ! Description of subroutine arguments:
212)
213) ! Residual - 1D array storing residual entries in units mol/sec
214) ! Jacobian - 2D array storing Jacobian entires in units kg water/sec
215) !
216) ! Jacobian [kg water/sec] * dc [mol/kg water] = -Res [mol/sec]
217) !
218) ! compute_derivative - Flag indicating whether analtical derivative should
219) ! be calculated. The user must provide either the analytical derivatives
220) ! or a numerical approximation unless always running with
221) ! NUMERICAL_JACOBIAN_RXN defined in input deck. If the use of
222) ! NUMERICAL_JACOBIAN_RXN is assumed, the user should provide an error
223) ! message when compute_derivative is true. E.g.
224) !
225) ! option%io_buffer = 'NUMERICAL_JACOBIAN_RXN must always be used ' // &
226) ! 'due to assumptions in Example'
227) ! call printErrMsg(option)
228) !
229) ! rt_auxvar - Object holding chemistry information (e.g. concentrations,
230) ! activity coefficients, mineral volume fractions, etc.). See
231) ! reactive_transport_aux.F90.
232) !
233) ! Useful variables:
234) ! rt_auxvar%total(:,iphase) - total component concentrations
235) ! [mol/L water] for phase
236) ! rt_auxvar%pri_molal(:) - free ion concentrations [mol/kg water]
237) ! rt_auxvar%pri_act_coef(:) - activity coefficients for primary species
238) ! rt_auxvar%aqueous%dtotal(:,iphase) - derivative of total component
239) ! concentration with respect to free ion [kg water/L water]
240) !
241) ! global_auxvar - Object holding information on flow (e.g. saturation,
242) ! density, viscosity, temperature, etc)
243) !
244) ! Useful variables:
245) ! global_auxvar%den(iphase) - fluid density [mol/m^3]
246) ! global_auxvar%den_kg(iphase) - fluid density [kg/m^3]
247) ! global_auxvar%sat(iphase) - saturation
248) ! global_auxvar%temp - temperature [C]
249) !
250) ! porosity - effective porosity of grid cell [m^3 pore/m^3 bulk]
251) ! volume - volume of grid cell [m^3]
252) ! reaction - Provides access to variable describing chemistry. E.g.
253) ! reaction%ncomp - # chemical degrees of freedom (mobile and immobile)
254) ! reaction%naqcomp - # chemical degrees of freedom on water
255) ! reaction%primary_species_names(:) - names of primary species
256) !
257) ! option - Provides handle for controlling simulation, catching and
258) ! reporting errors.
259)
260) ! Unit of the residual must be in moles/second
261) ! global_auxvar%sat(iphase) = saturation of cell
262) ! 1.d3 converts m^3 water -> L water
263) L_water = material_auxvar%porosity*global_auxvar%sat(iphase)* &
264) material_auxvar%volume*1.d3
265) ! always subtract contribution from residual
266) Residual(this%species_id) = Residual(this%species_id) - &
267) this%rate_constant * & ! 1/sec
268) L_water * & ! L water
269) ! rt_auxvar%total(this%species_id,iphase) = species total component
270) ! concentration
271) rt_auxvar%total(this%species_id,iphase) ! mol/L water
272)
273)
274)
275) if (compute_derivative) then
276)
277) ! 11. If using an analytical Jacobian, add code for Jacobian evaluation
278)
279) ! always add contribution to Jacobian
280) ! units = (mol/sec)*(kg water/mol) = kg water/sec
281) Jacobian(this%species_id,this%species_id) = &
282) Jacobian(this%species_id,this%species_id) + &
283) this%rate_constant * & ! 1/sec
284) L_water * & ! L water
285) ! kg water/L water
286) ! rt_auxvar%aqueous%dtotal(this%species_id,this%species_id,iphase) =
287) ! derivative of total component concentration with respect to the
288) ! free ion concentration of the same species.
289) rt_auxvar%aqueous%dtotal(this%species_id,this%species_id,iphase)
290)
291) endif
292)
293) end subroutine ExampleReact
294)
295) ! ************************************************************************** !
296)
297) subroutine ExampleDestroy(this)
298) !
299) ! Destroys allocatable or pointer objects created in this
300) ! module
301) !
302) ! Author: John Doe
303) ! Date: 00/00/00
304) !
305)
306) implicit none
307)
308) class(reaction_sandbox_example_type) :: this
309)
310) ! 12. Add code to deallocate contents of the example object
311)
312) end subroutine ExampleDestroy
313)
314) end module Reaction_Sandbox_Example_class