reaction_sandbox_clm_cn.F90       coverage:  100.00 %func     82.69 %block


     1) module Reaction_Sandbox_CLM_CN_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)                           ! 14.00674d0 / 12.011d0
    17)   PetscReal, parameter :: CN_ratio_mass_to_mol = 1.16616d0 
    18)   PetscInt, parameter :: CARBON_INDEX = 1
    19)   PetscInt, parameter :: NITROGEN_INDEX = 2
    20)   PetscInt, parameter :: SOM_INDEX = 1
    21) 
    22)   type, public, &
    23)     extends(reaction_sandbox_base_type) :: reaction_sandbox_clm_cn_type
    24)     PetscInt :: nrxn
    25)     PetscInt :: npool
    26)     PetscReal, pointer :: CN_ratio(:)
    27)     PetscReal, pointer :: rate_constant(:)
    28)     PetscReal, pointer :: respiration_fraction(:)
    29)     PetscReal, pointer :: inhibition_constant(:)
    30)     PetscInt, pointer :: upstream_pool_id(:)
    31)     PetscInt, pointer :: downstream_pool_id(:)
    32)     PetscInt, pointer :: pool_id_to_species_id(:,:)
    33)     PetscInt :: C_species_id
    34)     PetscInt :: N_species_id
    35)     type(pool_type), pointer :: pools
    36)     type(clm_cn_reaction_type), pointer :: reactions
    37)   contains
    38)     procedure, public :: ReadInput => CLM_CN_Read
    39)     procedure, public :: Setup => CLM_CN_Setup
    40)     procedure, public :: Evaluate => CLM_CN_React
    41)     procedure, public :: Destroy => CLM_CN_Destroy
    42)   end type reaction_sandbox_clm_cn_type
    43)   
    44)   type :: pool_type
    45)     character(len=MAXWORDLENGTH) :: name
    46)     PetscReal :: CN_ratio
    47)     type(pool_type), pointer :: next
    48)   end type pool_type
    49)   
    50)   type :: clm_cn_reaction_type
    51)     character(len=MAXWORDLENGTH) :: upstream_pool_name
    52)     character(len=MAXWORDLENGTH) :: downstream_pool_name
    53)     PetscReal :: rate_constant
    54)     PetscReal :: respiration_fraction
    55)     PetscReal :: inhibition_constant
    56)     type(clm_cn_reaction_type), pointer :: next
    57)   end type clm_cn_reaction_type
    58)   
    59)   public :: CLM_CN_Create
    60) 
    61) contains
    62) 
    63) ! ************************************************************************** !
    64) 
    65) function CLM_CN_Create()
    66)   ! 
    67)   ! Allocates CLM-CN reaction object.
    68)   ! 
    69)   ! Author: Glenn Hammond
    70)   ! Date: 02/04/13
    71)   ! 
    72) 
    73)   implicit none
    74)   
    75)   type(reaction_sandbox_clm_cn_type), pointer :: CLM_CN_Create
    76)   
    77)   allocate(CLM_CN_Create)
    78)   CLM_CN_Create%nrxn = 0
    79)   CLM_CN_Create%npool = 0
    80)   nullify(CLM_CN_Create%CN_ratio)
    81)   nullify(CLM_CN_Create%rate_constant)
    82)   nullify(CLM_CN_Create%respiration_fraction)
    83)   nullify(CLM_CN_Create%inhibition_constant)
    84)   nullify(CLM_CN_Create%pool_id_to_species_id)
    85)   nullify(CLM_CN_Create%upstream_pool_id)
    86)   nullify(CLM_CN_Create%downstream_pool_id)
    87)   CLM_CN_Create%C_species_id = 0
    88)   CLM_CN_Create%N_species_id = 0
    89)   nullify(CLM_CN_Create%next)
    90)   nullify(CLM_CN_Create%pools)
    91)   nullify(CLM_CN_Create%reactions)
    92) 
    93) end function CLM_CN_Create
    94) 
    95) ! ************************************************************************** !
    96) 
    97) subroutine CLM_CN_Read(this,input,option)
    98)   ! 
    99)   ! Reads input deck for reaction sandbox parameters
   100)   ! 
   101)   ! Author: Glenn Hammond
   102)   ! Date: 02/04/13
   103)   ! 
   104) 
   105)   use Option_module
   106)   use String_module
   107)   use Input_Aux_module
   108)   use Utility_module
   109)   use Units_module, only : UnitsConvertToInternal
   110)   
   111)   implicit none
   112)   
   113)   class(reaction_sandbox_clm_cn_type) :: this
   114)   type(input_type), pointer :: input
   115)   type(option_type) :: option
   116)   
   117)   character(len=MAXWORDLENGTH) :: word, internal_units
   118)   
   119)   type(pool_type), pointer :: new_pool, prev_pool
   120)   type(clm_cn_reaction_type), pointer :: new_reaction, prev_reaction
   121)   
   122)   PetscReal :: rate_constant, turnover_time
   123)   
   124)   nullify(new_pool)
   125)   nullify(prev_pool)
   126)   
   127)   nullify(new_reaction)
   128)   nullify(prev_reaction)
   129)   
   130)   do 
   131)     call InputReadPflotranString(input,option)
   132)     if (InputError(input)) exit
   133)     if (InputCheckExit(input,option)) exit
   134) 
   135)     call InputReadWord(input,option,word,PETSC_TRUE)
   136)     call InputErrorMsg(input,option,'keyword', &
   137)                        'CHEMISTRY,REACTION_SANDBOX,CLM-CN')
   138)     call StringToUpper(word)   
   139) 
   140)     select case(trim(word))
   141)       case('POOLS')
   142)         do
   143)           call InputReadPflotranString(input,option)
   144)           if (InputError(input)) exit
   145)           if (InputCheckExit(input,option)) exit   
   146) 
   147)           allocate(new_pool)
   148)           new_pool%name = ''
   149)           new_pool%CN_ratio = 0.d0
   150)           nullify(new_pool%next)
   151) 
   152)           call InputReadWord(input,option,new_pool%name,PETSC_TRUE)
   153)           call InputErrorMsg(input,option,'pool name', &
   154)             'CHEMISTRY,REACTION_SANDBOX,CLM-CN,POOLS')
   155)           call InputReadDouble(input,option,new_pool%CN_ratio)
   156)           if (InputError(input)) then
   157)             new_pool%CN_ratio = UNINITIALIZED_DOUBLE
   158)           else
   159)             ! convert CN ratio from mass C/mass N to mol C/mol N
   160)             new_pool%CN_ratio = new_pool%CN_ratio * CN_ratio_mass_to_mol
   161)           endif
   162)           if (associated(this%pools)) then
   163)             prev_pool%next => new_pool
   164)           else
   165)             this%pools => new_pool
   166)           endif
   167)           prev_pool => new_pool
   168)           nullify(new_pool)
   169)         enddo
   170)       case('REACTION')
   171)       
   172)         allocate(new_reaction)
   173)         new_reaction%upstream_pool_name = ''
   174)         new_reaction%downstream_pool_name = ''
   175)         new_reaction%rate_constant = UNINITIALIZED_DOUBLE
   176)         new_reaction%respiration_fraction = UNINITIALIZED_DOUBLE
   177)         new_reaction%inhibition_constant = 0.d0
   178)         nullify(new_reaction%next)
   179)         
   180)         ! need to set these temporarily in order to check that they
   181)         ! are not both set.
   182)         turnover_time = 0.d0
   183)         rate_constant = 0.d0
   184)         
   185)         do 
   186)           call InputReadPflotranString(input,option)
   187)           if (InputError(input)) exit
   188)           if (InputCheckExit(input,option)) exit
   189) 
   190)           call InputReadWord(input,option,word,PETSC_TRUE)
   191)           call InputErrorMsg(input,option,'keyword', &
   192)                              'CHEMISTRY,REACTION_SANDBOX,CLM-CN,REACTION')
   193)           call StringToUpper(word)   
   194) 
   195)           select case(trim(word))
   196)             case('UPSTREAM_POOL')
   197)               call InputReadWord(input,option, &
   198)                                  new_reaction%upstream_pool_name,PETSC_TRUE)
   199)               call InputErrorMsg(input,option,'upstream pool name', &
   200)                      'CHEMISTRY,REACTION_SANDBOX,CLM-CN,REACTION')
   201)             case('DOWNSTREAM_POOL')
   202)               call InputReadWord(input,option, &
   203)                                  new_reaction%downstream_pool_name,PETSC_TRUE)
   204)               call InputErrorMsg(input,option,'downstream pool name', &
   205)                      'CHEMISTRY,REACTION_SANDBOX,CLM-CN,REACTION')
   206)             case('RATE_CONSTANT')
   207)               call InputReadDouble(input,option,rate_constant)
   208)               call InputErrorMsg(input,option,'rate constant', &
   209)                      'CHEMISTRY,REACTION_SANDBOX,CLM-CN,REACTION')
   210)               call InputReadWord(input,option,word,PETSC_TRUE)
   211)               internal_units = 'unitless/sec'
   212)               if (InputError(input)) then
   213)                 input%err_buf = 'CLM-CN RATE CONSTANT UNITS'
   214)                 call InputDefaultMsg(input,option)
   215)               else              
   216)                 rate_constant = rate_constant * &
   217)                   UnitsConvertToInternal(word,internal_units,option)
   218)               endif
   219)             case('TURNOVER_TIME')
   220)               call InputReadDouble(input,option,turnover_time)
   221)               call InputErrorMsg(input,option,'turnover time', &
   222)                      'CHEMISTRY,REACTION_SANDBOX,CLM-CN,REACTION')
   223)               call InputReadWord(input,option,word,PETSC_TRUE)
   224)               internal_units = 'sec'
   225)               if (InputError(input)) then
   226)                 input%err_buf = 'CLM-CN TURNOVER TIME UNITS'
   227)                 call InputDefaultMsg(input,option)
   228)               else              
   229)                 turnover_time = turnover_time * &
   230)                   UnitsConvertToInternal(word,internal_units,option)
   231)               endif
   232)             case('RESPIRATION_FRACTION')
   233)               call InputReadDouble(input,option, &
   234)                                    new_reaction%respiration_fraction)
   235)               call InputErrorMsg(input,option,'respiration fraction', &
   236)                      'CHEMISTRY,REACTION_SANDBOX,CLM-CN,REACTION')
   237)             case('N_INHIBITION')
   238)               call InputReadDouble(input,option, &
   239)                                    new_reaction%inhibition_constant)
   240)               call InputErrorMsg(input,option,'inhibition constant', &
   241)                      'CHEMISTRY,REACTION_SANDBOX,CLM-CN,REACTION')
   242)             case default
   243)               call InputKeywordUnrecognized(word, &
   244)                      'CHEMISTRY,REACTION_SANDBOX,CLM-CN,REACTION',option)
   245)           end select
   246)         enddo
   247)         
   248)         ! check to ensure that one of turnover time or rate constant is set.
   249)         if (turnover_time > 0.d0 .and. rate_constant > 0.d0) then
   250)           option%io_buffer = 'Only TURNOVER_TIME or RATE_CONSTANT may ' // &
   251)             'be included in a CLM-CN reaction definition, but not both. ' // &
   252)             'See reaction with upstream pool "' // &
   253)             trim(new_reaction%upstream_pool_name) // '".'
   254)           call printErrMsg(option)
   255)         else if (turnover_time > 0.d0) then
   256)           new_reaction%rate_constant = 1.d0 / turnover_time
   257)         else
   258)           new_reaction%rate_constant = rate_constant
   259)         endif
   260)         ! ensure that respiration fraction is 0-1.
   261)         if (new_reaction%respiration_fraction < 0.d0 .or. &
   262)                   new_reaction%respiration_fraction > 1.d0) then
   263)           option%io_buffer = 'Respiratory fraction (rf) must be between ' // &
   264)             'zero and one (i.e. 0. <= rf <= 1.) in a CLM-CN reaction ' // &
   265)             'definition. See reaction with upstream pool "' // &
   266)             trim(new_reaction%upstream_pool_name) // '".'
   267)           call printErrMsg(option)
   268)         endif
   269)         ! If no downstream pool exists, ensure that respiration fraction = 1
   270)         if (len_trim(new_reaction%downstream_pool_name) < 1 .and. &
   271)             (1.d0 - new_reaction%respiration_fraction) > 1.d-40) then
   272)           option%io_buffer = 'Respiratory fraction (rf) must be set to ' // &
   273)             '1.0 if no downstream pool is specified in a CLM-CN reaction ' // &
   274)             'definition. See reaction with upstream pool "' // &
   275)             trim(new_reaction%upstream_pool_name) // '".'
   276)           call printErrMsg(option)
   277)         endif
   278)         if (associated(this%reactions)) then
   279)           prev_reaction%next => new_reaction
   280)         else
   281)           this%reactions => new_reaction
   282)         endif
   283)         prev_reaction => new_reaction
   284)         nullify(new_reaction)        
   285)       case default
   286)         call InputKeywordUnrecognized(word, &
   287)                      'CHEMISTRY,REACTION_SANDBOX,CLM-CN',option)
   288)     end select
   289)   enddo
   290)   
   291) end subroutine CLM_CN_Read
   292) 
   293) ! ************************************************************************** !
   294) 
   295) subroutine CLM_CN_Setup(this,reaction,option)
   296)   ! 
   297)   ! Sets up CLM-CN reaction after it has been read from input
   298)   ! 
   299)   ! Author: Glenn Hammond
   300)   ! Date: 02/04/13
   301)   ! 
   302) 
   303)   use Reaction_Aux_module, only : reaction_type
   304)   use Option_module
   305)   
   306)   implicit none
   307)   
   308)   class(reaction_sandbox_clm_cn_type) :: this
   309)   type(reaction_type) :: reaction  
   310)   type(option_type) :: option
   311)   
   312)   call CLM_CN_Map(this,reaction,option)
   313) 
   314) end subroutine CLM_CN_Setup
   315) 
   316) ! ************************************************************************** !
   317) 
   318) subroutine CLM_CN_Map(this,reaction,option)
   319)   ! 
   320)   ! Maps coefficients to primary dependent variables
   321)   ! 
   322)   ! Author: Glenn Hammond
   323)   ! Date: 02/04/13
   324)   ! 
   325) 
   326)   use Reaction_Aux_module, only : reaction_type
   327)   use Option_module
   328)   use String_module
   329)   use Reaction_Immobile_Aux_module
   330)   
   331)   implicit none
   332) 
   333)   class(reaction_sandbox_clm_cn_type) :: this
   334)   type(option_type) :: option
   335)   type(reaction_type) :: reaction
   336)   
   337)   character(len=MAXWORDLENGTH), allocatable :: pool_names(:)
   338)   character(len=MAXWORDLENGTH) :: word
   339)   
   340)   PetscInt :: icount
   341) 
   342)   type(pool_type), pointer :: cur_pool
   343)   type(clm_cn_reaction_type), pointer :: cur_rxn
   344)   
   345)   ! count # pools
   346)   icount = 0
   347)   cur_pool => this%pools
   348)   do
   349)     if (.not.associated(cur_pool)) exit
   350)     icount = icount + 1
   351)     cur_pool => cur_pool%next
   352)   enddo
   353)   this%npool = icount
   354)   
   355)   ! count # reactions
   356)   icount = 0
   357)   cur_rxn => this%reactions
   358)   do
   359)     if (.not.associated(cur_rxn)) exit
   360)     icount = icount + 1
   361)     cur_rxn => cur_rxn%next
   362)   enddo
   363)   this%nrxn = icount
   364)   
   365)   ! allocate and initialize arrays
   366)   allocate(this%CN_ratio(this%npool))
   367)   allocate(this%pool_id_to_species_id(0:2,this%npool))
   368)   allocate(this%upstream_pool_id(this%nrxn))
   369)   allocate(this%downstream_pool_id(this%nrxn))
   370)   allocate(this%rate_constant(this%nrxn))
   371)   allocate(this%respiration_fraction(this%nrxn))
   372)   allocate(this%inhibition_constant(this%nrxn))
   373)   this%CN_ratio = 0.d0
   374)   this%pool_id_to_species_id = 0
   375)   this%upstream_pool_id = 0
   376)   this%downstream_pool_id = 0
   377)   this%rate_constant = 0.d0
   378)   this%respiration_fraction = 0.d0
   379)   this%inhibition_constant = 0.d0
   380)   
   381)   ! temporary array for mapping pools in reactions
   382)   allocate(pool_names(this%npool))
   383)   
   384)   icount = 0
   385)   cur_pool => this%pools
   386)   do
   387)     if (.not.associated(cur_pool)) exit
   388)     icount = icount + 1
   389)     this%CN_ratio(icount) = cur_pool%CN_ratio
   390)     pool_names(icount) = cur_pool%name
   391)     if (cur_pool%CN_ratio < 0.d0) then
   392)       ! Since no CN ratio provided, must provide two species with the
   393)       ! same name as the pool with C or N appended.
   394)       word = trim(cur_pool%name) // 'C'
   395)       this%pool_id_to_species_id(CARBON_INDEX,icount) = &
   396)         GetImmobileSpeciesIDFromName(word,reaction%immobile, &
   397)                                      PETSC_FALSE,option)
   398)       word = trim(cur_pool%name) // 'N'
   399)       this%pool_id_to_species_id(NITROGEN_INDEX,icount) = &
   400)         GetImmobileSpeciesIDFromName(word,reaction%immobile, &
   401)                                      PETSC_FALSE,option)
   402)       this%pool_id_to_species_id(0,icount) = 2
   403)       if (minval(this%pool_id_to_species_id(:,icount)) <= 0) then
   404)         option%io_buffer = 'For CLM-CN pools with no CN ratio defined, ' // &
   405)           'the user must define two immobile species with the same root ' // &
   406)           'name as the pool with "C" or "N" appended, respectively.'
   407)         call printErrMsg(option)
   408)       endif
   409)     else ! only one species (e.g. SOMX)
   410)       this%pool_id_to_species_id(SOM_INDEX,icount) = &
   411)         GetImmobileSpeciesIDFromName(cur_pool%name,reaction%immobile, &
   412)                                      PETSC_TRUE,option)
   413)       this%pool_id_to_species_id(0,icount) = 1
   414)     endif
   415)     cur_pool => cur_pool%next
   416)   enddo
   417)   
   418)   ! map C and N species (solid phase for now)
   419)   word = 'C'
   420)   this%C_species_id = &
   421)       GetImmobileSpeciesIDFromName(word,reaction%immobile, &
   422)                                    PETSC_TRUE,option)
   423)   word = 'N'
   424)   this%N_species_id = &
   425)       GetImmobileSpeciesIDFromName(word,reaction%immobile, &
   426)                                    PETSC_TRUE,option)
   427)   
   428)   icount = 0
   429)   cur_rxn => this%reactions
   430)   do
   431)     if (.not.associated(cur_rxn)) exit
   432)     icount = icount + 1
   433)     this%upstream_pool_id(icount) = &
   434)       StringFindEntryInList(cur_rxn%upstream_pool_name,pool_names)
   435)     if (this%upstream_pool_id(icount) == 0) then
   436)       option%io_buffer = 'Upstream pool "' // &
   437)         trim(cur_rxn%upstream_pool_name) // &
   438)         '" in reaction with downstream pool "' // &
   439)         trim(cur_rxn%downstream_pool_name) // '" not found in list of pools.'
   440)       call printErrMsg(option)
   441)     endif
   442)     if (len_trim(cur_rxn%downstream_pool_name) > 0) then
   443)       this%downstream_pool_id(icount) = &
   444)         StringFindEntryInList(cur_rxn%downstream_pool_name,pool_names)
   445)       if (this%downstream_pool_id(icount) == 0) then
   446)         option%io_buffer = 'Downstream pool "' // &
   447)           trim(cur_rxn%downstream_pool_name) // &
   448)           '" in reaction with upstream pool "' // &
   449)           trim(cur_rxn%upstream_pool_name) // '" not found in list of pools.'
   450)         call printErrMsg(option)
   451)       endif
   452)       if (this%CN_ratio(this%downstream_pool_id(icount)) < 0.d0) then
   453)         option%io_buffer = 'For CLM-CN reactions, downstream pools ' // &
   454)           'must have a constant C:N ratio (i.e. C and N are not tracked ' // &
   455)           ' individually.  Therefore, pool "' // &
   456)           trim(cur_rxn%downstream_pool_name) // &
   457)           '" may not be used as a downstream pool.'
   458)         call printErrMsg(option)
   459)       endif
   460)     endif
   461)     this%rate_constant(icount) = cur_rxn%rate_constant
   462)     this%respiration_fraction(icount) = cur_rxn%respiration_fraction
   463)     this%inhibition_constant(icount) = cur_rxn%inhibition_constant
   464)     cur_rxn => cur_rxn%next
   465)   enddo 
   466)   
   467)   deallocate(pool_names)
   468)   
   469) end subroutine CLM_CN_Map
   470) 
   471) ! ************************************************************************** !
   472) 
   473) subroutine CLM_CN_React(this,Residual,Jacobian,compute_derivative,rt_auxvar, &
   474)                         global_auxvar,material_auxvar,reaction,option)
   475)   ! 
   476)   ! Evaluates reaction storing residual and/or Jacobian
   477)   ! 
   478)   ! Author: Glenn Hammond
   479)   ! Date: 02/04/13
   480)   ! 
   481) 
   482)   use Option_module
   483)   use Reaction_Aux_module, only : reaction_type
   484)   use Material_Aux_class, only : material_auxvar_type
   485)   
   486)   implicit none
   487)   
   488)   class(reaction_sandbox_clm_cn_type) :: this
   489)   type(option_type) :: option
   490)   type(reaction_type) :: reaction
   491)   PetscBool :: compute_derivative
   492)   ! the following arrays must be declared after reaction
   493)   PetscReal :: Residual(reaction%ncomp)
   494)   PetscReal :: Jacobian(reaction%ncomp,reaction%ncomp)
   495)   type(reactive_transport_auxvar_type) :: rt_auxvar
   496)   type(global_auxvar_type) :: global_auxvar
   497)   class(material_auxvar_type) :: material_auxvar
   498) 
   499)   PetscInt, parameter :: iphase = 1
   500)   PetscInt :: ipool_up, ipool_down
   501)   PetscInt :: ispec_pool_down
   502)   PetscInt :: ispecC_pool_up, ispecN_pool_up
   503)   PetscInt :: ires_pool_down, ires_C, ires_N
   504)   PetscInt :: iresC_pool_up, iresN_pool_up
   505)   PetscInt :: ispec_N
   506)   PetscReal :: drate, scaled_rate_const, rate
   507)   PetscInt :: irxn
   508)   
   509)   ! inhibition variables
   510)   PetscReal :: F_t
   511)   PetscReal :: F_theta
   512)   PetscReal :: constant_inhibition
   513)   PetscReal :: temp_K
   514)   PetscReal, parameter :: one_over_71_02 = 1.408054069d-2
   515)   PetscReal, parameter :: theta_min = 0.01d0     ! 1/nat log(0.01d0)
   516)   PetscReal, parameter :: one_over_log_theta_min = -2.17147241d-1
   517)   PetscReal, parameter :: twelve_over_14 = 0.857142857143d0
   518) 
   519)   PetscReal :: CN_ratio_up, CN_ratio_down
   520)   PetscBool :: constant_CN_ratio_up
   521)   PetscReal :: resp_frac
   522)   PetscReal :: stoich_N
   523)   PetscReal :: stoich_C
   524)   PetscReal :: stoich_downstreamC_pool
   525)   PetscReal :: stoich_upstreamC_pool, stoich_upstreamN_pool
   526)   
   527)   PetscReal :: N_inhibition, d_N_inhibition
   528)   PetscReal :: drate_dN_inhibition
   529)   PetscBool :: use_N_inhibition
   530)   PetscReal :: temp_real
   531)   
   532)   PetscReal :: dCN_ratio_up_dC_pool_up, dCN_ratio_up_dN_pool_up
   533)   PetscReal :: dstoich_upstreamN_pool_dC_pool_up
   534)   PetscReal :: dstoich_upstreamN_pool_dN_pool_up
   535)   PetscReal :: dstoichN_dC_pool_up
   536)   PetscReal :: dstoichN_dN_pool_up
   537)   
   538)   PetscReal :: sumC
   539)   PetscReal :: sumN
   540)   
   541)   ! inhibition due to temperature
   542)   ! Equation: F_t = exp(308.56*(1/17.02 - 1/(T - 227.13)))
   543)   temp_K = global_auxvar%temp + 273.15d0
   544) 
   545)   if (temp_K > 227.15d0) then
   546)     F_t = exp(308.56d0*(one_over_71_02 - 1.d0/(temp_K - 227.13d0)))
   547)   else
   548)     F_t = 0.0
   549)     return
   550)   endif
   551)   
   552)   ! inhibition due to moisture content
   553)   ! Equation: F_theta = log(theta_min/theta) / log(theta_min/theta_max)
   554)   ! Assumptions: theta is saturation
   555)   !              theta_min = 0.01, theta_max = 1.
   556)   F_theta = log(theta_min/max(theta_min,global_auxvar%sat(1))) * one_over_log_theta_min 
   557)   
   558)   constant_inhibition = F_t * F_theta
   559)   
   560)   ! indices for C and N species
   561)   ires_C = reaction%offset_immobile + this%C_species_id
   562)   ispec_N = this%N_species_id
   563)   ires_N = reaction%offset_immobile + ispec_N
   564) 
   565)   do irxn = 1, this%nrxn
   566)   
   567)     sumC = 0.d0
   568)     sumN = 0.d0
   569)   
   570)     ! scaled_rate_const units: (m^3 bulk / s) = (1/s) * (m^3 bulk)
   571)     scaled_rate_const = this%rate_constant(irxn)*material_auxvar%volume* &
   572)                         constant_inhibition
   573)     resp_frac = this%respiration_fraction(irxn)
   574)     
   575)     ipool_up = this%upstream_pool_id(irxn)
   576)     constant_CN_ratio_up = (this%pool_id_to_species_id(0,ipool_up) == 1)
   577) 
   578)     if (.not. constant_CN_ratio_up) then
   579)       ! upstream pool is Litter pool with two species (C,N)
   580)       !
   581)       !  a LitC_i + b LitN_i -> c SOM_i + d C + e N
   582)       !
   583)       ispecC_pool_up = this%pool_id_to_species_id(CARBON_INDEX,ipool_up)
   584)       ispecN_pool_up = this%pool_id_to_species_id(NITROGEN_INDEX,ipool_up)
   585)       CN_ratio_up = rt_auxvar%immobile(ispecC_pool_up) / &
   586)                     rt_auxvar%immobile(ispecN_pool_up)
   587)     else
   588)       ! upstream pool is an SOM pool with one species
   589)       !
   590)       !  a SOM_i -> c SOM_(i+1) + d C + e N
   591)       !
   592)       ispecC_pool_up = this%pool_id_to_species_id(SOM_INDEX,ipool_up)
   593)       CN_ratio_up = this%CN_ratio(ipool_up)
   594)     endif
   595) 
   596)     ! a = fraction_C_up = 1.
   597)     stoich_upstreamC_pool = 1.d0
   598)     ! b = a / CN_ratio_up
   599)     stoich_upstreamN_pool = stoich_upstreamC_pool / CN_ratio_up
   600) 
   601)     ! downstream pool
   602)     ipool_down = this%downstream_pool_id(irxn)
   603)     ! a downstream pool need not exist if last in succession and
   604)     ! respiration fraction = 1.
   605)     if (ipool_down > 0) then
   606)       ! downstream pool is always SOM (i.e. not split between C and N
   607)       ! as a litter would be).
   608)       ispec_pool_down = this%pool_id_to_species_id(SOM_INDEX,ipool_down)
   609)       CN_ratio_down = this%CN_ratio(ipool_down)
   610)       ! c = (1-resp_frac) * a
   611)       stoich_downstreamC_pool = (1.d0-resp_frac) * stoich_upstreamC_pool
   612)     else    
   613)       ispec_pool_down = 0
   614)       stoich_downstreamC_pool = 0.d0
   615)       CN_ratio_down = 1.d0 ! to prevent divide by zero below.
   616)     endif
   617)       
   618)     ! d = resp_frac * a
   619)     stoich_C = resp_frac * stoich_upstreamC_pool
   620)     ! e = b - c / CN_ratio_dn
   621)     stoich_N = stoich_upstreamN_pool - stoich_downstreamC_pool / CN_ratio_down
   622)  
   623)     ! Inhibition by nitrogen (inhibition concentration > 0 and N is a reactant)
   624)     ! must be calculated here as the sign on the stoichiometry for N is 
   625)     ! required.
   626)     if (this%inhibition_constant(irxn) > 1.d-40 .and. stoich_N < 0.d0) then
   627)       use_N_inhibition = PETSC_TRUE
   628)       temp_real = rt_auxvar%immobile(ispec_N) + &
   629)                   this%inhibition_constant(irxn)
   630)       N_inhibition = rt_auxvar%immobile(ispec_N) / temp_real
   631)       d_N_inhibition = this%inhibition_constant(irxn) / &
   632)                              (temp_real * temp_real)
   633)     else 
   634)       use_N_inhibition = PETSC_FALSE
   635)       N_inhibition = 1.d0
   636)       d_N_inhibition = 0.d0
   637)     endif
   638)     
   639)     ! residual units: (mol/sec) = (m^3 bulk/s) * (mol/m^3 bulk)
   640)     rate = rt_auxvar%immobile(ispecC_pool_up) * &
   641)            scaled_rate_const * N_inhibition
   642) 
   643)     ! calculation of residual
   644)     
   645)     ! carbon
   646)     Residual(ires_C) = Residual(ires_C) - stoich_C * rate
   647)     sumC = sumC + stoich_C * rate
   648)     
   649)     ! nitrogen
   650)     Residual(ires_N) = Residual(ires_N) - stoich_N * rate
   651)     sumN = sumN + stoich_N * rate
   652) 
   653)     ! C species in upstream pool (litter or SOM)
   654)     iresC_pool_up = reaction%offset_immobile + ispecC_pool_up
   655)     ! scaled by negative one since it is a reactant 
   656)     Residual(iresC_pool_up) = Residual(iresC_pool_up) - &
   657)       (-1.d0) * stoich_upstreamC_pool * rate
   658)     sumC = sumC - stoich_upstreamC_pool * rate
   659)     
   660)     ! N species in upstream pool (litter only)
   661)     if (.not.constant_CN_ratio_up) then
   662)       ! N species in upstream pool
   663)       iresN_pool_up = reaction%offset_immobile + ispecN_pool_up
   664)       ! scaled by negative one since it is a reactant 
   665)       Residual(iresN_pool_up) = Residual(iresN_pool_up) - &
   666)         (-1.d0) * stoich_upstreamN_pool * rate
   667)     endif
   668)     sumN = sumN - stoich_upstreamN_pool * rate
   669)     
   670)     if (ispec_pool_down > 0) then
   671)       ! downstream pool
   672)       ires_pool_down = reaction%offset_immobile + ispec_pool_down
   673)       Residual(ires_pool_down) = Residual(ires_pool_down) - &
   674)         stoich_downstreamC_pool * rate
   675)       sumC = sumC + stoich_downstreamC_pool * rate
   676)       sumN = sumN + stoich_downstreamC_pool / CN_ratio_down * rate
   677)     endif
   678)     
   679)     !for debugging
   680) !    if (dabs(sumC) > 1.d-40 .or. dabs(sumN) > 1.d-40) then
   681) !      print *, "sum C: ", sumC
   682) !      print *, "sum N: ", sumN
   683) !    endif
   684)     
   685)     if (compute_derivative) then
   686)     
   687)       drate = scaled_rate_const * N_inhibition
   688)       
   689)       ! upstream C pool
   690)       Jacobian(iresC_pool_up,iresC_pool_up) = &
   691)         Jacobian(iresC_pool_up,iresC_pool_up) - &
   692)         ! scaled by negative one since it is a reactant 
   693)         (-1.d0) * stoich_upstreamC_pool * drate
   694)       if (use_N_inhibition) then
   695) !       revision to avoid division by 0 when N -> 0 (N_inhibition -> 0)
   696)         drate_dN_inhibition = rt_auxvar%immobile(ispecC_pool_up) * &
   697)            scaled_rate_const * d_N_inhibition
   698) 
   699)         ! scaled by negative one since it is a reactant 
   700)         Jacobian(iresC_pool_up,ires_N) = &
   701)           Jacobian(iresC_pool_up,ires_N) - &
   702)           (-1.d0) * stoich_upstreamC_pool * drate_dN_inhibition
   703)       endif
   704)       
   705)       ! downstream pool
   706)       if (ispec_pool_down > 0) then
   707)         Jacobian(ires_pool_down,iresC_pool_up) = &
   708)           Jacobian(ires_pool_down,iresC_pool_up) - &
   709)           stoich_downstreamC_pool * drate
   710)         if (use_N_inhibition) then
   711)           Jacobian(ires_pool_down,ires_N) = &
   712)             Jacobian(ires_pool_down,ires_N) - &
   713)             stoich_downstreamC_pool * drate_dN_inhibition
   714)         endif
   715)       endif
   716)       
   717)       ! variable upstream N pool (for litter pools only!!!)
   718)       if (.not.constant_CN_ratio_up) then
   719) 
   720)         ! derivative of upstream N pool with respect to upstream C pool
   721)         ! scaled by negative one since it is a reactant 
   722)         Jacobian(iresN_pool_up,iresC_pool_up) = &
   723)           Jacobian(iresN_pool_up,iresC_pool_up) - &
   724)           (-1.d0) * stoich_upstreamN_pool * drate
   725)         if (use_N_inhibition) then
   726)           ! scaled by negative one since it is a reactant 
   727)           Jacobian(iresN_pool_up,ires_N) = &
   728)             Jacobian(iresN_pool_up,ires_N) - &
   729)             (-1.d0) * stoich_upstreamN_pool * drate_dN_inhibition
   730)         endif
   731) 
   732)         ! Since CN_ratio is a function of unknowns, other derivatives to 
   733)         ! calculate
   734)         dCN_ratio_up_dC_pool_up = 1.d0 / rt_auxvar%immobile(ispecN_pool_up)
   735)         dCN_ratio_up_dN_pool_up = -1.d0 * CN_ratio_up / &
   736)                                   rt_auxvar%immobile(ispecN_pool_up)
   737)         
   738)         ! stoich_upstreamC_pool = 1.d0 (for upstream litter pool)
   739)         ! dstoich_upstreamC_pool_dC_up = 0.
   740)         ! stoich_upstreamN_pool = stoich_upstreamC_pool / CN_ratio_up
   741)         temp_real = -1.d0 * stoich_upstreamN_pool / CN_ratio_up
   742)         dstoich_upstreamN_pool_dC_pool_up = temp_real * dCN_ratio_up_dC_pool_up
   743)         dstoich_upstreamN_pool_dN_pool_up = temp_real * dCN_ratio_up_dN_pool_up        
   744) 
   745)         Jacobian(iresN_pool_up,iresC_pool_up) = &
   746)           Jacobian(iresN_pool_up,iresC_pool_up) - &
   747)           ! scaled by negative one since it is a reactant 
   748) !     revision to avoid division by 0 when upstream N -> 0 (line 727, 734, 735)
   749)           (-1.d0) * (-1.d0) * rt_auxvar%immobile(ispecN_pool_up)/ &
   750)                               rt_auxvar%immobile(ispecC_pool_up)* &
   751)                               scaled_rate_const * N_inhibition
   752) 
   753)         Jacobian(iresN_pool_up,iresN_pool_up) = &
   754)           Jacobian(iresN_pool_up,iresN_pool_up) - &
   755)           ! scaled by negative one since it is a reactant 
   756) !     revision to avoid division by 0 when upstream N -> 0 (line 727, 734, 735)
   757)           (-1.d0) * scaled_rate_const * N_inhibition
   758) 
   759)         ! stoich_C = resp_frac * stoich_upstreamC_pool
   760)         ! dstoichC_dC_pool_up = 0.
   761)         ! dstoichC_dN_pool_up = 0
   762) 
   763)         ! nitrogen (stoichiometry a function of upstream C/N)
   764)         ! stoich_N = stoich_upstreamN_pool - stoich_downstreamC_pool / &
   765)         !                                    CN_ratio_down
   766)         ! latter half is constant
   767)         dstoichN_dC_pool_up = dstoich_upstreamN_pool_dC_pool_up
   768)         dstoichN_dN_pool_up = dstoich_upstreamN_pool_dN_pool_up
   769)         Jacobian(ires_N,iresC_pool_up) = Jacobian(ires_N,iresC_pool_up) - &
   770) !     revision to avoid division by 0 when upstream N -> 0 (line 727, 734, 735)
   771)           (-1.d0) * rt_auxvar%immobile(ispecN_pool_up)/ &
   772)                     rt_auxvar%immobile(ispecC_pool_up)* &
   773)                     scaled_rate_const * N_inhibition
   774) 
   775)         Jacobian(ires_N,iresN_pool_up) = Jacobian(ires_N,iresN_pool_up) - &
   776) !     revision to avoid division by 0 when upstream N -> 0 (line 727, 734, 735)
   777)           scaled_rate_const * N_inhibition
   778)       endif
   779)       
   780)       ! carbon
   781)       Jacobian(ires_C,iresC_pool_up) = Jacobian(ires_C,iresC_pool_up) - &
   782)         stoich_C * drate
   783)       ! nitrogen
   784)       Jacobian(ires_N,iresC_pool_up) = Jacobian(ires_N,iresC_pool_up) - &
   785)         stoich_N * drate
   786)       if (use_N_inhibition) then
   787)         Jacobian(ires_C,ires_N) = Jacobian(ires_C,ires_N) - & 
   788)           stoich_C * drate_dN_inhibition
   789)         Jacobian(ires_N,ires_N) = Jacobian(ires_N,ires_N) - &
   790)           stoich_N * drate_dN_inhibition
   791)       endif
   792)     endif
   793)   enddo
   794)   
   795) end subroutine CLM_CN_React
   796) 
   797) ! ************************************************************************** !
   798) 
   799) subroutine CLM_CN_Destroy(this)
   800)   ! 
   801)   ! Destroys allocatable or pointer objects created in this
   802)   ! module
   803)   ! 
   804)   ! Author: Glenn Hammond
   805)   ! Date: 02/04/13
   806)   ! 
   807) 
   808)   use Utility_module, only : DeallocateArray
   809) 
   810)   implicit none
   811)   
   812)   class(reaction_sandbox_clm_cn_type) :: this
   813)   
   814)   type(pool_type), pointer :: cur_pool, prev_pool
   815)   type(clm_cn_reaction_type), pointer :: cur_reaction, prev_reaction
   816)   
   817)   cur_pool => this%pools
   818)   do
   819)     if (.not.associated(cur_pool)) exit
   820)     prev_pool => cur_pool
   821)     cur_pool => cur_pool%next
   822)     deallocate(prev_pool)
   823)     nullify(prev_pool)
   824)   enddo
   825)   
   826)   cur_reaction => this%reactions
   827)   do
   828)     if (.not.associated(cur_reaction)) exit
   829)     prev_reaction => cur_reaction
   830)     cur_reaction => cur_reaction%next
   831)     deallocate(prev_reaction)
   832)     nullify(prev_reaction)
   833)   enddo
   834)   
   835)   call DeallocateArray(this%CN_ratio)
   836)   call DeallocateArray(this%rate_constant)
   837)   call DeallocateArray(this%respiration_fraction)
   838)   call DeallocateArray(this%inhibition_constant)
   839)   call DeallocateArray(this%upstream_pool_id)
   840)   call DeallocateArray(this%downstream_pool_id)
   841)   call DeallocateArray(this%pool_id_to_species_id)
   842)   
   843) end subroutine CLM_CN_Destroy
   844) 
   845) end module Reaction_Sandbox_CLM_CN_class

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