reaction_microbial.F90       coverage:  100.00 %func     86.04 %block


     1) module Reaction_Microbial_module
     2) 
     3)   use Reaction_Microbial_Aux_module
     4)   
     5)   use PFLOTRAN_Constants_module
     6) 
     7)   implicit none
     8)   
     9)   private 
    10) 
    11) #include "petsc/finclude/petscsys.h"
    12) 
    13)   public :: MicrobialRead, &
    14)             RMicrobial
    15) 
    16) contains
    17) 
    18) ! ************************************************************************** !
    19) 
    20) subroutine MicrobialRead(microbial,input,option)
    21)   ! 
    22)   ! Reads chemical species
    23)   ! 
    24)   ! Author: Glenn Hammond
    25)   ! Date: 08/16/12
    26)   ! 
    27) 
    28)   use Option_module
    29)   use String_module
    30)   use Input_Aux_module
    31)   use Utility_module
    32)   
    33)   implicit none
    34)   
    35)   type(microbial_type) :: microbial
    36)   type(input_type), pointer :: input
    37)   type(option_type) :: option
    38)   
    39)   character(len=MAXWORDLENGTH) :: word
    40)   type(microbial_rxn_type), pointer :: microbial_rxn, cur_microbial_rxn
    41)   type(monod_type), pointer :: monod, prev_monod
    42)   type(microbial_biomass_type), pointer :: microbial_biomass
    43)   type(inhibition_type), pointer :: inhibition, prev_inhibition
    44)   
    45)   microbial%nrxn = microbial%nrxn + 1
    46)         
    47)   microbial_rxn => MicrobialRxnCreate()
    48)   nullify(prev_monod)
    49)   nullify(prev_inhibition)
    50)   nullify(microbial_biomass)
    51)   do 
    52)     call InputReadPflotranString(input,option)
    53)     if (InputError(input)) exit
    54)     if (InputCheckExit(input,option)) exit
    55) 
    56)     call InputReadWord(input,option,word,PETSC_TRUE)
    57)     call InputErrorMsg(input,option,'keyword','CHEMISTRY,MICROBIAL_REACTION')
    58)     call StringToUpper(word)   
    59) 
    60)     select case(trim(word))
    61)       case('REACTION')
    62)         ! remainder of string should be the reaction equation
    63)         microbial_rxn%reaction = trim(adjustl(input%buf))
    64)         ! set flag for error message
    65)         if (len_trim(microbial_rxn%reaction) < 2) input%ierr = 1
    66)         call InputErrorMsg(input,option,'reaction string', &
    67)                             'CHEMISTRY,MICROBIAL_REACTION,REACTION')     
    68)       case('RATE_CONSTANT')
    69)         call InputReadDouble(input,option,microbial_rxn%rate_constant)  
    70)         call InputErrorMsg(input,option,'rate constant', &
    71)                            'CHEMISTRY,MICROBIAL_REACTION') 
    72)       case('ACTIVATION_ENERGY')
    73)         call InputReadDouble(input,option,microbial_rxn%activation_energy)  
    74)         call InputErrorMsg(input,option,'activation energy', &
    75)                            'CHEMISTRY,MICROBIAL_REACTION') 
    76)       case('MONOD')
    77)         monod => MicrobialMonodCreate()
    78)         do 
    79)           call InputReadPflotranString(input,option)
    80)           if (InputError(input)) exit
    81)           if (InputCheckExit(input,option)) exit
    82)           call InputReadWord(input,option,word,PETSC_TRUE)
    83)           call InputErrorMsg(input,option,'keyword', &
    84)                              'CHEMISTRY,MICROBIAL_REACTION,MONOD')
    85)           call StringToUpper(word)   
    86)           select case(trim(word))
    87)             case('SPECIES_NAME')
    88)               call InputReadWord(input,option,word,PETSC_TRUE)
    89)               call InputErrorMsg(input,option,'species name', &
    90)                            'CHEMISTRY,MICROBIAL_REACTION,MONOD')
    91)               monod%species_name = word
    92)             case('HALF_SATURATION_CONSTANT')
    93)               call InputReadDouble(input,option,monod%half_saturation_constant)
    94)               call InputErrorMsg(input,option,'half saturation constant', &
    95)                                  'CHEMISTRY,MICROBIAL_REACTION,MONOD')
    96)             case('THRESHOLD_CONCENTRATION')
    97)               call InputReadDouble(input,option,monod%threshold_concentration)  
    98)               call InputErrorMsg(input,option,'threshold concdntration', &
    99)                                  'CHEMISTRY,MICROBIAL_REACTION,MONOD')
   100)             case default
   101)               call InputKeywordUnrecognized(word, &
   102)                                             'CHEMISTRY,MICROBIAL_REACTION,MONOD', &
   103)                                             option)
   104)           end select
   105)         enddo
   106)         ! append to list
   107)         if (.not.associated(microbial_rxn%monod)) then
   108)           microbial_rxn%monod => monod
   109)         else
   110)           prev_monod%next => monod
   111)         endif
   112)         prev_monod => monod
   113)         nullify(monod)
   114)       case('INHIBITION')
   115)         inhibition => MicrobialInhibitionCreate()
   116)         do 
   117)           call InputReadPflotranString(input,option)
   118)           if (InputError(input)) exit
   119)           if (InputCheckExit(input,option)) exit
   120)           call InputReadWord(input,option,word,PETSC_TRUE)
   121)           call InputErrorMsg(input,option,'keyword', &
   122)                              'CHEMISTRY,MICROBIAL_REACTION,INHIBITION')
   123)           call StringToUpper(word)   
   124)           select case(trim(word))
   125)             case('SPECIES_NAME')
   126)               call InputReadWord(input,option,word,PETSC_TRUE)
   127)               call InputErrorMsg(input,option,'species name', &
   128)                                  'CHEMISTRY,MICROBIAL_REACTION,INHIBITION')
   129)               inhibition%species_name = word
   130)             case('TYPE')
   131)               call InputReadWord(input,option,word,PETSC_TRUE)
   132)               call InputErrorMsg(input,option,'inhibition type', &
   133)                                  'CHEMISTRY,MICROBIAL_REACTION,INHIBITION')
   134)               call StringToUpper(word)   
   135)               select case(word)
   136)                 case('MONOD')
   137)                   inhibition%itype = INHIBITION_MONOD
   138)                 case('INVERSE_MONOD')
   139)                   inhibition%itype = INHIBITION_INVERSE_MONOD
   140)                 case('THRESHOLD')
   141)                   inhibition%itype = INHIBITION_THRESHOLD
   142)                   call InputReadDouble(input,option, &
   143)                                        inhibition%inhibition_constant2)  
   144)                   call InputErrorMsg(input,option,'scaling factor', &
   145)                                      'CHEMISTRY,MICROBIAL_REACTION,&
   146)                                      &INHIBITION,THRESHOLD_INHIBITION')
   147)                 case default
   148)                   call InputKeywordUnrecognized(word, &
   149)                          'CHEMISTRY,MICROBIAL_REACTION,INHIBITION,TYPE',option)
   150)               end select
   151)             case('INHIBITION_CONSTANT')
   152)               call InputReadDouble(input,option,inhibition%inhibition_constant)  
   153)               call InputErrorMsg(input,option,'inhibition constant', &
   154)                                  'CHEMISTRY,MICROBIAL_REACTION,INHIBITION')
   155)             case default
   156)               call InputKeywordUnrecognized(word, &
   157)                       'CHEMISTRY,MICROBIAL_REACTION,INHIBITION',option)
   158)           end select        
   159)         enddo
   160)         if (len_trim(inhibition%species_name) < 2 .or. &
   161)             inhibition%itype == 0 .or. &
   162)             Uninitialized(inhibition%inhibition_constant)) then
   163)           option%io_buffer = 'A SPECIES_NAME, TYPE, and INHIBITION_CON' // &
   164)             'STANT must be defined for INHIBITION in MICROBIAL_REACTION ' // &
   165)             'with REACTION "' // trim(microbial_rxn%reaction) // '".'
   166)           call printErrMsg(option)
   167)         endif
   168)         ! append to list
   169)         if (.not.associated(microbial_rxn%inhibition)) then
   170)           microbial_rxn%inhibition => inhibition
   171)         else
   172)           prev_inhibition%next => inhibition
   173)         endif
   174)         prev_inhibition => inhibition
   175)         nullify(inhibition)
   176)       case('BIOMASS')
   177)         if (associated(microbial_biomass)) then
   178)           call MicrobialBiomassDestroy(microbial_biomass)
   179)         endif
   180)         microbial_biomass => MicrobialBiomassCreate()
   181)         do 
   182)           call InputReadPflotranString(input,option)
   183)           if (InputError(input)) exit
   184)           if (InputCheckExit(input,option)) exit
   185)           call InputReadWord(input,option,word,PETSC_TRUE)
   186)           call InputErrorMsg(input,option,'keyword', &
   187)                              'CHEMISTRY,MICROBIAL_REACTION,BIOMASS')
   188)           call StringToUpper(word)   
   189)           select case(trim(word))
   190)             case('SPECIES_NAME')
   191)               call InputReadWord(input,option,word,PETSC_TRUE)
   192)               call InputErrorMsg(input,option,'species name', &
   193)                                  'CHEMISTRY,MICROBIAL_REACTION,BIOMASS')
   194)               microbial_biomass%species_name = word
   195)             case('YIELD')
   196)               call InputReadDouble(input,option,microbial_biomass%yield)
   197)               call InputErrorMsg(input,option,'yield', &
   198)                                  'CHEMISTRY,MICROBIAL_REACTION,BIOMASS')
   199)             case default
   200)               call InputKeywordUnrecognized(word, &
   201)                                       'CHEMISTRY,MICROBIAL_REACTION,BIOMASS', &
   202)                                             option)
   203)           end select
   204)         enddo
   205)       case default
   206)         call InputKeywordUnrecognized(word,'CHEMISTRY,MICROBIAL_REACTION', &
   207)                                       option)
   208)     end select
   209)   enddo
   210)   
   211)   ! add linkage to biomass if exists
   212)   microbial_rxn%biomass => microbial_biomass
   213)   
   214)   ! add to list
   215)   if (.not.associated(microbial%microbial_rxn_list)) then
   216)     microbial%microbial_rxn_list => microbial_rxn
   217)     microbial_rxn%id = 1
   218)   else
   219)     cur_microbial_rxn => microbial%microbial_rxn_list
   220)     do
   221)       if (.not.associated(cur_microbial_rxn%next)) then
   222)         cur_microbial_rxn%next => microbial_rxn
   223)         microbial_rxn%id = cur_microbial_rxn%id + 1
   224)         exit
   225)       endif
   226)       cur_microbial_rxn => cur_microbial_rxn%next
   227)     enddo
   228)   endif
   229)   nullify(microbial_rxn)
   230) 
   231) end subroutine MicrobialRead
   232) 
   233) ! ************************************************************************** !
   234) 
   235) subroutine RMicrobial(Res,Jac,compute_derivative,rt_auxvar, &
   236)                       global_auxvar,material_auxvar,reaction,option)
   237)   ! 
   238)   ! Computes the microbial reaction
   239)   ! 
   240)   ! Author: Glenn Hammond
   241)   ! Date: 10/31/12
   242)   ! 
   243) 
   244)   use Option_module, only : option_type, printErrMsg
   245)   use Reactive_Transport_Aux_module, only : reactive_transport_auxvar_type
   246)   use Global_Aux_module, only : global_auxvar_type
   247)   use Material_Aux_class, only : material_auxvar_type
   248)   use Reaction_Aux_module, only : reaction_type
   249)   use Reaction_Immobile_Aux_module, only : immobile_type
   250)   
   251)   implicit none
   252)   
   253)   type(option_type) :: option
   254)   type(reaction_type) :: reaction
   255)   PetscBool :: compute_derivative
   256)   PetscReal :: Res(reaction%ncomp)
   257)   PetscReal :: Jac(reaction%ncomp,reaction%ncomp)
   258)   type(reactive_transport_auxvar_type) :: rt_auxvar
   259)   type(global_auxvar_type) :: global_auxvar
   260)   class(material_auxvar_type) :: material_auxvar
   261)   
   262)   PetscInt, parameter :: iphase = 1
   263)   PetscReal :: por_sat_vol
   264)   PetscInt :: irxn, i, ii, icomp, jcomp, ncomp
   265)   PetscInt :: imonod, iinhibition, ibiomass
   266)   PetscReal :: Im
   267)   PetscReal :: rate_constant
   268)   PetscReal :: activity
   269)   PetscReal :: act_coef
   270)   PetscReal :: monod(10)
   271)   PetscReal :: inhibition(10)
   272)   PetscReal :: biomass_conc, yield
   273)   PetscInt :: immobile_id
   274)   PetscReal :: denominator, dR_dX, dX_dc, dR_dc, dR_dbiomass
   275)   PetscReal :: tempreal
   276)   type(microbial_type), pointer :: microbial
   277)   type(immobile_type), pointer :: immobile
   278)   
   279)   microbial => reaction%microbial
   280)   immobile => reaction%immobile
   281)   
   282)   ! units:
   283)   ! Residual: mol/sec
   284)   ! Jacobian: (mol/sec)*(kg water/mol) = kg water/sec
   285)   
   286)   do irxn = 1, microbial%nrxn
   287)   
   288)     ! units:
   289)     !   without biomass: mol/L-sec
   290)     !   with biomass: mol/L-sec * (m^3 bulk / mol biomass)
   291) 
   292)     ncomp = microbial%specid(0,irxn)
   293)     rate_constant = microbial%rate_constant(irxn)
   294)     Im = rate_constant
   295)     if (associated(microbial%activation_energy)) then
   296)       ! ideal gas constant units: J/mol-K
   297)       Im = Im * exp(microbial%activation_energy(irxn)/IDEAL_GAS_CONSTANT* &
   298)                     (1.d0/298.15d0-1.d0/(global_auxvar%temp+273.15d0)))
   299)     endif
   300)     yield = 0.d0
   301)     biomass_conc = 0.d0
   302) 
   303)     ! monod expressions
   304)     do ii = 1, microbial%monodid(0,irxn)
   305)       imonod = microbial%monodid(ii,irxn)
   306)       icomp = microbial%monod_specid(imonod)
   307)       activity = rt_auxvar%pri_molal(icomp)*rt_auxvar%pri_act_coef(icomp)
   308)       monod(ii) = (activity - microbial%monod_Cth(imonod)) / &
   309)                   (microbial%monod_K(imonod) + activity - &
   310)                    microbial%monod_Cth(imonod))
   311)       Im = Im*monod(ii)
   312)     enddo
   313) 
   314)     ! inhibition expressions
   315)     do ii = 1, microbial%inhibitionid(0,irxn)
   316)       iinhibition = microbial%inhibitionid(ii,irxn)
   317)       icomp = microbial%inhibition_specid(iinhibition)
   318)       activity = rt_auxvar%pri_molal(icomp)*rt_auxvar%pri_act_coef(icomp)
   319)       select case(microbial%inhibition_type(iinhibition))
   320)         case(INHIBITION_MONOD)
   321)           inhibition(ii) = microbial%inhibition_C(iinhibition) / &
   322)                           (microbial%inhibition_C(iinhibition) + activity)
   323)         case(INHIBITION_INVERSE_MONOD)
   324)           inhibition(ii) = activity / &
   325)                           (microbial%inhibition_C(iinhibition) + activity)
   326)         case(INHIBITION_THRESHOLD)
   327)           ! if microbial%inhibition_C2 is negative, inhibition kicks in 
   328)           ! above the concentration
   329)           inhibition(ii) = 0.5d0 + &
   330)                            atan((activity - &
   331)                                  microbial%inhibition_C(iinhibition)) * &
   332)                                 microbial%inhibition_C2(iinhibition)) / PI
   333)       end select        
   334)       Im = Im*inhibition(ii)
   335)     enddo
   336)     
   337)     ! biomass term
   338)     ibiomass = microbial%biomassid(irxn)
   339)     if (ibiomass > 0) then
   340)       immobile_id = reaction%offset_immobile + ibiomass
   341)       biomass_conc = rt_auxvar%immobile(ibiomass)
   342)       yield = microbial%biomass_yield(irxn)
   343)       Im = Im*biomass_conc
   344)     endif
   345)     
   346)     ! por_sat_vol units: m^3 water
   347)     por_sat_vol = material_auxvar%porosity*global_auxvar%sat(iphase)* &
   348)                   material_auxvar%volume
   349) 
   350)     ! Im units (before): mol/L-sec
   351)     Im = Im * 1.d3*por_sat_vol
   352)     ! Im units (after): mol/sec
   353)     
   354)     do i = 1, ncomp
   355)       icomp = microbial%specid(i,irxn)
   356)       Res(icomp) = Res(icomp) - microbial%stoich(i,irxn)*Im
   357)     enddo
   358) 
   359)     if (ibiomass > 0) then
   360)       Res(immobile_id) = Res(immobile_id) - yield*Im
   361)     endif
   362)     
   363)     if (.not. compute_derivative) cycle
   364)     
   365)     ! monod expressions
   366)     do ii = 1, microbial%monodid(0,irxn)
   367)       imonod = microbial%monodid(ii,irxn)
   368)       jcomp = microbial%monod_specid(imonod)
   369)       act_coef = rt_auxvar%pri_act_coef(jcomp)
   370)       activity = rt_auxvar%pri_molal(jcomp)*act_coef
   371)         
   372)       dR_dX = Im / monod(ii)
   373)         
   374)       denominator = microbial%monod_K(imonod) + activity - &
   375)                     microbial%monod_Cth(imonod)
   376)         
   377)       dX_dc = act_coef / denominator - &
   378)               act_coef * (activity - microbial%monod_Cth(imonod)) / &
   379)               (denominator*denominator)
   380)         
   381)       dR_dc = -1.d0*dR_dX*dX_dc
   382)       do i = 1, ncomp
   383)         icomp = microbial%specid(i,irxn)
   384)         ! units = (mol/sec)*(kg water/mol) = kg water/sec
   385)         Jac(icomp,jcomp) = Jac(icomp,jcomp) + &
   386)                             microbial%stoich(i,irxn)*dR_dc
   387)       enddo
   388)       if (ibiomass > 0) then
   389)         Jac(immobile_id,jcomp) = Jac(immobile_id,jcomp) + yield*dR_dc      
   390)       endif      
   391)     enddo
   392) 
   393)     ! inhibition expressions
   394)     do ii = 1, microbial%inhibitionid(0,irxn)
   395)       iinhibition = microbial%inhibitionid(ii,irxn)
   396)       jcomp = microbial%inhibition_specid(iinhibition)
   397)       act_coef = rt_auxvar%pri_act_coef(jcomp)
   398)       activity = rt_auxvar%pri_molal(jcomp)*act_coef
   399) 
   400)       dR_dX = Im / inhibition(ii)
   401)         
   402)       select case(microbial%inhibition_type(iinhibition))
   403)         case(INHIBITION_MONOD)
   404)           denominator = microbial%inhibition_C(iinhibition) + activity
   405)           dX_dc = -1.d0 * act_coef *microbial%inhibition_C(iinhibition) / &
   406)                   (denominator*denominator)
   407)         case(INHIBITION_INVERSE_MONOD)
   408)           denominator = microbial%inhibition_C(iinhibition) + activity
   409)           dX_dc = act_coef / denominator - &
   410)                   act_coef * activity / &
   411)                   (denominator*denominator)
   412)         case(INHIBITION_THRESHOLD)
   413)           ! derivative of atan(X) = 1 / (1 + X^2) dX
   414)           tempreal = (activity - microbial%inhibition_C(iinhibition)) * &
   415)                      microbial%inhibition_C2(iinhibition)
   416)           dX_dc = (microbial%inhibition_C2(iinhibition) * act_coef / &
   417)                    (1.d0 + tempreal*tempreal)) / PI
   418)       end select        
   419)       
   420)       dR_dc = -1.d0*dR_dX*dX_dc
   421)       do i = 1, ncomp
   422)         icomp = microbial%specid(i,irxn)
   423)         ! units = (mol/sec)*(kg water/mol) = kg water/sec
   424)         Jac(icomp,jcomp) = Jac(icomp,jcomp) + &
   425)                             microbial%stoich(i,irxn)*dR_dc
   426)       enddo
   427)       if (ibiomass > 0) then
   428)         Jac(immobile_id,jcomp) = Jac(immobile_id,jcomp) + yield*dR_dc        
   429)       endif
   430)     enddo
   431) 
   432)     ! biomass expression
   433)     if (ibiomass > 0) then
   434)       dR_dbiomass = -1.d0*Im / biomass_conc
   435)       do i = 1, ncomp
   436)         icomp = microbial%specid(i,irxn)
   437)         ! units = (mol/sec)*(m^3/mol) = m^3/sec
   438)         Jac(icomp,immobile_id) = Jac(icomp,immobile_id) + &
   439)                             microbial%stoich(i,irxn)*dR_dbiomass
   440)       enddo
   441)       Jac(immobile_id,immobile_id) = Jac(immobile_id,immobile_id) + &
   442)         yield*dR_dbiomass
   443)     endif
   444)     
   445)   enddo
   446)     
   447) end subroutine RMicrobial
   448) 
   449) end module Reaction_Microbial_module

generated by
Intel(R) C++/Fortran Compiler code-coverage tool
Web-Page Owner: Nobody