reaction_database_aux.F90 coverage: 100.00 %func 80.57 %block
1) module Reaction_Database_Aux_module
2)
3) use PFLOTRAN_Constants_module
4)
5) implicit none
6)
7) private
8)
9) #include "petsc/finclude/petscsys.h"
10)
11) type, public :: database_rxn_type
12) PetscInt :: nspec
13) character(len=MAXWORDLENGTH), pointer :: spec_name(:)
14) PetscReal, pointer :: stoich(:)
15) PetscInt, pointer :: spec_ids(:)
16) PetscReal, pointer :: logK(:)
17) PetscReal, pointer :: logKCoeff_hpt(:)
18) end type database_rxn_type
19)
20) public :: BasisAlignSpeciesInRxn, &
21) BasisSubSpeciesInGasOrSecRxn, &
22) BasisSubSpeciesinMineralRxn, &
23) DatabaseRxnCreate, &
24) DatabaseRxnCreateFromRxnString, &
25) DatabaseCheckLegitimateLogKs, &
26) DatabaseRxnDestroy
27)
28) contains
29)
30) ! ************************************************************************** !
31)
32) function DatabaseRxnCreate()
33) !
34) ! Allocate and initialize an equilibrium reaction
35) !
36) ! Author: Glenn Hammond
37) ! Date: 09/01/08
38) !
39)
40) implicit none
41)
42) type(database_rxn_type), pointer :: DatabaseRxnCreate
43)
44) type(database_rxn_type), pointer :: dbaserxn
45)
46) allocate(dbaserxn)
47) dbaserxn%nspec = 0
48) nullify(dbaserxn%spec_name)
49) nullify(dbaserxn%stoich)
50) nullify(dbaserxn%spec_ids)
51) nullify(dbaserxn%logK)
52)
53) DatabaseRxnCreate => dbaserxn
54)
55) end function DatabaseRxnCreate
56)
57) ! ************************************************************************** !
58)
59) function DatabaseRxnCreateFromRxnString(reaction_string, &
60) naqcomp, aq_offset, &
61) primary_aq_species_names, &
62) nimcomp, im_offset, &
63) primary_im_species_names, &
64) consider_immobile_species,&
65) option)
66) !
67) ! Creates a database reaction given a
68) ! reaction string
69) !
70) ! Author: Glenn Hammond
71) ! Date: 10/30/12
72) !
73)
74) use Option_module
75) use String_module
76) use Input_Aux_module
77)
78) implicit none
79)
80) character(len=MAXSTRINGLENGTH) :: reaction_string
81) PetscInt :: naqcomp ! mobile aqueoues species
82) PetscInt :: aq_offset ! offset for aqueous species
83) character(len=MAXWORDLENGTH) :: primary_aq_species_names(naqcomp)
84) PetscInt :: nimcomp ! immobile primary speces (e.g. biomass)
85) PetscInt :: im_offset ! offset for aqueous species
86) character(len=MAXWORDLENGTH), pointer :: primary_im_species_names(:)
87) PetscBool :: consider_immobile_species
88) type(option_type) :: option
89)
90) type(database_rxn_type), pointer :: DatabaseRxnCreateFromRxnString
91)
92) character(len=MAXSTRINGLENGTH) :: string, string2
93) character(len=MAXWORDLENGTH) :: word, word2
94) PetscInt :: icount
95) PetscInt :: midpoint
96) PetscInt :: i, j, idum
97) PetscReal :: value
98) PetscBool :: negative_flag
99) PetscBool :: found
100) PetscErrorCode :: ierr
101) type(database_rxn_type), pointer :: dbaserxn
102)
103)
104) dbaserxn => DatabaseRxnCreate()
105)
106) icount = 0
107) ! Be sure to copy as words are removed when read. Need to full string for
108) ! later below
109) string = reaction_string
110) do
111) ierr = 0
112) call InputReadWord(string,word,PETSC_TRUE,ierr)
113) if (InputError(ierr)) exit
114)
115) select case(word)
116) case('+')
117) case('-')
118) case('=','<=>','<->')
119) case default
120) ! try reading as double precision
121) string2 = word
122) if (.not.StringStartsWithAlpha(string2)) then
123) ! the word is the stoichiometry value
124) else
125) ! check water
126) word2 = 'H2O'
127) if (.not.StringCompareIgnoreCase(word,word2)) then
128) ! the word is the species name
129) icount = icount + 1
130) endif
131) endif
132) end select
133)
134) enddo
135)
136) ! load species into database format
137) dbaserxn%nspec = icount
138) allocate(dbaserxn%spec_name(icount))
139) dbaserxn%spec_name = ''
140) allocate(dbaserxn%stoich(icount))
141) dbaserxn%stoich = UNINITIALIZED_DOUBLE
142) allocate(dbaserxn%spec_ids(icount))
143) dbaserxn%spec_ids = 0
144)
145) string = reaction_string
146) icount = 1
147) ! midpoint points to the first product species, as in
148) ! reactant1 + reactant2 <-> product1 + product2
149) midpoint = 0
150) negative_flag = PETSC_FALSE
151) do
152) !geh: This conditional ensures that if water is at the end of
153) ! the reaction expression, it is skipped.
154) if (icount > dbaserxn%nspec) exit
155)
156) ierr = 0
157) call InputReadWord(string,word,PETSC_TRUE,ierr)
158) if (InputError(ierr)) exit
159)
160) select case(word)
161) case('+')
162) case('-')
163) ! toggle negative flag
164) if (negative_flag) then
165) negative_flag = PETSC_FALSE
166) else
167) negative_flag = PETSC_TRUE
168) endif
169) case('=','<=>','<->')
170) midpoint = icount
171) case default
172) ! try reading as double precision
173) string2 = word
174) if (.not.StringStartsWithAlpha(string2)) then
175) ! negate if a product
176) call InputReadDouble(string2,option,value,ierr)
177) if (ierr /= 0) then
178) option%io_buffer = 'Keyword "' // trim(word) // &
179) '" not recognized in reaction string "' // &
180) trim(reaction_string) // '".'
181) call printErrMsg(option)
182) endif
183) ! negate if negative stoichiometry
184) if (negative_flag) value = -1.0*value
185) dbaserxn%stoich(icount) = value
186) else
187) dbaserxn%spec_name(icount) = word
188) if (negative_flag .and. &
189) (dbaserxn%stoich(icount) + 999.d0) < 1.d-10) then
190) dbaserxn%stoich(icount) = -1.d0
191) endif
192)
193) ! set the primary aqueous species id
194) found = PETSC_FALSE
195) do i = 1, naqcomp
196) if (StringCompare(word,primary_aq_species_names(i), &
197) MAXWORDLENGTH)) then
198) dbaserxn%spec_ids(icount) = i + aq_offset
199) found = PETSC_TRUE
200) exit
201) endif
202) enddo
203) ! set the primary immobile species id
204) if (.not.found .and. consider_immobile_species) then
205) do i = 1, nimcomp
206) if (StringCompare(word,primary_im_species_names(i), &
207) MAXWORDLENGTH)) then
208) dbaserxn%spec_ids(icount) = i + im_offset
209) found = PETSC_TRUE
210) exit
211) endif
212) enddo
213) endif
214)
215) ! check water
216) word2 = 'H2O'
217) if (StringCompareIgnoreCase(word,word2)) then
218) ! set stoichiometry back to uninitialized
219) dbaserxn%stoich(icount) = UNINITIALIZED_DOUBLE
220) ! don't increment icount
221) else if (.not.found) then
222) if (consider_immobile_species) then
223) option%io_buffer = 'Species ' // trim(word) // &
224) ' in reaction not found among primary species list.'
225) else
226) option%io_buffer = 'Species ' // trim(word) // &
227) ' in reaction not found among primary aqueous species list.'
228) endif
229) call printErrMsg(option)
230) else
231) icount = icount + 1
232) endif
233) endif
234) negative_flag = PETSC_FALSE
235) end select
236) enddo
237)
238) ! if no stoichiometry specified, default = 1.
239) do i = 1, dbaserxn%nspec
240) if ((dbaserxn%stoich(i) + 999.d0) < 1.d-10) dbaserxn%stoich(i) = 1.d0
241) enddo
242) if (midpoint > 0) then
243) ! negate stoichiometries after midpoint
244) do i = midpoint, dbaserxn%nspec
245) dbaserxn%stoich(i) = -1.d0*dbaserxn%stoich(i)
246) enddo
247) endif
248) ! now negate all stoichiometries to have - for reactants; + for products
249) do i = 1, dbaserxn%nspec
250) dbaserxn%stoich(i) = -1.d0*dbaserxn%stoich(i)
251) enddo
252) ! reorder species ids in ascending order
253) do i = 1, dbaserxn%nspec
254) do j = i+1, dbaserxn%nspec
255) if (dbaserxn%spec_ids(i) > dbaserxn%spec_ids(j)) then
256) ! swap ids
257) idum = dbaserxn%spec_ids(j)
258) dbaserxn%spec_ids(j) = dbaserxn%spec_ids(i)
259) dbaserxn%spec_ids(i) = idum
260) ! swap stoichiometry
261) value = dbaserxn%stoich(j)
262) dbaserxn%stoich(j) = dbaserxn%stoich(i)
263) dbaserxn%stoich(i) = value
264) ! swap names
265) word = dbaserxn%spec_name(j)
266) dbaserxn%spec_name(j) = dbaserxn%spec_name(i)
267) dbaserxn%spec_name(i) = word
268) endif
269) enddo
270) enddo
271)
272) DatabaseRxnCreateFromRxnString => dbaserxn
273)
274) end function DatabaseRxnCreateFromRxnString
275)
276) ! ************************************************************************** !
277)
278) subroutine BasisAlignSpeciesInRxn(num_basis_species,basis_names, &
279) num_rxn_species,rxn_species_names, &
280) rxn_stoich,rxn_species_ids,species_name, &
281) option)
282) !
283) ! Aligns the ordering of species in reaction with
284) ! the current basis
285) !
286) ! Author: Glenn Hammond
287) ! Date: 10/07/08
288) !
289)
290) use Option_module
291) use String_module
292)
293) implicit none
294)
295) PetscInt :: num_basis_species
296) character(len=MAXWORDLENGTH) :: basis_names(num_basis_species), species_name
297) PetscInt :: num_rxn_species
298) character(len=MAXWORDLENGTH) :: rxn_species_names(num_rxn_species)
299) PetscReal :: rxn_stoich(num_rxn_species)
300) PetscInt :: rxn_species_ids(num_rxn_species)
301) type(option_type) :: option
302)
303) PetscInt :: i_rxn_species
304) PetscInt :: i_basis_species
305) PetscReal :: stoich_new(num_basis_species)
306) PetscBool :: found
307)
308) stoich_new = 0.d0
309) do i_rxn_species = 1, num_rxn_species
310) found = PETSC_FALSE
311) do i_basis_species = 1, num_basis_species
312) if (StringCompare(rxn_species_names(i_rxn_species), &
313) basis_names(i_basis_species), &
314) MAXWORDLENGTH)) then
315) stoich_new(i_basis_species) = rxn_stoich(i_rxn_species)
316) found = PETSC_TRUE
317) exit
318) endif
319) enddo
320) if (.not.found) then
321) option%io_buffer = trim(rxn_species_names(i_rxn_species)) // &
322) ' not found in basis (BasisAlignSpeciesInRxn) for species ' // &
323) trim(species_name)
324) call printErrMsg(option)
325) endif
326) enddo
327)
328) ! zero everthing out
329) rxn_species_names = ''
330) rxn_stoich = 0.d0
331) rxn_species_ids = 0
332)
333) ! fill in
334) i_rxn_species = 0
335) do i_basis_species = 1, num_basis_species
336) if (dabs(stoich_new(i_basis_species)) > 1.d-40) then
337) i_rxn_species = i_rxn_species + 1
338) rxn_species_names(i_rxn_species) = basis_names(i_basis_species)
339) rxn_stoich(i_rxn_species) = stoich_new(i_basis_species)
340) rxn_species_ids(i_rxn_species) = i_basis_species
341) endif
342) enddo
343)
344) if (i_rxn_species /= num_rxn_species) then
345) write(option%io_buffer,*) &
346) 'Number of reaction species does not match original:', &
347) i_rxn_species, num_rxn_species
348) call printErrMsg(option)
349) endif
350)
351) end subroutine BasisAlignSpeciesInRxn
352)
353) ! ************************************************************************** !
354)
355) subroutine BasisSubSpeciesInGasOrSecRxn(name1,dbaserxn1,dbaserxn2,scale)
356) !
357) ! Swaps out a chemical species in a chemical
358) ! reaction, replacing it with the species in a
359) ! secondary reaction (swaps 1 into 2)
360) !
361) ! Author: Glenn Hammond
362) ! Date: 10/06/08
363) !
364)
365) use String_module
366)
367) implicit none
368)
369) character(len=MAXWORDLENGTH) :: name1
370) type(database_rxn_type) :: dbaserxn1
371) type(database_rxn_type) :: dbaserxn2
372) PetscReal :: scale
373)
374) PetscInt :: i, j, tempcount, prevcount
375) character(len=MAXWORDLENGTH) :: tempnames(20)
376) PetscReal :: tempstoich(20)
377) PetscBool :: found
378)
379) tempnames = ''
380) tempstoich = 0.d0
381)
382) ! load species in reaction other than species 1 into new arrays
383) scale = 1.d0
384) tempcount = 0
385) do i=1,dbaserxn2%nspec
386) if (.not.StringCompare(name1, &
387) dbaserxn2%spec_name(i), &
388) MAXWORDLENGTH)) then
389) tempcount = tempcount + 1
390) tempnames(tempcount) = dbaserxn2%spec_name(i)
391) tempstoich(tempcount) = dbaserxn2%stoich(i)
392) else
393) scale = dbaserxn2%stoich(i)
394) endif
395) enddo
396)
397) ! search for duplicate species and add stoichs or add new species
398) ! if not duplicated
399) do j=1,dbaserxn1%nspec
400) found = PETSC_FALSE
401) do i=1,tempcount
402) if (StringCompare(tempnames(i), &
403) dbaserxn1%spec_name(j), &
404) MAXWORDLENGTH)) then
405) tempstoich(i) = tempstoich(i) + scale*dbaserxn1%stoich(j)
406) found = PETSC_TRUE
407) exit
408) endif
409) enddo
410) if (.not.found) then
411) tempcount = tempcount + 1
412) tempnames(tempcount) = dbaserxn1%spec_name(j)
413) tempstoich(tempcount) = scale*dbaserxn1%stoich(j)
414) endif
415) enddo
416)
417) ! deallocate arrays
418) deallocate(dbaserxn2%spec_name)
419) deallocate(dbaserxn2%stoich)
420)
421) ! check for zero stoichiometries due to cancelation
422) prevcount = tempcount
423) tempcount = 0
424) do i=1,prevcount
425) if (dabs(tempstoich(i)) > 1.d-10) then
426) tempcount = tempcount + 1
427) tempnames(tempcount) = tempnames(i)
428) tempstoich(tempcount) = tempstoich(i)
429) endif
430) enddo
431)
432) tempnames(tempcount+1:) = ''
433) tempstoich(tempcount+1:) = 0.d0
434)
435) ! reallocate
436) allocate(dbaserxn2%spec_name(tempcount))
437) allocate(dbaserxn2%stoich(tempcount))
438)
439) ! fill arrays in dbaserxn
440) dbaserxn2%nspec = tempcount
441) do i=1,tempcount
442) dbaserxn2%spec_name(i) = tempnames(i)
443) dbaserxn2%stoich(i) = tempstoich(i)
444) enddo
445)
446) dbaserxn2%logK = dbaserxn2%logK + scale*dbaserxn1%logK
447)
448) end subroutine BasisSubSpeciesInGasOrSecRxn
449)
450) ! ************************************************************************** !
451)
452) subroutine BasisSubSpeciesInMineralRxn(name,sec_dbaserxn,mnrl_dbaserxn,scale)
453) !
454) ! Swaps out a chemical species in a chemical
455) ! reaction, replacing it with the species in a
456) ! secondary reaction (swaps 1 into 2)
457) !
458) ! Author: Glenn Hammond
459) ! Date: 10/06/08
460) !
461)
462) use String_module
463)
464) implicit none
465)
466) character(len=MAXWORDLENGTH) :: name
467) type(database_rxn_type) :: sec_dbaserxn
468) type(database_rxn_type) :: mnrl_dbaserxn
469) PetscReal :: scale
470)
471) PetscInt :: i, j, tempcount, prevcount
472) character(len=MAXWORDLENGTH) :: tempnames(20)
473) PetscReal :: tempstoich(20)
474) PetscBool :: found
475)
476) tempnames = ''
477) tempstoich = 0.d0
478)
479) ! load species in reaction other than species 1 into new arrays
480) scale = 1.d0
481) tempcount = 0
482) do i=1,mnrl_dbaserxn%nspec
483) if (.not.StringCompare(name, &
484) mnrl_dbaserxn%spec_name(i), &
485) MAXWORDLENGTH)) then
486) tempcount = tempcount + 1
487) tempnames(tempcount) = mnrl_dbaserxn%spec_name(i)
488) tempstoich(tempcount) = mnrl_dbaserxn%stoich(i)
489) else
490) scale = mnrl_dbaserxn%stoich(i)
491) endif
492) enddo
493)
494) ! search for duplicate species and add stoichs or add new species
495) ! if not duplicated
496) do j=1,sec_dbaserxn%nspec
497) found = PETSC_FALSE
498) do i=1,mnrl_dbaserxn%nspec
499) if (StringCompare(tempnames(i), &
500) sec_dbaserxn%spec_name(j), &
501) MAXWORDLENGTH)) then
502) tempstoich(i) = tempstoich(i) + scale*sec_dbaserxn%stoich(j)
503) found = PETSC_TRUE
504) exit
505) endif
506) enddo
507) if (.not.found) then
508) tempcount = tempcount + 1
509) tempnames(tempcount) = sec_dbaserxn%spec_name(j)
510) tempstoich(tempcount) = scale*sec_dbaserxn%stoich(j)
511) endif
512) enddo
513)
514) ! deallocate arrays
515) deallocate(mnrl_dbaserxn%spec_name)
516) deallocate(mnrl_dbaserxn%stoich)
517)
518) ! check for zero stoichiometries due to cancelation
519) prevcount = tempcount
520) tempcount = 0
521) do i=1,prevcount
522) if (dabs(tempstoich(i)) > 1.d-10) then
523) tempcount = tempcount + 1
524) tempnames(tempcount) = tempnames(i)
525) tempstoich(tempcount) = tempstoich(i)
526) endif
527) enddo
528)
529) tempnames(tempcount+1:) = ''
530) tempstoich(tempcount+1:) = 0.d0
531)
532) ! reallocate
533) allocate(mnrl_dbaserxn%spec_name(tempcount))
534) allocate(mnrl_dbaserxn%stoich(tempcount))
535)
536) ! fill arrays in dbaserxn
537) mnrl_dbaserxn%nspec = tempcount
538) do i=1,tempcount
539) mnrl_dbaserxn%spec_name(i) = tempnames(i)
540) mnrl_dbaserxn%stoich(i) = tempstoich(i)
541) enddo
542)
543) mnrl_dbaserxn%logK = mnrl_dbaserxn%logK + scale*sec_dbaserxn%logK
544)
545) end subroutine BasisSubSpeciesInMineralRxn
546)
547) ! ************************************************************************** !
548)
549) function DatabaseCheckLegitimateLogKs(dbaserxn,species_name,temperatures, &
550) option)
551) !
552) ! Checks whether legitimate log Ks exist for
553) ! all database temperatures if running
554) ! non-isothermal
555) !
556) ! Author: Glenn Hammond
557) ! Date: 01/07/13
558) !
559)
560) use Option_module
561)
562) implicit none
563)
564) type(database_rxn_type), pointer :: dbaserxn
565) character(len=MAXWORDLENGTH) :: species_name
566) PetscReal :: temperatures(:)
567) type(option_type) :: option
568)
569) PetscBool :: DatabaseCheckLegitimateLogKs
570)
571) PetscInt :: itemp
572) character(len=MAXSTRINGLENGTH) :: string
573) character(len=MAXWORDLENGTH) :: word
574)
575) DatabaseCheckLegitimateLogKs = PETSC_TRUE
576)
577) if (.not.associated(dbaserxn) .or. option%use_isothermal) return
578)
579) string = ''
580) do itemp = 1, size(dbaserxn%logK)
581) if (dabs(dbaserxn%logK(itemp) - 500.) < 1.d-10) then
582) write(word,'(f5.1)') temperatures(itemp)
583) string = trim(string) // ' ' // word
584) DatabaseCheckLegitimateLogKs = PETSC_FALSE
585) endif
586) enddo
587)
588) if (.not.DatabaseCheckLegitimateLogKs) then
589) option%io_buffer = 'Undefined log Ks for temperatures (' // &
590) trim(adjustl(string)) // ') for species "' // &
591) trim(species_name) // '" in database.'
592) call printWrnMsg(option)
593) endif
594)
595) end function DatabaseCheckLegitimateLogKs
596)
597) ! ************************************************************************** !
598)
599) subroutine DatabaseRxnDestroy(dbaserxn)
600) !
601) ! Deallocates a database reaction
602) !
603) ! Author: Glenn Hammond
604) ! Date: 05/29/08
605) !
606) use Utility_module, only : DeallocateArray
607)
608) implicit none
609)
610) type(database_rxn_type), pointer :: dbaserxn
611)
612) if (.not.associated(dbaserxn)) return
613)
614) if (associated(dbaserxn%spec_name)) deallocate(dbaserxn%spec_name)
615) nullify(dbaserxn%spec_name)
616) call DeallocateArray(dbaserxn%spec_ids)
617) call DeallocateArray(dbaserxn%stoich)
618) call DeallocateArray(dbaserxn%logK)
619)
620) deallocate(dbaserxn)
621) nullify(dbaserxn)
622)
623) end subroutine DatabaseRxnDestroy
624)
625) end module Reaction_Database_Aux_module