reaction_mineral.F90 coverage: 87.50 %func 61.11 %block
1) module Reaction_Mineral_module
2)
3) use Reaction_Mineral_Aux_module
4) use Reaction_Aux_module
5) use Reactive_Transport_Aux_module
6) use Global_Aux_module
7)
8) use PFLOTRAN_Constants_module
9)
10) implicit none
11)
12) private
13)
14) #include "petsc/finclude/petscsys.h"
15)
16) PetscReal, parameter :: perturbation_tolerance = 1.d-5
17)
18) public :: MineralRead, &
19) MineralReadKinetics, &
20) MineralReadFromDatabase, &
21) MineralProcessConstraint, &
22) RKineticMineral, &
23) RMineralSaturationIndex, &
24) MineralUpdateTempDepCoefs
25)
26) contains
27)
28) ! ************************************************************************** !
29)
30) subroutine MineralRead(mineral,input,option)
31) !
32) ! Reads chemical species
33) !
34) ! Author: Glenn Hammond
35) ! Date: 08/16/12
36) !
37)
38) use Option_module
39) use String_module
40) use Input_Aux_module
41) use Utility_module
42)
43) implicit none
44)
45) type(mineral_type) :: mineral
46) type(input_type), pointer :: input
47) type(option_type) :: option
48)
49) type(mineral_rxn_type), pointer :: cur_mineral, prev_mineral
50)
51) nullify(prev_mineral)
52) do
53) call InputReadPflotranString(input,option)
54) if (InputError(input)) exit
55) if (InputCheckExit(input,option)) exit
56)
57) mineral%nmnrl = mineral%nmnrl + 1
58)
59) cur_mineral => MineralRxnCreate()
60) call InputReadWord(input,option,cur_mineral%name,PETSC_TRUE)
61) call InputErrorMsg(input,option,'keyword','CHEMISTRY,MINERALS')
62) if (.not.associated(mineral%mineral_list)) then
63) mineral%mineral_list => cur_mineral
64) cur_mineral%id = 1
65) endif
66) if (associated(prev_mineral)) then
67) prev_mineral%next => cur_mineral
68) cur_mineral%id = prev_mineral%id + 1
69) endif
70) prev_mineral => cur_mineral
71) nullify(cur_mineral)
72) enddo
73)
74) end subroutine MineralRead
75)
76) ! ************************************************************************** !
77)
78) subroutine MineralReadKinetics(mineral,input,option)
79) !
80) ! Reads mineral kinetics
81) !
82) ! Author: Glenn Hammond
83) ! Date: 10/16/08
84) !
85)
86) use Input_Aux_module
87) use String_module
88) use Option_module
89) use Units_module
90)
91) implicit none
92)
93) type(mineral_type) :: mineral
94) type(input_type), pointer :: input
95) type(option_type) :: option
96)
97) character(len=MAXSTRINGLENGTH) :: string
98) character(len=MAXSTRINGLENGTH) :: error_string
99) character(len=MAXWORDLENGTH) :: word
100) character(len=MAXWORDLENGTH) :: name
101) character(len=MAXWORDLENGTH) :: card
102) character(len=MAXWORDLENGTH) :: internal_units
103)
104) type(mineral_rxn_type), pointer :: cur_mineral
105) type(transition_state_rxn_type), pointer :: tstrxn, cur_tstrxn
106) type(transition_state_prefactor_type), pointer :: prefactor, &
107) cur_prefactor
108) type(ts_prefactor_species_type), pointer :: prefactor_species, &
109) cur_prefactor_species
110) PetscBool :: found
111) PetscInt :: imnrl,icount
112) PetscReal :: temp_real
113)
114) cur_mineral => mineral%mineral_list
115) do
116) if (.not.associated(cur_mineral)) exit
117) cur_mineral%id = -1*abs(cur_mineral%id)
118) cur_mineral => cur_mineral%next
119) enddo
120)
121) input%ierr = 0
122) icount = 0
123) do
124)
125) call InputReadPflotranString(input,option)
126) if (InputError(input)) exit
127) if (InputCheckExit(input,option)) exit
128)
129) call InputReadWord(input,option,name,PETSC_TRUE)
130) call InputErrorMsg(input,option,'keyword','CHEMISTRY,MINERAL_KINETICS')
131)
132) cur_mineral => mineral%mineral_list
133) found = PETSC_FALSE
134) do
135) if (.not.associated(cur_mineral)) exit
136) if (StringCompare(cur_mineral%name,name,MAXWORDLENGTH)) then
137) found = PETSC_TRUE
138) cur_mineral%itype = MINERAL_KINETIC
139) tstrxn => TransitionStateTheoryRxnCreate()
140) ! initialize to UNINITIALIZED_INTEGER to ensure that it is set
141) tstrxn%rate = UNINITIALIZED_DOUBLE
142) do
143) call InputReadPflotranString(input,option)
144) call InputReadStringErrorMsg(input,option,card)
145) if (InputCheckExit(input,option)) exit
146) call InputReadWord(input,option,word,PETSC_TRUE)
147) error_string = 'CHEMISTRY,MINERAL_KINETICS'
148) call InputErrorMsg(input,option,'word',error_string)
149)
150) select case(trim(word))
151) case('RATE_CONSTANT')
152) ! read rate constant
153) call InputReadDouble(input,option,tstrxn%rate)
154) if (tstrxn%rate < 0.d0) then
155) tstrxn%rate = 10.d0**tstrxn%rate
156) endif
157) call InputErrorMsg(input,option,'rate',error_string)
158) ! read units if they exist
159) internal_units = 'mol/m^2-sec'
160) call InputReadWord(input,option,word,PETSC_TRUE)
161) if (InputError(input)) then
162) input%err_buf = trim(cur_mineral%name) // ' RATE UNITS'
163) call InputDefaultMsg(input,option)
164) else
165) tstrxn%rate = tstrxn%rate * &
166) UnitsConvertToInternal(word,internal_units,option)
167) endif
168) case('ACTIVATION_ENERGY')
169) ! read activation energy for Arrhenius law
170) call InputReadDouble(input,option,tstrxn%activation_energy)
171) call InputErrorMsg(input,option,'activation',error_string)
172) case('AFFINITY_THRESHOLD')
173) ! read affinity threshold for precipitation
174) call InputReadDouble(input,option,tstrxn%affinity_threshold)
175) call InputErrorMsg(input,option,'affinity threshold', &
176) error_string)
177) case('AFFINITY_POWER')
178) ! reads exponent on affinity term
179) call InputReadDouble(input,option,tstrxn%affinity_factor_beta)
180) call InputErrorMsg(input,option,'affinity power',error_string)
181) case('MINERAL_SCALE_FACTOR')
182) ! read mineral scale factor term
183) call InputReadDouble(input,option,tstrxn%min_scale_factor)
184) call InputErrorMsg(input,option,"Mineral scale fac", &
185) error_string)
186) case('TEMKIN_CONSTANT')
187) ! reads exponent on affinity term
188) call InputReadDouble(input,option,tstrxn%affinity_factor_sigma)
189) call InputErrorMsg(input,option,"Temkin's constant", &
190) error_string)
191) case('SURFACE_AREA_POROSITY_POWER')
192) call InputReadDouble(input,option,tstrxn%surf_area_porosity_pwr)
193) call InputErrorMsg(input,option,'surface area porosity power', &
194) error_string)
195) case('SURFACE_AREA_VOL_FRAC_POWER')
196) call InputReadDouble(input,option,tstrxn%surf_area_vol_frac_pwr)
197) call InputErrorMsg(input,option, &
198) 'surface area volume fraction power', &
199) error_string)
200) case('RATE_LIMITER')
201) ! read rate limiter for precipitation
202) call InputReadDouble(input,option,tstrxn%rate_limiter)
203) call InputErrorMsg(input,option,'rate_limiter',error_string)
204) case('IRREVERSIBLE')
205) ! read flag for irreversible reaction
206) option%io_buffer = 'IRREVERSIBLE mineral precipitation/' // &
207) 'dissolution no longer supported. The code is commented out.'
208) call printErrMsg(option)
209) tstrxn%irreversible = 1
210) call InputErrorMsg(input,option,'irreversible',error_string)
211) case('ARMOR_MINERAL')
212) ! read mineral name
213) call InputReadWord(input,option,tstrxn%armor_min_name,PETSC_TRUE)
214) call InputErrorMsg(input,option,'name',error_string)
215) case('ARMOR_PWR')
216) ! read power law exponent
217) call InputReadDouble(input,option,tstrxn%armor_pwr)
218) call InputErrorMsg(input,option,'armor_pwr',error_string)
219) case('ARMOR_CRIT_VOL_FRAC')
220) ! read critical volume fraction
221) call InputReadDouble(input,option,tstrxn%armor_crit_vol_frac)
222) call InputErrorMsg(input,option,'armor_crit_vol_frac',error_string)
223)
224) case('PREFACTOR')
225) error_string = 'CHEMISTRY,MINERAL_KINETICS,PREFACTOR'
226) prefactor => TransitionStatePrefactorCreate()
227) ! Initialize to UNINITIALIZED_DOUBLE to check later whether they were set
228) prefactor%rate = UNINITIALIZED_DOUBLE
229) prefactor%activation_energy = UNINITIALIZED_DOUBLE
230) do
231) call InputReadPflotranString(input,option)
232) call InputReadStringErrorMsg(input,option,card)
233) if (InputCheckExit(input,option)) exit
234) call InputReadWord(input,option,word,PETSC_TRUE)
235) call InputErrorMsg(input,option,'word',error_string)
236) select case(trim(word))
237) case('RATE_CONSTANT')
238) ! read rate constant
239) call InputReadDouble(input,option,prefactor%rate)
240) call InputErrorMsg(input,option,'rate',error_string)
241) if (prefactor%rate < 0.d0) then
242) prefactor%rate = 10.d0**prefactor%rate
243) endif
244) ! read units if they exist
245) internal_units = 'mol/m^2-sec'
246) call InputReadWord(input,option,word,PETSC_TRUE)
247) if (InputError(input)) then
248) input%err_buf = trim(cur_mineral%name) // &
249) 'PREFACTOR RATE UNITS'
250) call InputDefaultMsg(input,option)
251) else
252) prefactor%rate = prefactor%rate * &
253) UnitsConvertToInternal(word,internal_units,option)
254) endif
255) case('ACTIVATION_ENERGY')
256) ! read activation energy for Arrhenius law
257) call InputReadDouble(input,option, &
258) prefactor%activation_energy)
259) call InputErrorMsg(input,option,'activation',error_string)
260) case('PREFACTOR_SPECIES')
261) error_string = 'CHEMISTRY,MINERAL_KINETICS,PREFACTOR,&
262) &SPECIES'
263) prefactor_species => TSPrefactorSpeciesCreate()
264) call InputReadWord(input,option,prefactor_species%name, &
265) PETSC_TRUE)
266) call InputErrorMsg(input,option,'name',error_string)
267) do
268) call InputReadPflotranString(input,option)
269) call InputReadStringErrorMsg(input,option,card)
270) if (InputCheckExit(input,option)) exit
271) call InputReadWord(input,option,word,PETSC_TRUE)
272) call InputErrorMsg(input,option,'keyword',error_string)
273) select case(trim(word))
274) case('ALPHA')
275) call InputReadDouble(input,option, &
276) prefactor_species%alpha)
277) call InputErrorMsg(input,option,'alpha',error_string)
278) case('BETA')
279) call InputReadDouble(input,option, &
280) prefactor_species%beta)
281) call InputErrorMsg(input,option,'beta',error_string)
282) case('ATTENUATION_COEF')
283) call InputReadDouble(input,option, &
284) prefactor_species%attenuation_coef)
285) call InputErrorMsg(input,option, &
286) 'attenuation coefficient', &
287) error_string)
288) case default
289) call InputKeywordUnrecognized(word, &
290) 'CHEMISTRY,MINERAL_KINETICS,PREFACTOR,SPECIES', &
291) option)
292) end select
293) enddo
294) ! add prefactor species
295) if (.not.associated(prefactor%species)) then
296) prefactor%species => prefactor_species
297) else ! append to end of list
298) cur_prefactor_species => prefactor%species
299) do
300) if (.not.associated(cur_prefactor_species%next)) then
301) cur_prefactor_species%next => prefactor_species
302) exit
303) else
304) cur_prefactor_species => cur_prefactor_species%next
305) endif
306) enddo
307) endif
308) error_string = 'CHEMISTRY,MINERAL_KINETICS,PREFACTOR'
309) case default
310) call InputKeywordUnrecognized(word, &
311) 'CHEMISTRY,MINERAL_KINETICS,PREFACTOR',option)
312) end select
313) enddo
314) ! add prefactor
315) if (.not.associated(tstrxn%prefactor)) then
316) tstrxn%prefactor => prefactor
317) else ! append to end of list
318) cur_prefactor => tstrxn%prefactor
319) do
320) if (.not.associated(cur_prefactor%next)) then
321) cur_prefactor%next => prefactor
322) exit
323) else
324) cur_prefactor => cur_prefactor%next
325) endif
326) enddo
327) endif
328) error_string = 'CHEMISTRY,MINERAL_KINETICS'
329) case default
330) call InputKeywordUnrecognized(word, &
331) 'CHEMISTRY,MINERAL_KINETICS',option)
332) end select
333) enddo
334) ! Loop over prefactors and set kinetic rates and activation energies
335) ! equal to the "outer" values if zero.
336) cur_prefactor => tstrxn%prefactor
337) do
338) if (.not.associated(cur_prefactor)) exit
339) ! if not initialized
340) if (Uninitialized(cur_prefactor%rate)) then
341) cur_prefactor%rate = tstrxn%rate
342) if (Uninitialized(cur_prefactor%rate)) then
343) option%io_buffer = 'Both outer and inner prefactor rate ' // &
344) 'constants uninitialized for kinetic mineral ' // &
345) cur_mineral%name // '.'
346) call printErrMsg(option)
347) endif
348) endif
349) if (Uninitialized(cur_prefactor%activation_energy)) then
350) cur_prefactor%activation_energy = tstrxn%activation_energy
351) endif
352) cur_prefactor => cur_prefactor%next
353) enddo
354) ! add tst rxn
355) if (.not.associated(cur_mineral%tstrxn)) then
356) cur_mineral%tstrxn => tstrxn
357) else ! append to end of list
358) cur_tstrxn => cur_mineral%tstrxn
359) do
360) if (.not.associated(cur_tstrxn%next)) then
361) cur_tstrxn%next => tstrxn
362) exit
363) else
364) cur_tstrxn => cur_tstrxn%next
365) endif
366) enddo
367) endif
368) cur_mineral%id = abs(cur_mineral%id)
369) exit
370) endif
371) cur_mineral => cur_mineral%next
372) enddo
373) if (.not.found) then
374) option%io_buffer = 'Mineral "' // trim(name) // '" specified under ' // &
375) 'CHEMISTRY,MINERAL_KINETICS not found in list of available minerals.'
376) call printErrMsg(option)
377) endif
378) enddo
379)
380) cur_mineral => mineral%mineral_list
381) imnrl = 0
382) do
383) if (.not.associated(cur_mineral)) exit
384) if (cur_mineral%id < 0 .and. &
385) cur_mineral%itype == MINERAL_KINETIC) then
386) option%io_buffer = 'No rate provided in input file for mineral: ' // &
387) trim(cur_mineral%name) // '.'
388) call printErrMsg(option)
389) endif
390) if (associated(cur_mineral%tstrxn)) then
391) imnrl = imnrl + 1
392) !geh reaction%kinmnrl_names(imnrl) = cur_mineral%name
393) endif
394) cur_mineral => cur_mineral%next
395) enddo
396)
397) cur_mineral => mineral%mineral_list
398) do
399) if (.not.associated(cur_mineral)) exit
400) cur_mineral%id = abs(cur_mineral%id)
401) cur_mineral => cur_mineral%next
402) enddo
403)
404) end subroutine MineralReadKinetics
405)
406) ! ************************************************************************** !
407)
408) subroutine MineralReadFromDatabase(mineral,num_dbase_temperatures,input, &
409) option)
410) !
411) ! Reads mineral from database
412) !
413) ! Author: Glenn Hammond
414) ! Date: 10/16/08
415) !
416)
417) use Input_Aux_module
418) use String_module
419) use Option_module
420) use Reaction_Database_Aux_module
421)
422) implicit none
423)
424) type(mineral_rxn_type) :: mineral
425) PetscInt :: num_dbase_temperatures
426) type(input_type), pointer :: input
427) type(option_type) :: option
428)
429) PetscInt :: ispec
430) PetscInt :: itemp
431)
432) ! read the molar volume
433) call InputReadDouble(input,option,mineral%molar_volume)
434) call InputErrorMsg(input,option,'MINERAL molar volume','DATABASE')
435) ! convert from cm^3/mol to m^3/mol
436) mineral%molar_volume = mineral%molar_volume*1.d-6
437) ! create mineral reaction
438) if (.not.associated(mineral%tstrxn)) then
439) mineral%tstrxn => TransitionStateTheoryRxnCreate()
440) endif
441) ! read the number of aqueous species in mineral rxn
442) mineral%dbaserxn => DatabaseRxnCreate()
443) call InputReadInt(input,option,mineral%dbaserxn%nspec)
444) call InputErrorMsg(input,option,'Number of species in mineral reaction', &
445) 'DATABASE')
446) ! allocate arrays for rxn
447) allocate(mineral%dbaserxn%spec_name(mineral%dbaserxn%nspec))
448) mineral%dbaserxn%spec_name = ''
449) allocate(mineral%dbaserxn%stoich(mineral%dbaserxn%nspec))
450) mineral%dbaserxn%stoich = 0.d0
451) allocate(mineral%dbaserxn%logK(num_dbase_temperatures))
452) mineral%dbaserxn%logK = 0.d0
453) ! read in species and stoichiometries
454) do ispec = 1, mineral%dbaserxn%nspec
455) call InputReadDouble(input,option,mineral%dbaserxn%stoich(ispec))
456) call InputErrorMsg(input,option,'MINERAL species stoichiometry','DATABASE')
457) call InputReadQuotedWord(input,option,mineral%dbaserxn% &
458) spec_name(ispec),PETSC_TRUE)
459) call InputErrorMsg(input,option,'MINERAL species name','DATABASE')
460) enddo
461) do itemp = 1, num_dbase_temperatures
462) call InputReadDouble(input,option,mineral%dbaserxn%logK(itemp))
463) call InputErrorMsg(input,option,'MINERAL logKs','DATABASE')
464) enddo
465) ! read the molar weight
466) call InputReadDouble(input,option,mineral%molar_weight)
467) call InputErrorMsg(input,option,'MINERAL molar weight','DATABASE')
468)
469) end subroutine MineralReadFromDatabase
470)
471) ! ************************************************************************** !
472)
473) subroutine MineralProcessConstraint(mineral,constraint_name,constraint,option)
474) !
475) ! Initializes constraints based on mineral
476) ! species in system
477) !
478) ! Author: Glenn Hammond
479) ! Date: 01/07/13
480) !
481)
482) use Option_module
483) use Input_Aux_module
484) use String_module
485) use Utility_module
486)
487) implicit none
488)
489) type(mineral_type), pointer :: mineral
490) character(len=MAXWORDLENGTH) :: constraint_name
491) type(mineral_constraint_type), pointer :: constraint
492) type(option_type) :: option
493)
494) PetscBool :: found
495) PetscInt :: imnrl, jmnrl
496)
497) character(len=MAXWORDLENGTH) :: mineral_name(mineral%nkinmnrl)
498) character(len=MAXWORDLENGTH) :: constraint_vol_frac_string(mineral%nkinmnrl)
499) character(len=MAXWORDLENGTH) :: constraint_area_string(mineral%nkinmnrl)
500) PetscReal :: constraint_vol_frac(mineral%nkinmnrl)
501) PetscReal :: constraint_area(mineral%nkinmnrl)
502) PetscBool :: external_vol_frac_dataset(mineral%nkinmnrl)
503) PetscBool :: external_area_dataset(mineral%nkinmnrl)
504)
505) if (.not.associated(constraint)) return
506)
507) mineral_name = ''
508) constraint_vol_frac_string = ''
509) constraint_area_string = ''
510) external_vol_frac_dataset = PETSC_FALSE
511) external_area_dataset = PETSC_FALSE
512) do imnrl = 1, mineral%nkinmnrl
513) found = PETSC_FALSE
514) do jmnrl = 1, mineral%nkinmnrl
515) if (StringCompare(constraint%names(imnrl), &
516) mineral%kinmnrl_names(jmnrl), &
517) MAXWORDLENGTH)) then
518) found = PETSC_TRUE
519) exit
520) endif
521) enddo
522) if (.not.found) then
523) option%io_buffer = &
524) 'Mineral ' // trim(constraint%names(imnrl)) // &
525) 'from CONSTRAINT ' // trim(constraint_name) // &
526) ' not found among kinetic minerals.'
527) call printErrMsg(option)
528) else
529) constraint_vol_frac(jmnrl) = &
530) constraint%constraint_vol_frac(imnrl)
531) constraint_area(jmnrl) = &
532) constraint%constraint_area(imnrl)
533) mineral_name(jmnrl) = constraint%names(imnrl)
534) constraint_vol_frac_string(jmnrl) = &
535) constraint%constraint_vol_frac_string(imnrl)
536) constraint_area_string(jmnrl) = &
537) constraint%constraint_area_string(imnrl)
538) external_vol_frac_dataset(jmnrl) = &
539) constraint%external_vol_frac_dataset(imnrl)
540) external_area_dataset(jmnrl) = &
541) constraint%external_area_dataset(imnrl)
542) endif
543) enddo
544) constraint%names = mineral_name
545) constraint%constraint_vol_frac = constraint_vol_frac
546) constraint%constraint_area = constraint_area
547) constraint%constraint_vol_frac_string = constraint_vol_frac_string
548) constraint%constraint_area_string = constraint_area_string
549) constraint%external_vol_frac_dataset = external_vol_frac_dataset
550) constraint%external_area_dataset = external_area_dataset
551)
552) end subroutine MineralProcessConstraint
553)
554) ! ************************************************************************** !
555)
556) subroutine RKineticMineral(Res,Jac,compute_derivative,rt_auxvar, &
557) global_auxvar,material_auxvar,reaction,option)
558) !
559) ! Computes the kinetic mineral precipitation/dissolution
560) ! rates
561) !
562) ! Author: Glenn Hammond
563) ! Date: 09/04/08
564) !
565)
566) use Option_module
567) use Material_Aux_class
568) #ifdef SOLID_SOLUTION
569) use Reaction_Solid_Soln_Aux_module
570) #endif
571)
572) implicit none
573)
574) type(option_type) :: option
575) type(reaction_type) :: reaction
576) PetscBool :: compute_derivative
577) PetscReal :: Res(reaction%ncomp)
578) PetscReal :: Jac(reaction%ncomp,reaction%ncomp)
579) type(reactive_transport_auxvar_type) :: rt_auxvar
580) type(global_auxvar_type) :: global_auxvar
581) class(material_auxvar_type) :: material_auxvar
582)
583) PetscInt :: i, j, k, imnrl, icomp, jcomp, kcplx, iphase, ncomp
584) PetscInt :: ipref, ipref_species
585) ! I am assuming a maximum of 10 prefactors and 5 species per prefactor
586) PetscReal :: tempreal, tempreal2
587) PetscReal :: affinity_factor, sign_
588) PetscReal :: Im, Im_const, dIm_dQK
589) PetscReal :: ln_conc(reaction%naqcomp)
590) PetscReal :: ln_sec(reaction%neqcplx)
591) PetscReal :: ln_act(reaction%naqcomp)
592) PetscReal :: ln_sec_act(reaction%neqcplx)
593) PetscReal :: QK, lnQK, dQK_dCj, dQK_dmj, den
594)
595) PetscReal :: ln_spec_act, spec_act_coef
596) PetscReal :: ln_prefactor, ln_numerator, ln_denominator
597) PetscReal :: prefactor(10), ln_prefactor_spec(5,10)
598) PetscReal :: sum_prefactor_rate
599) PetscReal :: dIm_dsum_prefactor_rate, dIm_dspec
600) PetscReal :: dprefactor_dprefactor_spec, dprefactor_spec_dspec
601) PetscReal :: dprefactor_spec_dspec_numerator
602) PetscReal :: dprefactor_spec_dspec_denominator
603) PetscReal :: denominator
604) PetscInt :: icplx
605) PetscReal :: ln_gam_m_beta
606)
607) #ifdef SOLID_SOLUTION
608) PetscBool :: cycle_
609) PetscReal :: rate_scale(reaction%mineral%nkinmnrl)
610) type(solid_solution_type), pointer :: cur_solid_soln
611) PetscInt :: istoich_solid
612) #endif
613)
614) type(mineral_type), pointer :: mineral
615)
616) PetscInt, parameter :: needs_to_be_fixed = 1
617)
618) PetscReal :: arrhenius_factor
619)
620) iphase = 1
621) mineral => reaction%mineral
622)
623) ln_conc = log(rt_auxvar%pri_molal)
624) ln_act = ln_conc+log(rt_auxvar%pri_act_coef)
625)
626) if (reaction%neqcplx > 0) then
627) ln_sec = log(rt_auxvar%sec_molal)
628) ln_sec_act = ln_sec+log(rt_auxvar%sec_act_coef)
629) endif
630)
631) #ifdef SOLID_SOLUTION
632) rate_scale = 1.d0
633) if (associated(reaction%solid_solution_list)) then
634) do imnrl = 1, mineral%nkinmnrl ! for each mineral
635) call RMineralRate(imnrl,ln_act,ln_sec_act,rt_auxvar,global_auxvar, &
636) QK,Im,Im_const,sum_prefactor_rate,affinity_factor, &
637) prefactor,ln_prefactor_spec,cycle_, &
638) reaction,mineral,option)
639) enddo
640)
641) cur_solid_soln => reaction%solid_solution_list
642) do
643) if (.not.associated(cur_solid_soln)) exit
644) tempreal = 0.d0
645) do istoich_solid = 1, cur_solid_soln%num_stoich_solid
646) imnrl = cur_solid_soln%stoich_solid_ids(istoich_solid)
647) ! do something with mineral ikinmnrl rate, e.g.
648) tempreal = tempreal + rt_auxvar%mnrl_rate(imnrl)
649) !tempreal = max(tempreal,dabs(rt_auxvar%mnrl_rate(ikinmnrl))
650) enddo
651) tempreal = tempreal / dble(cur_solid_soln%num_stoich_solid)
652) do istoich_solid = 1, cur_solid_soln%num_stoich_solid
653) imnrl = cur_solid_soln%stoich_solid_ids(istoich_solid)
654) rate_scale(imnrl) = tempreal
655) enddo
656) cur_solid_soln => cur_solid_soln%next
657) enddo
658) endif
659) #endif
660)
661) ! Zero all rates as default
662) rt_auxvar%mnrl_rate(:) = 0.d0
663)
664) do imnrl = 1, mineral%nkinmnrl ! for each mineral
665)
666) #ifdef SOLID_SOLUTION
667) call RMineralRate(imnrl,ln_act,ln_sec_act,rt_auxvar,global_auxvar, &
668) QK,Im,Im_const,sum_prefactor_rate,affinity_factor, &
669) prefactor,ln_prefactor_spec,cycle_, &
670) reaction,mineral,option)
671) if (cycle_) cycle
672)
673) Im = rate_scale(imnrl)*Im
674) Im_const = rate_scale(imnrl)*Im_const
675) #else
676) ! compute ion activity product
677) lnQK = -mineral%kinmnrl_logK(imnrl)*LOG_TO_LN
678)
679) ! activity of water
680) if (mineral%kinmnrlh2oid(imnrl) > 0) then
681) lnQK = lnQK + mineral%kinmnrlh2ostoich(imnrl)* &
682) rt_auxvar%ln_act_h2o
683) endif
684)
685) ncomp = mineral%kinmnrlspecid(0,imnrl)
686) do i = 1, ncomp
687) icomp = mineral%kinmnrlspecid(i,imnrl)
688) lnQK = lnQK + mineral%kinmnrlstoich(i,imnrl)*ln_act(icomp)
689) enddo
690)
691) if (lnQK <= 6.90776d0) then
692) QK = exp(lnQK)
693) else
694) QK = 1.d3
695) endif
696)
697) if (associated(mineral%kinmnrl_Temkin_const)) then
698) if (associated(mineral%kinmnrl_min_scale_factor)) then
699) affinity_factor = 1.d0-QK**(1.d0/ &
700) (mineral%kinmnrl_min_scale_factor(imnrl)* &
701) mineral%kinmnrl_Temkin_const(imnrl)))
702) else
703) affinity_factor = 1.d0-QK**(1.d0/ &
704) mineral%kinmnrl_Temkin_const(imnrl))
705) endif
706) else if (associated(mineral%kinmnrl_min_scale_factor)) then
707) affinity_factor = 1.d0-QK**(1.d0/ &
708) mineral%kinmnrl_min_scale_factor(imnrl))
709) else
710) affinity_factor = 1.d0-QK
711) endif
712)
713) sign_ = sign(1.d0,affinity_factor)
714)
715) if (rt_auxvar%mnrl_volfrac(imnrl) > 0 .or. sign_ < 0.d0) then
716)
717) ! if ((mineral%kinmnrl_irreversible(imnrl) == 0 &
718) ! .and. (rt_auxvar%mnrl_volfrac(imnrl) > 0 .or. sign_ < 0.d0)) &
719) ! .or. (mineral%kinmnrl_irreversible(imnrl) == 1 .and. sign_ < 0.d0)) then
720)
721) ! check for supersaturation threshold for precipitation
722) ! if (associated(mineral%kinmnrl_affinity_threshold)) then
723) if (mineral%kinmnrl_affinity_threshold(imnrl) > 0.d0) then
724) if (sign_ < 0.d0 .and. &
725) QK < mineral%kinmnrl_affinity_threshold(imnrl)) cycle
726) endif
727)
728) ! check for rate limiter for precipitation
729) if (mineral%kinmnrl_rate_limiter(imnrl) > 0.d0) then
730) affinity_factor = affinity_factor/(1.d0+(1.d0-affinity_factor) &
731) /mineral%kinmnrl_rate_limiter(imnrl))
732) endif
733)
734) ! compute prefactor
735) if (mineral%kinmnrl_num_prefactors(imnrl) > 0) then
736) sum_prefactor_rate = 0.d0
737) prefactor = 0.d0
738) ln_prefactor_spec = 0.d0
739) ! sum over parallel prefactors
740) do ipref = 1, mineral%kinmnrl_num_prefactors(imnrl)
741) ln_prefactor = 0.d0
742) ! product of "monod" equations
743) do ipref_species = 1, mineral%kinmnrl_prefactor_id(0,ipref,imnrl)
744) icomp = mineral%kinmnrl_prefactor_id(ipref_species,ipref,imnrl)
745) if (icomp > 0) then ! primary species
746) ln_spec_act = ln_act(icomp)
747) else ! secondary species (given a negative id to differentiate)
748) ln_spec_act = ln_sec_act(-icomp)
749) endif
750) ln_numerator = &
751) mineral%kinmnrl_pref_alpha(ipref_species,ipref,imnrl)* &
752) ln_spec_act
753) ln_denominator = log(1.d0 + &
754) exp(log(mineral%kinmnrl_pref_atten_coef(ipref_species,ipref,imnrl)) + &
755) mineral%kinmnrl_pref_beta(ipref_species,ipref,imnrl)* &
756) ln_spec_act))
757) ln_prefactor = ln_prefactor + ln_numerator
758) ln_prefactor = ln_prefactor - ln_denominator
759) ln_prefactor_spec(ipref_species,ipref) = ln_numerator - ln_denominator
760) enddo
761) prefactor(ipref) = exp(ln_prefactor)
762) ! Arrhenius factor
763) arrhenius_factor = 1.d0
764) if (mineral%kinmnrl_pref_activation_energy(ipref,imnrl) > 0.d0) then
765) arrhenius_factor = &
766) exp(mineral%kinmnrl_pref_activation_energy(ipref,imnrl)/ &
767) IDEAL_GAS_CONSTANT &
768) *(1.d0/(25.d0+273.15d0)-1.d0/(global_auxvar%temp+ &
769) 273.15d0)))
770) endif
771) sum_prefactor_rate = sum_prefactor_rate + prefactor(ipref)* &
772) mineral%kinmnrl_pref_rate(ipref,imnrl)* &
773) arrhenius_factor
774) enddo
775) else
776) ! Arrhenius factor
777) arrhenius_factor = 1.d0
778) if (mineral%kinmnrl_activation_energy(imnrl) > 0.d0) then
779) arrhenius_factor = exp(mineral%kinmnrl_activation_energy(imnrl)/ &
780) IDEAL_GAS_CONSTANT &
781) *(1.d0/(25.d0+273.15d0)-1.d0/(global_auxvar%temp+273.15d0)))
782) endif
783) sum_prefactor_rate = mineral%kinmnrl_rate_constant(imnrl)* &
784) arrhenius_factor
785) endif
786)
787) ! compute rate
788) ! rate: mol/m^2 mnrl/sec
789) ! area: m^2 mnrl/m^3 bulk
790) ! volume: m^3 bulk
791) Im_const = -rt_auxvar%mnrl_area(imnrl)
792) if (associated(mineral%kinmnrl_min_scale_factor)) then
793) Im_const = Im_const/mineral%kinmnrl_min_scale_factor(imnrl)
794) endif
795)
796) ! units: mol/sec/m^3 bulk
797) if (associated(mineral%kinmnrl_affinity_power)) then
798) ! Im_const: m^2 mnrl/m^3 bulk
799) ! sum_prefactor_rate: mol/m^2 mnrl/sec
800) Im = Im_const*sign_* &
801) abs(affinity_factor)**mineral%kinmnrl_affinity_power(imnrl)* &
802) sum_prefactor_rate
803) else
804) Im = Im_const*sign_*abs(affinity_factor)*sum_prefactor_rate
805) endif
806) ! store volumetric rate to be used in updating mineral volume fractions
807) ! at end of time step
808) rt_auxvar%mnrl_rate(imnrl) = Im ! mol/sec/m^3 bulk
809)
810) else ! rate is already zero by default; move on to next mineral
811) cycle
812) endif
813) #endif
814)
815) ! scale Im_const by volume for calculating derivatives below
816) ! units: m^2 mnrl
817) Im_const = Im_const*material_auxvar%volume
818)
819) ! convert rate from volumetric (mol/sec/m^3 bulk) to mol/sec
820) ! units: mol/sec
821) Im = Im*material_auxvar%volume
822)
823) ncomp = mineral%kinmnrlspecid(0,imnrl)
824) do i = 1, ncomp
825) icomp = mineral%kinmnrlspecid(i,imnrl)
826) Res(icomp) = Res(icomp) + mineral%kinmnrlstoich(i,imnrl)*Im
827) enddo
828)
829) if (.not. compute_derivative) cycle
830)
831) ! calculate derivatives of rate with respect to free
832) ! units = mol/sec
833) if (associated(mineral%kinmnrl_affinity_power)) then
834) dIm_dQK = -Im*mineral%kinmnrl_affinity_power(imnrl)/abs(affinity_factor)
835) else
836) dIm_dQK = -Im_const*sum_prefactor_rate
837) endif
838)
839) if (associated(mineral%kinmnrl_Temkin_const)) then
840) if (associated(mineral%kinmnrl_min_scale_factor)) then
841) dIm_dQK = dIm_dQK*(1.d0/(mineral%kinmnrl_min_scale_factor(imnrl)* &
842) mineral%kinmnrl_Temkin_const(imnrl)))/QK*(1.d0-affinity_factor)
843) else
844) dIm_dQK = dIm_dQK*(1.d0/mineral%kinmnrl_Temkin_const(imnrl))/QK* &
845) (1.d0-affinity_factor)
846) endif
847) else if (associated(mineral%kinmnrl_min_scale_factor)) then
848) dIm_dQK = dIm_dQK*(1.d0/mineral%kinmnrl_min_scale_factor(imnrl))/QK* &
849) (1.d0-affinity_factor)
850) endif
851)
852) ! derivatives with respect to primary species in reaction quotient
853) if (mineral%kinmnrl_rate_limiter(imnrl) <= 0.d0) then
854) do j = 1, ncomp
855) jcomp = mineral%kinmnrlspecid(j,imnrl)
856) ! unit = L water/mol
857) dQK_dCj = mineral%kinmnrlstoich(j,imnrl)*QK*exp(-ln_conc(jcomp))
858) ! units = (L water/mol)*(kg water/m^3 water)*(m^3 water/1000 L water) = kg water/mol
859) dQK_dmj = dQK_dCj*global_auxvar%den_kg(iphase)*1.d-3 ! the multiplication by density could be moved
860) ! outside the loop
861) do i = 1, ncomp
862) icomp = mineral%kinmnrlspecid(i,imnrl)
863) ! units = (mol/sec)*(kg water/mol) = kg water/sec
864) Jac(icomp,jcomp) = Jac(icomp,jcomp) + &
865) mineral%kinmnrlstoich(i,imnrl)*dIm_dQK*dQK_dmj
866) enddo
867) enddo
868)
869) else
870)
871) den = 1.d0+(1.d0-affinity_factor)/mineral%kinmnrl_rate_limiter(imnrl)
872) do j = 1, ncomp
873) jcomp = mineral%kinmnrlspecid(j,imnrl)
874) ! unit = L water/mol
875) dQK_dCj = mineral%kinmnrlstoich(j,imnrl)*QK*exp(-ln_conc(jcomp))
876) ! units = (L water/mol)*(kg water/m^3 water)*(m^3 water/1000 L water) = kg water/mol
877) dQK_dmj = dQK_dCj*global_auxvar%den_kg(iphase)*1.d-3 ! the multiplication by density could be moved
878) ! outside the loop
879) do i = 1, ncomp
880) icomp = mineral%kinmnrlspecid(i,imnrl)
881) ! units = (mol/sec)*(kg water/mol) = kg water/sec
882) Jac(icomp,jcomp) = Jac(icomp,jcomp) + &
883) mineral%kinmnrlstoich(i,imnrl)*dIm_dQK &
884) *(1.d0 + QK/mineral%kinmnrl_rate_limiter(imnrl)/den)*dQK_dmj/den
885) enddo
886) enddo
887) endif
888)
889) if (mineral%kinmnrl_num_prefactors(imnrl) > 0) then ! add contribution of derivative in prefactor - messy
890) #if 1
891) dIm_dsum_prefactor_rate = Im/sum_prefactor_rate
892) ! summation over parallel reactions (prefactors)
893) do ipref = 1, mineral%kinmnrl_num_prefactors(imnrl)
894) arrhenius_factor = 1.d0
895) if (mineral%kinmnrl_pref_activation_energy(ipref,imnrl) > 0.d0) then
896) arrhenius_factor = &
897) exp(mineral%kinmnrl_pref_activation_energy(ipref,imnrl)/ &
898) IDEAL_GAS_CONSTANT &
899) *(1.d0/(25.d0+273.15d0)-1.d0/(global_auxvar%temp+ &
900) 273.15d0)))
901) endif
902) ! prefactor() saved in residual calc above
903) ln_prefactor = log(prefactor(ipref))
904) ! product of "monod" equations
905) do ipref_species = 1, mineral%kinmnrl_prefactor_id(0,ipref,imnrl)
906) ! derivative of 54 with respect to a single "monod" equation
907) ! ln_prefactor_spec(,) saved in residual calc above
908) dprefactor_dprefactor_spec = exp(ln_prefactor- &
909) ln_prefactor_spec(ipref_species,ipref))
910) icomp = mineral%kinmnrl_prefactor_id(ipref_species,ipref,imnrl)
911) if (icomp > 0) then ! primary species
912) ln_spec_act = ln_act(icomp)
913) spec_act_coef = rt_auxvar%pri_act_coef(icomp)
914) else ! secondary species
915) ln_spec_act = ln_sec_act(-icomp)
916) spec_act_coef = rt_auxvar%sec_act_coef(-icomp)
917) endif
918) ! derivative of numerator in eq. 54 with respect to species activity
919) dprefactor_spec_dspec_numerator = &
920) mineral%kinmnrl_pref_alpha(ipref_species,ipref,imnrl) * &
921) exp(ln_prefactor_spec(ipref_species,ipref) - ln_spec_act)
922) ln_gam_m_beta = mineral%kinmnrl_pref_beta(ipref_species,ipref,imnrl)* &
923) ln_spec_act
924) ! denominator
925) denominator = 1.d0 + &
926) exp(log(mineral%kinmnrl_pref_atten_coef(ipref_species,ipref,imnrl)) + &
927) ln_gam_m_beta)
928) ! derivative of denominator in eq. 54 with respect to species activity
929) dprefactor_spec_dspec_denominator = -1.d0 * &
930) exp(ln_prefactor_spec(ipref_species,ipref)) / denominator * &
931) mineral%kinmnrl_pref_atten_coef(ipref_species,ipref,imnrl) * &
932) mineral%kinmnrl_pref_beta(ipref_species,ipref,imnrl) * &
933) exp(ln_gam_m_beta - ln_spec_act)
934)
935) ! chain rule for derivative of "monod" equation
936) dprefactor_spec_dspec = dprefactor_spec_dspec_numerator + &
937) dprefactor_spec_dspec_denominator
938)
939) ! thus far the derivative is with respect to the activity, convert to with
940) ! respect to molality
941) dprefactor_spec_dspec = dprefactor_spec_dspec * spec_act_coef
942)
943) dIm_dspec = dIm_dsum_prefactor_rate * dprefactor_dprefactor_spec * &
944) dprefactor_spec_dspec * &
945) mineral%kinmnrl_pref_rate(ipref,imnrl)* &
946) arrhenius_factor
947)
948)
949) if (icomp > 0) then
950) ! add derivative for primary species
951) do i = 1, ncomp
952) jcomp = mineral%kinmnrlspecid(i,imnrl)
953) ! units = (mol/sec)*(kg water/mol) = kg water/sec
954) Jac(jcomp,icomp) = Jac(jcomp,icomp) + &
955) mineral%kinmnrlstoich(i,imnrl)*dIm_dspec
956) enddo
957) else ! secondary species -- have to calculate the derivative
958) ! have to recalculate the reaction quotient (QK) for secondary species
959) icplx = -icomp
960)
961) ! compute secondary species concentration
962) lnQK = -reaction%eqcplx_logK(icplx)*LOG_TO_LN
963)
964) ! activity of water
965) if (reaction%eqcplxh2oid(icplx) > 0) then
966) lnQK = lnQK + reaction%eqcplxh2ostoich(icplx)*rt_auxvar%ln_act_h2o
967) endif
968)
969) ncomp = reaction%eqcplxspecid(0,icplx)
970) do i = 1, ncomp
971) icomp = reaction%eqcplxspecid(i,icplx)
972) lnQK = lnQK + reaction%eqcplxstoich(i,icplx)*ln_act(icomp)
973) enddo
974) ! add contribution to derivatives secondary prefactor with respect to free
975) do j = 1, ncomp
976) jcomp = reaction%eqcplxspecid(j,icplx)
977) tempreal = reaction%eqcplxstoich(j,icplx)*exp(lnQK-ln_conc(jcomp))/ &
978) rt_auxvar%sec_act_coef(icplx)
979) do i = 1, ncomp
980) icomp = reaction%eqcplxspecid(i,icplx)
981) Jac(icomp,jcomp) = Jac(icomp,jcomp) + &
982) reaction%eqcplxstoich(i,icplx)*tempreal*dIm_dspec
983) enddo
984) enddo
985) endif
986) enddo
987) enddo ! loop over prefactors
988) #endif
989) endif
990) enddo ! loop over minerals
991)
992) end subroutine RKineticMineral
993)
994) ! ************************************************************************** !
995)
996) subroutine RMineralRate(imnrl,ln_act,ln_sec_act,rt_auxvar,global_auxvar, &
997) QK,Im,Im_const,sum_prefactor_rate,affinity_factor, &
998) prefactor,ln_prefactor_spec,cycle_, &
999) reaction,mineral,option)
1000) !
1001) ! Calculates the mineral saturation index
1002) !
1003) ! Author: Glenn Hammond
1004) ! Date: 08/29/11
1005) !
1006)
1007) use Option_module
1008)
1009) implicit none
1010)
1011) type(option_type) :: option
1012) type(reaction_type) :: reaction
1013) type(mineral_type) :: mineral
1014) type(reactive_transport_auxvar_type) :: rt_auxvar
1015) type(global_auxvar_type) :: global_auxvar
1016) PetscReal :: ln_act(reaction%naqcomp)
1017) PetscReal :: ln_sec_act(reaction%neqcplx)
1018) PetscReal :: QK
1019) PetscReal :: Im, Im_const
1020) PetscReal :: sum_prefactor_rate
1021) PetscReal :: affinity_factor
1022) PetscReal :: prefactor(10), ln_prefactor_spec(5,10)
1023) PetscBool :: cycle_
1024)
1025) PetscReal :: lnQK
1026) PetscInt :: i, imnrl, icomp, ncomp, ipref, ipref_species
1027) PetscReal :: sign_
1028)
1029) PetscReal :: ln_spec_act
1030) PetscReal :: ln_prefactor, ln_numerator, ln_denominator
1031)
1032) PetscReal :: arrhenius_factor
1033) PetscInt, parameter :: iphase = 1
1034)
1035) cycle_ = PETSC_FALSE
1036)
1037) ! compute ion activity product
1038) lnQK = -mineral%kinmnrl_logK(imnrl)*LOG_TO_LN
1039)
1040) ! activity of water
1041) if (mineral%kinmnrlh2oid(imnrl) > 0) then
1042) lnQK = lnQK + mineral%kinmnrlh2ostoich(imnrl)* &
1043) rt_auxvar%ln_act_h2o
1044) endif
1045)
1046) ncomp = mineral%kinmnrlspecid(0,imnrl)
1047) do i = 1, ncomp
1048) icomp = mineral%kinmnrlspecid(i,imnrl)
1049) lnQK = lnQK + mineral%kinmnrlstoich(i,imnrl)*ln_act(icomp)
1050) enddo
1051)
1052) if (lnQK <= 6.90776d0) then
1053) QK = exp(lnQK)
1054) else
1055) QK = 1.d3
1056) endif
1057)
1058) if (associated(mineral%kinmnrl_Temkin_const)) then
1059) affinity_factor = 1.d0-QK**(1.d0/ &
1060) mineral%kinmnrl_Temkin_const(imnrl))
1061) else
1062) affinity_factor = 1.d0-QK
1063) endif
1064)
1065) sign_ = sign(1.d0,affinity_factor)
1066)
1067) if (rt_auxvar%mnrl_volfrac(imnrl) > 0 .or. sign_ < 0.d0) then
1068)
1069) ! if ((mineral%kinmnrl_irreversible(imnrl) == 0 &
1070) ! .and. (rt_auxvar%mnrl_volfrac(imnrl) > 0 .or. sign_ < 0.d0)) &
1071) ! .or. (mineral%kinmnrl_irreversible(imnrl) == 1 .and. sign_ < 0.d0)) then
1072)
1073) ! check for supersaturation threshold for precipitation
1074) ! if (associated(mineral%kinmnrl_affinity_threshold)) then
1075) if (mineral%kinmnrl_affinity_threshold(imnrl) > 0.d0) then
1076) if (sign_ < 0.d0 .and. &
1077) QK < mineral%kinmnrl_affinity_threshold(imnrl)) then
1078) cycle_ = PETSC_TRUE
1079) return
1080) endif
1081) endif
1082)
1083) ! check for rate limiter for precipitation
1084) if (mineral%kinmnrl_rate_limiter(imnrl) > 0.d0) then
1085) affinity_factor = affinity_factor/(1.d0+(1.d0-affinity_factor) &
1086) /mineral%kinmnrl_rate_limiter(imnrl))
1087) endif
1088)
1089) ! compute prefactor
1090) if (mineral%kinmnrl_num_prefactors(imnrl) > 0) then
1091) sum_prefactor_rate = 0.d0
1092) prefactor = 0.d0
1093) ln_prefactor_spec = 0.d0
1094) ! sum over parallel prefactors
1095) do ipref = 1, mineral%kinmnrl_num_prefactors(imnrl)
1096) ln_prefactor = 0.d0
1097) ! product of "monod" equations
1098) do ipref_species = 1, mineral%kinmnrl_prefactor_id(0,ipref,imnrl)
1099) icomp = mineral%kinmnrl_prefactor_id(ipref_species,ipref,imnrl)
1100) if (icomp > 0) then ! primary species
1101) ln_spec_act = ln_act(icomp)
1102) else ! secondary species (given a negative id to differentiate)
1103) ln_spec_act = ln_sec_act(-icomp)
1104) endif
1105) ln_numerator = &
1106) mineral%kinmnrl_pref_alpha(ipref_species,ipref,imnrl)* &
1107) ln_spec_act
1108) ln_denominator = log(1.d0 + &
1109) exp(log(mineral%kinmnrl_pref_atten_coef(ipref_species,ipref,imnrl)) + &
1110) mineral%kinmnrl_pref_beta(ipref_species,ipref,imnrl)* &
1111) ln_spec_act))
1112) ln_prefactor = ln_prefactor + ln_numerator
1113) ln_prefactor = ln_prefactor - ln_denominator
1114) ln_prefactor_spec(ipref_species,ipref) = ln_numerator - ln_denominator
1115) enddo
1116) prefactor(ipref) = exp(ln_prefactor)
1117) ! Arrhenius factor
1118) arrhenius_factor = 1.d0
1119) if (mineral%kinmnrl_pref_activation_energy(ipref,imnrl) > 0.d0) then
1120) arrhenius_factor = &
1121) exp(mineral%kinmnrl_pref_activation_energy(ipref,imnrl)/ &
1122) IDEAL_GAS_CONSTANT &
1123) *(1.d0/(25.d0+273.15d0)-1.d0/(global_auxvar%temp+ &
1124) 273.15d0)))
1125) endif
1126) sum_prefactor_rate = sum_prefactor_rate + prefactor(ipref)* &
1127) mineral%kinmnrl_pref_rate(ipref,imnrl)* &
1128) arrhenius_factor
1129) enddo
1130) else
1131) ! Arrhenius factor
1132) arrhenius_factor = 1.d0
1133) if (mineral%kinmnrl_activation_energy(imnrl) > 0.d0) then
1134) arrhenius_factor = exp(mineral%kinmnrl_activation_energy(imnrl)/ &
1135) IDEAL_GAS_CONSTANT &
1136) *(1.d0/(25.d0+273.15d0)-1.d0/(global_auxvar%temp+273.15d0)))
1137) endif
1138) sum_prefactor_rate = mineral%kinmnrl_rate_constant(imnrl)* &
1139) arrhenius_factor
1140) endif
1141)
1142) ! compute rate
1143) ! rate: mol/m^2 mnrl/sec
1144) ! area: m^2 mnrl/m^3 bulk
1145) ! volume: m^3 bulk
1146) Im_const = -rt_auxvar%mnrl_area(imnrl)
1147)
1148) ! units: mol/sec/m^3 bulk
1149) if (associated(mineral%kinmnrl_affinity_power)) then
1150) ! Im_const: m^2 mnrl/m^3 bulk
1151) ! sum_prefactor_rate: mol/m^2 mnrl/sec
1152) Im = Im_const*sign_* &
1153) abs(affinity_factor)**mineral%kinmnrl_affinity_power(imnrl)* &
1154) sum_prefactor_rate
1155) else
1156) Im = Im_const*sign_*abs(affinity_factor)*sum_prefactor_rate
1157) endif
1158) ! store volumetric rate to be used in updating mineral volume fractions
1159) ! at end of time step
1160) rt_auxvar%mnrl_rate(imnrl) = Im ! mol/sec/m^3 bulk
1161) else ! rate is already zero by default; move on to next mineral
1162) cycle_ = PETSC_TRUE
1163) endif
1164)
1165) end subroutine RMineralRate
1166)
1167) ! ************************************************************************** !
1168)
1169) function RMineralSaturationIndex(imnrl,rt_auxvar,global_auxvar,reaction,option)
1170) !
1171) ! Calculates the mineral saturation index
1172) !
1173) ! Author: Glenn Hammond
1174) ! Date: 08/29/11
1175) !
1176)
1177) use Option_module
1178)
1179) implicit none
1180)
1181) type(option_type) :: option
1182) PetscInt :: imnrl
1183) type(reaction_type) :: reaction
1184) type(reactive_transport_auxvar_type) :: rt_auxvar
1185) type(global_auxvar_type) :: global_auxvar
1186)
1187) PetscReal :: RMineralSaturationIndex
1188) PetscInt :: i, icomp
1189) PetscReal :: lnQK
1190) PetscInt, parameter :: iphase = 1
1191) type(mineral_type), pointer :: mineral
1192)
1193) mineral => reaction%mineral
1194)
1195) if (.not.option%use_isothermal) then
1196) call MineralUpdateTempDepCoefs(global_auxvar%temp, &
1197) global_auxvar%pres(iphase), &
1198) reaction%mineral, &
1199) reaction%use_geothermal_hpt, &
1200) PETSC_TRUE,option)
1201) endif
1202)
1203) ! compute saturation
1204) lnQK = -mineral%mnrl_logK(imnrl)*LOG_TO_LN
1205) if (mineral%mnrlh2oid(imnrl) > 0) then
1206) lnQK = lnQK + mineral%mnrlh2ostoich(imnrl)*rt_auxvar%ln_act_h2o
1207) endif
1208) do i = 1, mineral%mnrlspecid(0,imnrl)
1209) icomp = mineral%mnrlspecid(i,imnrl)
1210) lnQK = lnQK + mineral%mnrlstoich(i,imnrl)* &
1211) log(rt_auxvar%pri_molal(icomp)*rt_auxvar%pri_act_coef(icomp))
1212) enddo
1213) RMineralSaturationIndex = exp(lnQK)
1214)
1215) end function RMineralSaturationIndex
1216)
1217) ! ************************************************************************** !
1218)
1219) subroutine MineralUpdateTempDepCoefs(temp,pres,mineral,use_geothermal_hpt, &
1220) update_mnrl,option)
1221) !
1222) ! Updates temperature dependent coefficients for
1223) ! anisothermal simulations
1224) !
1225) ! Author: Glenn Hammond
1226) ! Date: 01/25/13
1227) !
1228)
1229) use Option_module
1230)
1231) implicit none
1232)
1233) PetscReal :: temp
1234) PetscReal :: pres
1235) type(mineral_type) :: mineral
1236) PetscBool :: use_geothermal_hpt
1237) PetscBool :: update_mnrl
1238) type(option_type) :: option
1239)
1240) if (.not.use_geothermal_hpt) then
1241) if (associated(mineral%kinmnrl_logKcoef)) then
1242) call ReactionInterpolateLogK(mineral%kinmnrl_logKcoef, &
1243) mineral%kinmnrl_logK, &
1244) temp, &
1245) mineral%nkinmnrl)
1246) endif
1247) if (update_mnrl .and. associated(mineral%mnrl_logKcoef)) then
1248) call ReactionInterpolateLogK(mineral%mnrl_logKcoef, &
1249) mineral%mnrl_logK, &
1250) temp, &
1251) mineral%nmnrl)
1252) endif
1253) else
1254) if (associated(mineral%kinmnrl_logKcoef)) then
1255) call ReactionInterpolateLogK_hpt(mineral%kinmnrl_logKcoef, &
1256) mineral%kinmnrl_logK, &
1257) temp, &
1258) pres, &
1259) mineral%nkinmnrl)
1260) endif
1261) if (update_mnrl .and. associated(mineral%mnrl_logKcoef)) then
1262) call ReactionInterpolateLogK_hpt(mineral%mnrl_logKcoef, &
1263) mineral%mnrl_logK, &
1264) temp, &
1265) pres, &
1266) mineral%nmnrl)
1267) endif
1268) endif
1269)
1270) end subroutine MineralUpdateTempDepCoefs
1271)
1272) end module Reaction_Mineral_module