reaction_sandbox_pnnl_cyber.F90       coverage:  100.00 %func     92.26 %block


     1) module Reaction_Sandbox_Cyber_class
     2) 
     3)   use Reaction_Sandbox_Base_class
     4)   
     5)   use Global_Aux_module
     6)   use Reactive_Transport_Aux_module
     7)   
     8)   use PFLOTRAN_Constants_module
     9) 
    10)   implicit none
    11)   
    12)   private
    13)   
    14) #include "petsc/finclude/petscsys.h"
    15) 
    16) 
    17)   type, public, &
    18)     extends(reaction_sandbox_base_type) :: reaction_sandbox_cyber_type
    19)     PetscInt :: nh4_id
    20)     PetscInt :: o2_id
    21)     PetscInt :: no3_id
    22)     PetscInt :: no2_id
    23)     PetscInt :: n2_id
    24)     PetscInt :: doc_id
    25)     PetscInt :: biomass_id
    26)     PetscInt :: co2_id
    27)     PetscReal :: f1
    28)     PetscReal :: f2
    29)     PetscReal :: f3
    30)     PetscReal :: f_act   ! fraction of active biomass
    31)     PetscReal :: k_deg   ! biomass degradation rate
    32)     PetscReal :: k1      ! nitrate rate constant
    33)     PetscReal :: k2      ! nitrite rate constant
    34)     PetscReal :: k3      ! oxygen rate constant
    35)     PetscReal :: Kd1
    36)     PetscReal :: Ka1
    37)     PetscReal :: Kd2
    38)     PetscReal :: Ka2
    39)     PetscReal :: Kd3
    40)     PetscReal :: Ka3
    41)     PetscReal :: stoich_1_doc
    42)     PetscReal :: stoich_1_nh4
    43)     PetscReal :: stoich_1_no3
    44)     PetscReal :: stoich_1_no2
    45)     PetscReal :: stoich_1_co2
    46)     PetscReal :: stoich_1_biomass
    47)     PetscReal :: stoich_2_doc
    48)     PetscReal :: stoich_2_nh4
    49)     PetscReal :: stoich_2_no2
    50)     PetscReal :: stoich_2_n2
    51)     PetscReal :: stoich_2_co2
    52)     PetscReal :: stoich_2_biomass
    53)     PetscReal :: stoich_3_doc
    54)     PetscReal :: stoich_3_nh4
    55)     PetscReal :: stoich_3_o2
    56)     PetscReal :: stoich_3_co2
    57)     PetscReal :: stoich_3_biomass
    58)     PetscInt :: nrxn
    59)     PetscInt, pointer :: nrow(:)
    60)     PetscInt, pointer :: ncol(:)
    61)     PetscInt, pointer :: irow(:,:)
    62)     PetscInt, pointer :: icol(:,:)
    63)     PetscReal, pointer :: stoich_row(:,:)
    64)   contains
    65)     procedure, public :: ReadInput => CyberRead
    66)     procedure, public :: Setup => CyberSetup
    67)     procedure, public :: Evaluate => CyberReact
    68)     procedure, public :: Destroy => CyberDestroy
    69)   end type reaction_sandbox_cyber_type
    70)   
    71)   public :: CyberCreate
    72) 
    73) contains
    74) 
    75) ! ************************************************************************** !
    76) 
    77) function CyberCreate()
    78)   ! 
    79)   ! Allocates PNNL N reaction object.
    80)   ! 
    81)   ! Author: Glenn Hammond
    82)   ! Date: 10/01/15
    83)   ! 
    84) 
    85)   implicit none
    86)   
    87)   class(reaction_sandbox_cyber_type), pointer :: CyberCreate
    88) 
    89)   allocate(CyberCreate)
    90)   CyberCreate%o2_id = UNINITIALIZED_INTEGER
    91)   CyberCreate%nh4_id = UNINITIALIZED_INTEGER
    92)   CyberCreate%no3_id = UNINITIALIZED_INTEGER
    93)   CyberCreate%no2_id = UNINITIALIZED_INTEGER
    94)   CyberCreate%n2_id = UNINITIALIZED_INTEGER
    95)   CyberCreate%doc_id = UNINITIALIZED_INTEGER
    96)   CyberCreate%biomass_id = UNINITIALIZED_INTEGER
    97)   CyberCreate%co2_id = UNINITIALIZED_INTEGER
    98)   CyberCreate%f1 = UNINITIALIZED_DOUBLE  
    99)   CyberCreate%f2 = UNINITIALIZED_DOUBLE  
   100)   CyberCreate%f3 = UNINITIALIZED_DOUBLE  
   101)   CyberCreate%f_act = UNINITIALIZED_DOUBLE  
   102)   CyberCreate%k_deg = UNINITIALIZED_DOUBLE  
   103)   CyberCreate%k1 = UNINITIALIZED_DOUBLE  
   104)   CyberCreate%k2 = UNINITIALIZED_DOUBLE  
   105)   CyberCreate%k3 = UNINITIALIZED_DOUBLE  
   106)   CyberCreate%Kd1 = UNINITIALIZED_DOUBLE  
   107)   CyberCreate%Ka1 = UNINITIALIZED_DOUBLE  
   108)   CyberCreate%Kd2 = UNINITIALIZED_DOUBLE  
   109)   CyberCreate%Ka2 = UNINITIALIZED_DOUBLE  
   110)   CyberCreate%Kd3 = UNINITIALIZED_DOUBLE  
   111)   CyberCreate%Ka3 = UNINITIALIZED_DOUBLE  
   112)   CyberCreate%stoich_1_doc = UNINITIALIZED_DOUBLE  
   113)   CyberCreate%stoich_1_nh4 = UNINITIALIZED_DOUBLE  
   114)   CyberCreate%stoich_1_no3 = UNINITIALIZED_DOUBLE  
   115)   CyberCreate%stoich_1_no2 = UNINITIALIZED_DOUBLE  
   116)   CyberCreate%stoich_1_co2 = UNINITIALIZED_DOUBLE  
   117)   CyberCreate%stoich_1_biomass = UNINITIALIZED_DOUBLE  
   118)   CyberCreate%stoich_2_doc = UNINITIALIZED_DOUBLE  
   119)   CyberCreate%stoich_2_nh4 = UNINITIALIZED_DOUBLE  
   120)   CyberCreate%stoich_2_no2 = UNINITIALIZED_DOUBLE  
   121)   CyberCreate%stoich_2_n2 = UNINITIALIZED_DOUBLE  
   122)   CyberCreate%stoich_2_co2 = UNINITIALIZED_DOUBLE  
   123)   CyberCreate%stoich_2_biomass = UNINITIALIZED_DOUBLE  
   124)   CyberCreate%stoich_3_doc = UNINITIALIZED_DOUBLE  
   125)   CyberCreate%stoich_3_nh4 = UNINITIALIZED_DOUBLE  
   126)   CyberCreate%stoich_3_o2 = UNINITIALIZED_DOUBLE  
   127)   CyberCreate%stoich_3_co2 = UNINITIALIZED_DOUBLE  
   128)   CyberCreate%stoich_3_biomass = UNINITIALIZED_DOUBLE  
   129)   CyberCreate%nrxn = UNINITIALIZED_INTEGER
   130)   nullify(CyberCreate%nrow)
   131)   nullify(CyberCreate%ncol)
   132)   nullify(CyberCreate%irow)
   133)   nullify(CyberCreate%icol)
   134)   nullify(CyberCreate%stoich_row)
   135) 
   136)   nullify(CyberCreate%next)  
   137)       
   138) end function CyberCreate
   139) 
   140) ! ************************************************************************** !
   141) 
   142) subroutine CyberRead(this,input,option)
   143)   ! 
   144)   ! Reads input deck for PNNL N reaction parameters (if any)
   145)   ! 
   146)   ! Author: Glenn Hammond
   147)   ! Date: 10/01/15
   148)   ! 
   149) 
   150)   use Option_module
   151)   use String_module
   152)   use Input_Aux_module
   153)   
   154)   implicit none
   155)   
   156)   class(reaction_sandbox_cyber_type) :: this
   157)   type(input_type), pointer :: input
   158)   type(option_type) :: option
   159) 
   160)   PetscInt :: i
   161)   character(len=MAXWORDLENGTH) :: word, internal_units, units
   162)   character(len=MAXSTRINGLENGTH) :: error_string 
   163) 
   164)   error_string = 'CHEMISTRY,REACTION_SANDBOX,CYBERNETIC'
   165)   
   166)   do 
   167)     call InputReadPflotranString(input,option)
   168)     if (InputError(input)) exit
   169)     if (InputCheckExit(input,option)) exit
   170) 
   171)     call InputReadWord(input,option,word,PETSC_TRUE)
   172)     call InputErrorMsg(input,option,'keyword',error_string)
   173)     call StringToUpper(word)   
   174) 
   175)     select case(trim(word))
   176)       case('F1')
   177)         call InputReadDouble(input,option,this%f1)  
   178)         call InputErrorMsg(input,option,'f1',error_string)
   179)       case('F2')
   180)         call InputReadDouble(input,option,this%f2)  
   181)         call InputErrorMsg(input,option,'f2',error_string)
   182)       case('F3')
   183)         call InputReadDouble(input,option,this%f3)  
   184)         call InputErrorMsg(input,option,'f3',error_string)
   185)       case('K1','K_NO3-')
   186)         call InputReadDouble(input,option,this%k1)  
   187)         call InputErrorMsg(input,option,'k1',error_string)
   188)         call InputReadAndConvertUnits(input,this%k1,'1/sec', &
   189)                                       error_string//',k1',option)
   190)       case('K2','K_NO2-')
   191)         call InputReadDouble(input,option,this%k2)  
   192)         call InputErrorMsg(input,option,'k2',error_string)
   193)         call InputReadAndConvertUnits(input,this%k2,'1/sec', &
   194)                                       error_string//',k2',option)
   195)       case('K3','K_O2(aq)')
   196)         call InputReadDouble(input,option,this%k3)  
   197)         call InputErrorMsg(input,option,'k3',error_string)
   198)         call InputReadAndConvertUnits(input,this%k3,'1/sec', &
   199)                                       error_string//',k3',option)
   200)       case('KA1','KA_NO3-')
   201)         call InputReadDouble(input,option,this%Ka1)  
   202)         call InputErrorMsg(input,option,'Ka1',error_string)
   203)         call InputReadAndConvertUnits(input,this%Ka1,'M', &
   204)                                       error_string//',Ka1',option)
   205)       case('KA2','KA_NO2-')
   206)         call InputReadDouble(input,option,this%Ka2)  
   207)         call InputErrorMsg(input,option,'Ka2',error_string)
   208)         call InputReadAndConvertUnits(input,this%Ka2,'M', &
   209)                                       error_string//',Ka2',option)
   210)       case('KA3','KA_O2(aq)')
   211)         call InputReadDouble(input,option,this%Ka3)  
   212)         call InputErrorMsg(input,option,'Ka3',error_string)
   213)         call InputReadAndConvertUnits(input,this%Ka3,'M', &
   214)                                       error_string//',Ka3',option)
   215)       case('KD1','KD_NO3-')
   216)         call InputReadDouble(input,option,this%Kd1)  
   217)         call InputErrorMsg(input,option,'Kd1',error_string)
   218)         call InputReadAndConvertUnits(input,this%Kd1,'M', &
   219)                                       error_string//',Kd1',option)
   220)       case('KD2','KD_NO2-')
   221)         call InputReadDouble(input,option,this%Kd2)  
   222)         call InputErrorMsg(input,option,'Kd2',error_string)
   223)         call InputReadAndConvertUnits(input,this%Kd2,'M', &
   224)                                       error_string//',Kd2',option)
   225)       case('KD3','KD_O2(aq)')
   226)         call InputReadDouble(input,option,this%Kd3)  
   227)         call InputErrorMsg(input,option,'Kd3',error_string)
   228)         call InputReadAndConvertUnits(input,this%Kd3,'M', &
   229)                                       error_string//',Kd3',option)
   230)       case('KDEG')
   231)         call InputReadDouble(input,option,this%k_deg)  
   232)         call InputErrorMsg(input,option,'kdeg',error_string)
   233)         call InputReadAndConvertUnits(input,this%k_deg,'1/sec', &
   234)                                       error_string//',kdeg',option)
   235)       case('F_ACT')
   236)         call InputReadDouble(input,option,this%f_act)  
   237)         call InputErrorMsg(input,option,'f_act',error_string)
   238)       case default
   239)         call InputKeywordUnrecognized(word,error_string,option)
   240)     end select
   241)   enddo
   242)   
   243) end subroutine CyberRead
   244) 
   245) ! ************************************************************************** !
   246) 
   247) subroutine CyberSetup(this,reaction,option)
   248)   ! 
   249)   ! Sets up the PNNL N reaction either with parameters either
   250)   ! read from the input deck or hardwired.
   251)   ! 
   252)   ! Author: Glenn Hammond
   253)   ! Date: 10/01/15
   254)   ! 
   255) 
   256)   use Reaction_Aux_module, only : reaction_type, GetPrimarySpeciesIDFromName
   257)   use Reaction_Immobile_Aux_module, only : GetImmobileSpeciesIDFromName
   258)   use Option_module
   259) 
   260)   implicit none
   261)   
   262)   class(reaction_sandbox_cyber_type) :: this
   263)   type(reaction_type) :: reaction
   264)   type(option_type) :: option
   265)   
   266)   character(len=MAXWORDLENGTH) :: word
   267)   PetscInt :: irxn
   268)   
   269)   PetscReal, parameter :: per_day_to_per_sec = 1.d0 / 24.d0 / 3600.d0
   270) 
   271)   word = 'NH4+'
   272)   this%nh4_id = &
   273)     GetPrimarySpeciesIDFromName(word,reaction,option)
   274)   word = 'O2(aq)'
   275)   this%o2_id = &
   276)     GetPrimarySpeciesIDFromName(word,reaction,option)
   277)   word = 'NO3-'
   278)   this%no3_id = &
   279)     GetPrimarySpeciesIDFromName(word,reaction,option)
   280)   word = 'NO2-'
   281)   this%no2_id = &
   282)     GetPrimarySpeciesIDFromName(word,reaction,option)
   283)   word = 'N2(aq)'
   284)   this%n2_id = &
   285)     GetPrimarySpeciesIDFromName(word,reaction,option)
   286)   word = 'CH2O(aq)'
   287)   this%doc_id = &
   288)     GetPrimarySpeciesIDFromName(word,reaction,option)
   289)   word = 'C5H7O2N(aq)'
   290)   this%biomass_id = &
   291)     GetPrimarySpeciesIDFromName(word,reaction,option)
   292) !    GetImmobileSpeciesIDFromName(word,reaction%immobile,option) + reaction%offset_immobile
   293)   word = 'CO2(aq)'
   294)   this%co2_id = &
   295)     GetPrimarySpeciesIDFromName(word,reaction,option)
   296)   
   297)   ! constants based on Hyun's writeup on 6/21/16 entitled "Mini-cybernetic 
   298)   ! model of batch denitrification process"
   299)   if (Uninitialized(this%f1)) this%f1 = 0.497d0
   300)   if (Uninitialized(this%f2)) this%f2 = 0.999d0
   301)   if (Uninitialized(this%f3)) this%f3 = 0.066d0
   302)   if (Uninitialized(this%k1)) this%k1 = 25.04d0 * per_day_to_per_sec
   303)   if (Uninitialized(this%k2)) this%k2 = 17.82d0 * per_day_to_per_sec
   304)   if (Uninitialized(this%k3)) this%k3 = 75.12d0 * per_day_to_per_sec
   305)   if (Uninitialized(this%Kd1)) this%Kd1 = 0.25d-3
   306)   if (Uninitialized(this%Kd2)) this%Kd2 = 0.25d-3
   307)   if (Uninitialized(this%Kd3)) this%Kd3 = 0.25d-3
   308)   if (Uninitialized(this%Ka1)) this%Ka1 = 0.001d-3
   309)   if (Uninitialized(this%Ka2)) this%Ka2 = 0.004d-3
   310)   if (Uninitialized(this%Ka3)) this%Ka3 = 0.001d-3
   311)   if (Uninitialized(this%f_act)) this%f_act = 0.126d0
   312)   if (Uninitialized(this%k_deg)) this%k_deg = 0.532d0 * per_day_to_per_sec
   313) 
   314) ! uncomment these to zero out reactions
   315) !  this%k1 = 0.d0
   316) !  this%k2 = 0.d0
   317) !  this%k3 = 0.d0
   318) !  this%Kd1 = 0.d0
   319) !  this%Ka1 = 0.d0
   320) !  this%f_act = 1.d40
   321) !  this%k_deg = 0.d0
   322) 
   323)   this%stoich_1_doc = -1.d0
   324)   this%stoich_1_nh4 = -0.2d0*(1.d0-this%f1)
   325)   this%stoich_1_no3 = -2.d0*this%f1
   326)   this%stoich_1_no2 = 2.d0*this%f1
   327)   this%stoich_1_co2 = this%f1
   328)   this%stoich_1_biomass = 0.2d0*(1.d0-this%f1)
   329) 
   330)   this%stoich_2_doc = -1.d0
   331)   this%stoich_2_nh4 = -0.2d0*(1.d0-this%f2)
   332)   this%stoich_2_no2 = -4.d0/3.d0*this%f2
   333)   this%stoich_2_n2 = 2.d0/3.d0*this%f2
   334)   this%stoich_2_co2 = this%f2
   335)   this%stoich_2_biomass = 0.2d0*(1.d0-this%f2)
   336) 
   337)   this%stoich_3_doc = -1.d0
   338)   this%stoich_3_nh4 = -0.2d0*(1.d0-this%f3)
   339)   this%stoich_3_o2 = -1.d0*this%f3
   340)   this%stoich_3_co2 = this%f3
   341)   this%stoich_3_biomass = 0.2d0*(1.d0-this%f3)
   342) 
   343)   if (this%k3 > 1.d-40) then
   344)   this%nrxn = 3
   345)   else
   346)     this%nrxn = 2
   347)   endif
   348)   allocate(this%nrow(this%nrxn))
   349)   this%nrow = UNINITIALIZED_INTEGER
   350)   allocate(this%ncol(this%nrxn))
   351)   this%ncol = UNINITIALIZED_INTEGER
   352)   allocate(this%irow(6,this%nrxn))
   353)   this%irow = UNINITIALIZED_INTEGER
   354)   allocate(this%icol(4,this%nrxn))
   355)   this%icol = UNINITIALIZED_INTEGER
   356)   allocate(this%stoich_row(6,this%nrxn))
   357)   this%stoich_row = UNINITIALIZED_DOUBLE
   358)   
   359)   ! NO3- -> NO2-
   360)   irxn = 1
   361)   this%nrow(irxn) = 6
   362)   this%irow(1,irxn) = this%doc_id
   363)   this%irow(2,irxn) = this%nh4_id
   364)   this%irow(3,irxn) = this%no3_id
   365)   this%irow(4,irxn) = this%no2_id
   366)   this%irow(5,irxn) = this%co2_id
   367)   this%irow(6,irxn) = this%biomass_id
   368)   this%stoich_row(1,irxn) = this%stoich_1_doc
   369)   this%stoich_row(2,irxn) = this%stoich_1_nh4
   370)   this%stoich_row(3,irxn) = this%stoich_1_no3
   371)   this%stoich_row(4,irxn) = this%stoich_1_no2
   372)   this%stoich_row(5,irxn) = this%stoich_1_co2
   373)   this%stoich_row(6,irxn) = this%stoich_1_biomass
   374)   this%ncol(irxn) = 4
   375)   this%icol(1,irxn) = this%doc_id
   376)   this%icol(2,irxn) = this%no3_id
   377)   this%icol(3,irxn) = this%no2_id
   378)   this%icol(4,irxn) = this%o2_id
   379)   ! NO2- -> N2
   380)   irxn = 2
   381)   this%nrow(irxn) = 6
   382)   this%irow(1,irxn) = this%doc_id
   383)   this%irow(2,irxn) = this%nh4_id
   384)   this%irow(3,irxn) = this%no2_id
   385)   this%irow(4,irxn) = this%n2_id
   386)   this%irow(5,irxn) = this%co2_id
   387)   this%irow(6,irxn) = this%biomass_id
   388)   this%stoich_row(1,irxn) = this%stoich_2_doc
   389)   this%stoich_row(2,irxn) = this%stoich_2_nh4
   390)   this%stoich_row(3,irxn) = this%stoich_2_no2
   391)   this%stoich_row(4,irxn) = this%stoich_2_n2
   392)   this%stoich_row(5,irxn) = this%stoich_2_co2
   393)   this%stoich_row(6,irxn) = this%stoich_2_biomass
   394)   this%ncol(irxn) = 4
   395)   this%icol(1,irxn) = this%doc_id
   396)   this%icol(2,irxn) = this%no3_id
   397)   this%icol(3,irxn) = this%no2_id
   398)   this%icol(4,irxn) = this%o2_id
   399)   if (this%nrxn > 2) then
   400)   ! O2 -> H2O + CO2
   401)   irxn = 3 
   402)   this%nrow(irxn) = 5
   403)   this%irow(1,irxn) = this%doc_id
   404)   this%irow(2,irxn) = this%nh4_id
   405)   this%irow(3,irxn) = this%o2_id
   406)   this%irow(4,irxn) = this%co2_id
   407)   this%irow(5,irxn) = this%biomass_id
   408)   this%stoich_row(1,irxn) = this%stoich_3_doc
   409)   this%stoich_row(2,irxn) = this%stoich_3_nh4
   410)   this%stoich_row(3,irxn) = this%stoich_3_o2
   411)   this%stoich_row(4,irxn) = this%stoich_3_co2
   412)   this%stoich_row(5,irxn) = this%stoich_3_biomass
   413)   this%ncol(irxn) = 4
   414)   this%icol(1,irxn) = this%doc_id
   415)   this%icol(2,irxn) = this%no3_id
   416)   this%icol(3,irxn) = this%no2_id
   417)   this%icol(4,irxn) = this%o2_id
   418)   endif
   419)   
   420) end subroutine CyberSetup
   421) 
   422) ! ************************************************************************** !
   423) 
   424) subroutine CyberReact(this,Residual,Jacobian,compute_derivative, &
   425)                          rt_auxvar,global_auxvar,material_auxvar,reaction, &
   426)                          option)
   427)   ! 
   428)   ! Evaluates reaction storing residual and/or Jacobian
   429)   ! 
   430)   ! Author: Glenn Hammond
   431)   ! Date: 10/01/15
   432)   ! 
   433) 
   434)   use Option_module
   435)   use Reaction_Aux_module
   436)   use Material_Aux_class
   437)   
   438)   implicit none
   439)   
   440)   class(reaction_sandbox_cyber_type) :: this  
   441)   type(option_type) :: option
   442)   type(reaction_type) :: reaction
   443)   PetscBool :: compute_derivative
   444)   ! the following arrays must be declared after reaction
   445)   PetscReal :: Residual(reaction%ncomp)
   446)   PetscReal :: Jacobian(reaction%ncomp,reaction%ncomp)
   447)   type(reactive_transport_auxvar_type) :: rt_auxvar
   448)   type(global_auxvar_type) :: global_auxvar
   449)   class(material_auxvar_type) :: material_auxvar
   450) 
   451)   PetscInt, parameter :: iphase = 1
   452)   PetscReal :: L_water
   453)   PetscReal :: kg_water
   454)   
   455)   PetscInt :: i, j, irxn
   456) 
   457)   PetscReal :: Co2, Cno3, Cno2, Cn2, Cdoc, Cco2, Cnh4, X
   458)   PetscReal :: r1docmonod, r1docmonod_denom
   459)   PetscReal :: r2docmonod, r2docmonod_denom
   460)   PetscReal :: r3docmonod, r3docmonod_denom
   461)   PetscReal :: r1no3monod, r1no3monod_denom
   462)   PetscReal :: r2no2monod, r2no2monod_denom
   463)   PetscReal :: r3o2monod, r3o2monod_denom
   464)   PetscReal :: r1kin, r2kin, r3kin
   465)   PetscReal :: sumkin, sumkinsq
   466)   PetscReal :: u1, u2, u3
   467)   PetscReal :: dr1kin_ddoc, dr1kin_dno3
   468)   PetscReal :: dr2kin_ddoc, dr2kin_dno2
   469)   PetscReal :: dr3kin_ddoc, dr3kin_do2
   470)   PetscReal :: du_denom_dr
   471)   PetscReal :: du1_dr1kin, du1_dr2kin, du1_dr3kin
   472)   PetscReal :: du2_dr1kin, du2_dr2kin, du2_dr3kin
   473)   PetscReal :: du3_dr1kin, du3_dr2kin, du3_dr3kin
   474)   PetscReal :: dr1_ddoc, dr1_dno3, dr1_dno2, dr1_do2
   475)   PetscReal :: dr2_ddoc, dr2_dno3, dr2_dno2, dr2_do2
   476)   PetscReal :: dr3_ddoc, dr3_dno3, dr3_dno2, dr3_do2
   477)   PetscReal :: du1_ddoc, du1_dno3, du1_dno2, du1_do2
   478)   PetscReal :: du2_ddoc, du2_dno3, du2_dno2, du2_do2
   479)   PetscReal :: du3_ddoc, du3_dno3, du3_dno2, du3_do2
   480)   PetscReal :: molality_to_molarity
   481) 
   482)   PetscReal :: rate(3), derivative_col(6,3)
   483)   
   484)   L_water = material_auxvar%porosity*global_auxvar%sat(iphase)* &
   485)             material_auxvar%volume*1.d3 ! m^3 -> L
   486)   kg_water = material_auxvar%porosity*global_auxvar%sat(iphase)* &
   487)              global_auxvar%den_kg(iphase)*material_auxvar%volume
   488) 
   489)   molality_to_molarity = global_auxvar%den_kg(iphase)*1.d-3
   490)     
   491)   if (reaction%act_coef_update_frequency /= ACT_COEF_FREQUENCY_OFF) then
   492)     option%io_buffer = 'Activity coefficients not currently supported in &
   493)       &CyberReact().'
   494)     call printErrMsg(option)
   495)   endif
   496)   
   497)   ! concentrations are molarities [M]
   498)   Co2 = rt_auxvar%pri_molal(this%o2_id)* &
   499)         rt_auxvar%pri_act_coef(this%o2_id)*molality_to_molarity
   500)   Cnh4 = rt_auxvar%pri_molal(this%nh4_id)* &
   501)          rt_auxvar%pri_act_coef(this%nh4_id)*molality_to_molarity
   502)   Cno3 = rt_auxvar%pri_molal(this%no3_id)* &
   503)          rt_auxvar%pri_act_coef(this%no3_id)*molality_to_molarity
   504)   Cno2 = rt_auxvar%pri_molal(this%no2_id)* &
   505)          rt_auxvar%pri_act_coef(this%no2_id)*molality_to_molarity
   506)   Cn2 = rt_auxvar%pri_molal(this%n2_id)* &
   507)         rt_auxvar%pri_act_coef(this%n2_id)*molality_to_molarity
   508)   Cdoc = rt_auxvar%pri_molal(this%doc_id)* &
   509)          rt_auxvar%pri_act_coef(this%doc_id)*molality_to_molarity
   510) !  X = rt_auxvar%immobile(this%biomass_id-reaction%offset_immobile)
   511)   X = rt_auxvar%pri_molal(this%biomass_id)* &
   512)       rt_auxvar%pri_act_coef(this%biomass_id)*molality_to_molarity
   513)   
   514)   ! NO3- -> NO2-
   515)   r1docmonod_denom = Cdoc+this%Kd1
   516)   r1docmonod = Cdoc/r1docmonod_denom
   517)   r1no3monod_denom = Cno3+this%Ka1
   518)   r1no3monod = Cno3/r1no3monod_denom
   519)   r1kin = this%k1*r1docmonod*r1no3monod
   520)                  
   521)   ! NO2- -> N2
   522)   r2docmonod_denom = Cdoc+this%Kd2 
   523)   r2docmonod = Cdoc/r2docmonod_denom
   524)   r2no2monod_denom = Cno2+this%Ka2
   525)   r2no2monod = Cno2/r2no2monod_denom
   526)   r2kin = this%k2*r2docmonod*r2no2monod
   527)                  
   528)   ! O2 -> 
   529)   r3docmonod_denom = Cdoc+this%Kd3
   530)   r3docmonod = Cdoc/r3docmonod_denom
   531)   r3o2monod_denom = Co2+this%Ka3
   532)   r3o2monod = Co2/r3o2monod_denom
   533)   r3kin = this%k3*r3docmonod*r3o2monod
   534)                  
   535)   sumkin = r1kin + r2kin + r3kin
   536)   sumkinsq = sumkin * sumkin
   537)                  
   538)   u1 = 0.d0
   539)   if (r1kin > 0.d0) u1 = r1kin/sumkin
   540)   u2 = 0.d0
   541)   if (r2kin > 0.d0) u2 = r2kin/sumkin
   542)   u3 = 0.d0
   543)   if (r3kin > 0.d0) u3 = r3kin/sumkin
   544)   
   545)   rate(1) = u1*r1kin  ! mol/mol biomass/sec
   546)   rate(2) = u2*r2kin
   547)   rate(3) = u3*r3kin
   548)   
   549)   do irxn = 1, this%nrxn
   550)     do i = 1, this%nrow(irxn)
   551)       ! mol/sec
   552)       ! X is in [M]
   553)       Residual(this%irow(i,irxn)) = Residual(this%irow(i,irxn)) - &
   554)         this%stoich_row(i,irxn) * rate(irxn) * X * &
   555)         ! if biomass is aqueous multiply by L_water
   556)         ! if biomass is immobile multiply by material_auxvar%volume
   557)         L_water
   558)     enddo
   559)   enddo
   560)   
   561)   ! decay of biomass
   562)   ! if biomass is aqueous multiply by L_water
   563)   ! if biomass is immobile multiply by material_auxvar%volume
   564)   Residual(this%biomass_id) = Residual(this%biomass_id) + &
   565)                               this%k_deg * X * L_water
   566) 
   567)   ! production of doc by biomass decay
   568)   ! note the addition
   569)   ! mol/sec
   570)   Residual(this%doc_id) = Residual(this%doc_id) - &
   571)                           this%k_deg/this%f_act * X * L_water
   572)                  
   573)   if (compute_derivative) then
   574)   
   575)     dr1kin_ddoc = this%k1 * &
   576)                   (r1docmonod/Cdoc - r1docmonod/r1docmonod_denom) * &
   577)                   r1no3monod * &
   578)                   rt_auxvar%pri_act_coef(this%doc_id)
   579)     dr1kin_dno3 = this%k1 * &
   580)                   r1docmonod * &
   581)                   (r1no3monod/Cno3 - r1no3monod/r1no3monod_denom) * &
   582)                   rt_auxvar%pri_act_coef(this%no3_id)
   583)     dr2kin_ddoc = this%k2 * &
   584)                   (r2docmonod/Cdoc - r2docmonod/r2docmonod_denom) * &
   585)                   r2no2monod * &
   586)                   rt_auxvar%pri_act_coef(this%doc_id)
   587)     dr2kin_dno2 = this%k2 * &
   588)                   r2docmonod * &
   589)                   (r2no2monod/Cno2 - r2no2monod/r2no2monod_denom) * &
   590)                   rt_auxvar%pri_act_coef(this%no2_id)
   591)     dr3kin_ddoc = this%k3 * &
   592)                   (r3docmonod/Cdoc - r3docmonod/r3docmonod_denom) * &
   593)                   r3o2monod * &
   594)                   rt_auxvar%pri_act_coef(this%doc_id)
   595)     dr3kin_do2 = this%k3 * &
   596)                   r3docmonod * &
   597)                   (r3o2monod/Co2 - r3o2monod/r3o2monod_denom) * &
   598)                   rt_auxvar%pri_act_coef(this%o2_id)
   599)                 
   600)     du_denom_dr = -1.d0/sumkinsq
   601)     du1_dr1kin = 1.d0/sumkin + r1kin*du_denom_dr
   602)     du1_dr2kin = r1kin*du_denom_dr
   603)     du1_dr3kin = r1kin*du_denom_dr
   604)     du2_dr1kin = r2kin*du_denom_dr
   605)     du2_dr2kin = 1.d0/sumkin + r2kin*du_denom_dr
   606)     du2_dr3kin = r2kin*du_denom_dr
   607)     du3_dr1kin = r3kin*du_denom_dr
   608)     du3_dr2kin = r3kin*du_denom_dr
   609)     du3_dr3kin = 1.d0/sumkin + r3kin*du_denom_dr
   610) 
   611)     du1_ddoc = du1_dr1kin * dr1kin_ddoc + du1_dr2kin * dr2kin_ddoc + &
   612)                du1_dr3kin * dr3kin_ddoc
   613)     du1_dno3 = du1_dr1kin * dr1kin_dno3
   614)     du1_dno2 = du1_dr2kin * dr2kin_dno2
   615)     du1_do2 = du1_dr3kin * dr3kin_do2
   616)     du2_ddoc = du2_dr1kin * dr1kin_ddoc + du2_dr2kin * dr2kin_ddoc + &
   617)                du2_dr3kin * dr3kin_ddoc
   618)     du2_dno3 = du2_dr1kin * dr1kin_dno3
   619)     du2_dno2 = du2_dr2kin * dr2kin_dno2
   620)     du2_do2 = du2_dr3kin * dr3kin_do2
   621)     du3_ddoc = du3_dr1kin * dr1kin_ddoc + du3_dr2kin * dr2kin_ddoc + &
   622)                du3_dr3kin * dr3kin_ddoc
   623)     du3_dno3 = du3_dr1kin * dr1kin_dno3
   624)     du3_dno2 = du3_dr2kin * dr2kin_dno2
   625)     du3_do2 = du3_dr3kin * dr3kin_do2
   626) 
   627)     dr1_ddoc = dr1kin_ddoc * u1 + du1_ddoc * r1kin
   628)     dr1_dno3 = dr1kin_dno3 * u1 + du1_dno3 * r1kin
   629)     dr1_dno2 = du1_dno2 * r1kin
   630)     dr1_do2 = du1_do2 * r1kin
   631)     dr2_ddoc = dr2kin_ddoc * u2 + du2_ddoc * r2kin
   632)     dr2_dno3 = du2_dno3 * r2kin
   633)     dr2_dno2 = dr2kin_dno2 * u2 + du2_dno2 * r2kin
   634)     dr2_do2 = du2_do2 * r2kin
   635)     dr3_ddoc = dr3kin_ddoc * u3 + du3_ddoc * r3kin
   636)     dr3_dno3 = du3_dno3 * r3kin
   637)     dr3_dno2 = du3_dno2 * r3kin
   638)     dr3_do2 = dr3kin_do2 * u3 + du3_do2 * r3kin
   639)   
   640)     irxn = 1
   641)     derivative_col(1,irxn) = dr1_ddoc
   642)     derivative_col(2,irxn) = dr1_dno3
   643)     derivative_col(3,irxn) = dr1_dno2
   644)     derivative_col(4,irxn) = dr1_do2  
   645)     irxn = 2
   646)     derivative_col(1,irxn) = dr2_ddoc
   647)     derivative_col(2,irxn) = dr2_dno3
   648)     derivative_col(3,irxn) = dr2_dno2
   649)     derivative_col(4,irxn) = dr2_do2  
   650)     irxn = 3
   651)     derivative_col(1,irxn) = dr3_ddoc
   652)     derivative_col(2,irxn) = dr3_dno3
   653)     derivative_col(3,irxn) = dr3_dno2
   654)     derivative_col(4,irxn) = dr3_do2      
   655)     
   656)     ! fill the Jacobian
   657)     ! units = kg water/sec. Multiply by kg_water
   658)     do irxn = 1, this%nrxn
   659)       do j = 1, this%ncol(irxn)
   660)         do i = 1, this%nrow(irxn)
   661)           Jacobian(this%irow(i,irxn),this%icol(j,irxn)) = &
   662)             Jacobian(this%irow(i,irxn),this%icol(j,irxn)) - &
   663)             this%stoich_row(i,irxn) * derivative_col(j,irxn) * X * kg_water
   664)         enddo
   665)       enddo
   666)       ! if biomass is aqueous, units = kg water/sec. Multiply by kg_water
   667)       ! if biomass is immobile, units = m^3 bulk/sec. Multiply by 
   668)       !   material_auxvar%volume
   669)       do i = 1, this%nrow(irxn)
   670)         Jacobian(this%irow(i,irxn),this%biomass_id) = &
   671)           Jacobian(this%irow(i,irxn),this%biomass_id) - &
   672)            this%stoich_row(i,irxn) * rate(irxn) * kg_water
   673)       enddo
   674)     enddo 
   675) 
   676)     ! decay of biomass
   677)     ! if biomass is aqueous, units = kg water/sec. Multiply by kg_water
   678)     ! if biomass is immobile, units = m^3 bulk/sec. Multiply by 
   679)     !   material_auxvar%volume
   680)     Jacobian(this%biomass_id,this%biomass_id) = &
   681)       Jacobian(this%biomass_id,this%biomass_id) + &
   682)       this%k_deg * kg_water
   683) 
   684)     ! production of doc by biomass decay
   685)     ! units = kg water/sec. Multiply by kg_water
   686)     Jacobian(this%doc_id,this%biomass_id) = &
   687)       Jacobian(this%doc_id,this%biomass_id) - &
   688)       this%k_deg/this%f_act * kg_water
   689)     
   690)   endif
   691)   
   692) end subroutine CyberReact
   693) 
   694) ! ************************************************************************** !
   695) 
   696) subroutine CyberDestroy(this)
   697)   ! 
   698)   ! Destroys allocatable or pointer objects created in this
   699)   ! module
   700)   ! 
   701)   ! Author: Glenn Hammond
   702)   ! Date: 10/01/15
   703)   ! 
   704)   use Utility_module
   705)   
   706)   implicit none
   707)   
   708)   class(reaction_sandbox_cyber_type) :: this  
   709) 
   710)   call DeallocateArray(this%nrow)
   711)   call DeallocateArray(this%irow)
   712)   call DeallocateArray(this%icol)
   713)   call DeallocateArray(this%stoich_row)
   714) 
   715) end subroutine CyberDestroy
   716) 
   717) end module Reaction_Sandbox_Cyber_class

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