reaction_surf_complex.F90       coverage:  85.71 %func     46.48 %block


     1) module Reaction_Surface_Complexation_module
     2) 
     3)   use Reaction_Surface_Complexation_Aux_module
     4)   use Reaction_Aux_module
     5)   use Reactive_Transport_Aux_module
     6)   use Global_Aux_module
     7)   use Material_Aux_class
     8)   
     9)   use PFLOTRAN_Constants_module
    10) 
    11)   implicit none
    12)   
    13)   private 
    14) 
    15) #include "petsc/finclude/petscsys.h"
    16) 
    17)   PetscReal, parameter :: perturbation_tolerance = 1.d-5
    18)   
    19)   public :: SurfaceComplexationRead, &
    20)             SrfCplxProcessConstraint, &
    21)             RTotalSorbEqSurfCplx, &
    22)             RMultiRateSorption, &
    23)             RKineticSurfCplx, &
    24)             RTotalSorbMultiRateAsEQ
    25)             
    26) contains
    27) 
    28) ! ************************************************************************** !
    29) 
    30) subroutine SurfaceComplexationRead(reaction,input,option)
    31)   ! 
    32)   ! Reads chemical species
    33)   ! 
    34)   ! Author: Glenn Hammond
    35)   ! Date: 05/02/08
    36)   ! 
    37) 
    38)   use Option_module
    39)   use String_module
    40)   use Input_Aux_module
    41)   use Utility_module
    42)   
    43)   implicit none
    44)   
    45)   type(reaction_type) :: reaction
    46)   type(input_type), pointer :: input
    47)   type(option_type) :: option
    48)   
    49)   character(len=MAXSTRINGLENGTH) :: string
    50)   character(len=MAXWORDLENGTH) :: word
    51)   character(len=MAXWORDLENGTH) :: name
    52)   character(len=MAXWORDLENGTH) :: card
    53)   type(surface_complexation_type), pointer :: surface_complexation
    54)   type(surface_complex_type), pointer :: srfcplx, cur_srfcplx, prev_srfcplx, &
    55)                                          cur_srfcplx_in_rxn
    56)   type(surface_complex_type), pointer :: rate_list, cur_srfcplx_rate, &
    57)                                          prev_srfcplx_rate  
    58)   type(surface_complexation_rxn_type), pointer :: srfcplx_rxn, cur_srfcplx_rxn
    59)   PetscInt :: temp_srfcplx_count
    60)   PetscBool :: found
    61)   PetscReal :: tempreal
    62)   PetscInt :: i
    63)   PetscInt :: num_times_surface_type_set
    64)   
    65)   nullify(srfcplx_rxn)
    66)   nullify(cur_srfcplx_rxn)
    67)   nullify(cur_srfcplx)
    68)   
    69)   surface_complexation => reaction%surface_complexation
    70)            
    71)   srfcplx_rxn => SurfaceComplexationRxnCreate()
    72)   ! default
    73)   srfcplx_rxn%itype = SRFCMPLX_RXN_EQUILIBRIUM
    74)   temp_srfcplx_count = 0
    75)   num_times_surface_type_set = 0
    76)   do
    77)     call InputReadPflotranString(input,option)
    78)     if (InputError(input)) exit
    79)     if (InputCheckExit(input,option)) exit
    80) 
    81)     call InputReadWord(input,option,word,PETSC_TRUE)
    82)     call InputErrorMsg(input,option,'keyword', &
    83)                         'CHEMISTRY,SURFACE_COMPLEXATION_RXN')
    84)     call StringToUpper(word)
    85)                 
    86)     select case(trim(word))
    87)       case('EQUILIBRIUM')
    88)         srfcplx_rxn%itype = SRFCMPLX_RXN_EQUILIBRIUM
    89)       case('MULTIRATE_KINETIC')
    90)         srfcplx_rxn%itype = SRFCMPLX_RXN_MULTIRATE_KINETIC
    91)       case('KINETIC')
    92)         srfcplx_rxn%itype = SRFCMPLX_RXN_KINETIC
    93)       case('COMPLEX_KINETICS')
    94)         nullify(prev_srfcplx)
    95)         do
    96)           call InputReadPflotranString(input,option)
    97)           if (InputError(input)) exit
    98)           if (InputCheckExit(input,option)) exit
    99)                       
   100)           srfcplx => SurfaceComplexCreate()
   101)           call InputReadWord(input,option,srfcplx%name,PETSC_TRUE)
   102)           call InputErrorMsg(input,option,'keyword', &
   103)             'CHEMISTRY,SURFACE_COMPLEXATION_RXN,COMPLEX_KINETIC_RATE')
   104)                         
   105)           do
   106)             call InputReadPflotranString(input,option)
   107)             call InputReadStringErrorMsg(input,option,card)
   108)             if (InputCheckExit(input,option)) exit
   109)             call InputReadWord(input,option,word,PETSC_TRUE)
   110)             call InputErrorMsg(input,option,'word', &
   111)                     'CHEMISTRY,SURFACE_COMPLEXATION_RXN,COMPLEX_KINETIC_RATE') 
   112)             select case(trim(word))
   113)               case('FORWARD_RATE_CONSTANT')
   114)                 call InputReadDouble(input,option,srfcplx%forward_rate)
   115)                 call InputErrorMsg(input,option,'forward_rate', &
   116)                         'CHEMISTRY,SURFACE_COMPLEXATION_RXN,COMPLEX_KINETIC_RATE')
   117)               case('BACKWARD_RATE_CONSTANT')
   118)                 call InputReadDouble(input,option,srfcplx%backward_rate)
   119)                 call InputErrorMsg(input,option,'backward_rate', &
   120)                         'CHEMISTRY,SURFACE_COMPLEXATION_RXN,COMPLEX_KINETIC_RATE')
   121)               case default
   122)                 call InputKeywordUnrecognized(word, &
   123)                        'CHEMISTRY,SURFACE_COMPLEXATION_RXN,COMPLEX_KINETIC_RATE',option)
   124)             end select
   125)           enddo
   126)                                       
   127)           if (.not.associated(rate_list)) then
   128)             rate_list => srfcplx
   129)           endif
   130)           if (associated(prev_srfcplx)) then
   131)             prev_srfcplx%next => srfcplx
   132)           endif
   133)           prev_srfcplx => srfcplx
   134)           nullify(srfcplx)
   135)         enddo
   136)         nullify(prev_srfcplx)
   137)       case('RATE','RATES') 
   138)         srfcplx_rxn%itype = SRFCMPLX_RXN_MULTIRATE_KINETIC
   139)         string = 'RATES inside SURFACE_COMPLEXATION_RXN'
   140)         call UtilityReadArray(srfcplx_rxn%rates,NEG_ONE_INTEGER,string,input, &
   141)                               option) 
   142)       case('SITE_FRACTION') 
   143)         string = 'SITE_FRACTION inside SURFACE_COMPLEXATION_RXN'
   144)         call UtilityReadArray(srfcplx_rxn%site_fractions,NEG_ONE_INTEGER, &
   145)                               string,input,option) 
   146)       case('MULTIRATE_SCALE_FACTOR')
   147)         call InputReadDouble(input,option,srfcplx_rxn%kinmr_scale_factor)
   148)         call InputErrorMsg(input,option,'keyword', &
   149)           'CHEMISTRY,SURFACE_COMPLEXATION_RXN,MULTIRATE_SCALE_FACTOR')
   150)       case('MINERAL')
   151)         srfcplx_rxn%surface_itype = MINERAL_SURFACE
   152)         num_times_surface_type_set = num_times_surface_type_set + 1
   153)         call InputReadWord(input,option,srfcplx_rxn%surface_name, &
   154)           PETSC_TRUE)
   155)         call InputErrorMsg(input,option,'keyword', &
   156)           'CHEMISTRY,SURFACE_COMPLEXATION_RXN,MINERAL_NAME')
   157)       case('ROCK_DENSITY')
   158)         srfcplx_rxn%surface_itype = ROCK_SURFACE
   159)         num_times_surface_type_set = num_times_surface_type_set + 1
   160)       case('COLLOID')
   161)         srfcplx_rxn%surface_itype = COLLOID_SURFACE
   162)         num_times_surface_type_set = num_times_surface_type_set + 1
   163)         call InputReadWord(input,option,srfcplx_rxn%surface_name, &
   164)           PETSC_TRUE)
   165)         call InputErrorMsg(input,option,'keyword', &
   166)           'CHEMISTRY,SURFACE_COMPLEXATION_RXN,COLLOID_NAME')
   167)       case('SITE')
   168)         call InputReadWord(input,option,srfcplx_rxn%free_site_name, &
   169)           PETSC_TRUE)
   170)         call InputErrorMsg(input,option,'keyword', &
   171)           'CHEMISTRY,SURFACE_COMPLEXATION_RXN,SITE_NAME')
   172)         ! site density in mol/m^3 bulk
   173)         call InputReadDouble(input,option,srfcplx_rxn%site_density)
   174)         call InputErrorMsg(input,option,'keyword', &
   175)           'CHEMISTRY,SURFACE_COMPLEXATION_RXN,SITE_DENSITY')                   
   176)       case('COMPLEXES')
   177)         nullify(prev_srfcplx)
   178)         do
   179)           call InputReadPflotranString(input,option)
   180)           if (InputError(input)) exit
   181)           if (InputCheckExit(input,option)) exit
   182)                       
   183)           temp_srfcplx_count = temp_srfcplx_count + 1
   184)           srfcplx => SurfaceComplexCreate()
   185)           srfcplx%id = temp_srfcplx_count
   186)           call InputReadWord(input,option,srfcplx%name,PETSC_TRUE)
   187)           call InputErrorMsg(input,option,'keyword', &
   188)             'CHEMISTRY,SURFACE_COMPLEXATION_RXN,COMPLEX_NAME')
   189)                 
   190)           if (.not.associated(srfcplx_rxn%complex_list)) then
   191)             srfcplx_rxn%complex_list => srfcplx
   192)           endif
   193)           if (associated(prev_srfcplx)) then
   194)             prev_srfcplx%next => srfcplx
   195)           endif
   196)           prev_srfcplx => srfcplx
   197)           nullify(srfcplx)
   198)         enddo
   199)         nullify(prev_srfcplx)
   200)       case default
   201)         call InputKeywordUnrecognized(word, &
   202)                 'CHEMISTRY,SURFACE_COMPLEXATION_RXN',option)
   203)     end select
   204) 
   205)   enddo
   206)   
   207)   if (num_times_surface_type_set > 1) then
   208)     option%io_buffer = 'Surface site type (e.g. MINERAL, ROCK_DENSITY, ' // &
   209)       'COLLOID) may only be set once under the SURFACE_COMPLEXATION_RXN card.'
   210)     call printErrMsg(option)
   211)   endif
   212)   
   213)   if (.not.associated(surface_complexation%rxn_list)) then
   214)     surface_complexation%rxn_list => srfcplx_rxn
   215)     srfcplx_rxn%id = 1
   216)   else
   217)     cur_srfcplx_rxn => surface_complexation%rxn_list
   218)     do
   219)       if (.not.associated(cur_srfcplx_rxn%next)) then
   220)         cur_srfcplx_rxn%next => srfcplx_rxn
   221)         srfcplx_rxn%id = cur_srfcplx_rxn%id + 1
   222)         exit
   223)       endif
   224)       cur_srfcplx_rxn => cur_srfcplx_rxn%next
   225)     enddo
   226)     nullify(cur_srfcplx_rxn)
   227)   endif
   228)   
   229)   ! Add surface complexes in reaction to master list, without duplicating.
   230)   ! Set the id of the surface complex to the one in the master list, even
   231)   ! if duplicated.
   232)   cur_srfcplx_in_rxn => srfcplx_rxn%complex_list
   233)   do
   234)     if (.not.associated(cur_srfcplx_in_rxn)) exit
   235)     cur_srfcplx => surface_complexation%complex_list
   236)     found = PETSC_FALSE
   237)     do
   238)       if (.not.associated(cur_srfcplx)) exit
   239)       if (StringCompare(cur_srfcplx_in_rxn%name,cur_srfcplx%name, &
   240)                         MAXWORDLENGTH)) then
   241)         ! set id to id in master list
   242)         cur_srfcplx_in_rxn%id = cur_srfcplx%id
   243)         cur_srfcplx_in_rxn%ptr => cur_srfcplx
   244)         found = PETSC_TRUE
   245)         exit
   246)       endif
   247)       prev_srfcplx => cur_srfcplx
   248)       cur_srfcplx => cur_srfcplx%next
   249)     enddo
   250)     if (.not.found) then
   251)       srfcplx => SurfaceComplexCreate()
   252)       srfcplx%name = cur_srfcplx_in_rxn%name
   253)       if (.not.associated(prev_srfcplx)) then
   254)         surface_complexation%complex_list => srfcplx
   255)         srfcplx%id = 1
   256)       else
   257)         prev_srfcplx%next => srfcplx
   258)         srfcplx%id = prev_srfcplx%id + 1
   259)       endif
   260)     endif
   261)     cur_srfcplx_in_rxn => cur_srfcplx_in_rxn%next
   262)   enddo
   263) 
   264)   surface_complexation%nsrfcplxrxn = &
   265)     surface_complexation%nsrfcplxrxn + 1
   266)   select case(srfcplx_rxn%itype)
   267)     ! default (NULL) to EQUILIBRIUM
   268)     case(SRFCMPLX_RXN_NULL,SRFCMPLX_RXN_EQUILIBRIUM)
   269)       surface_complexation%neqsrfcplx = surface_complexation%neqsrfcplx + &
   270)         temp_srfcplx_count
   271)       surface_complexation%neqsrfcplxrxn = &
   272)         surface_complexation%neqsrfcplxrxn + 1
   273)     case(SRFCMPLX_RXN_MULTIRATE_KINETIC)
   274)       surface_complexation%nkinmrsrfcplx = &
   275)         surface_complexation%nkinmrsrfcplx + temp_srfcplx_count
   276)       surface_complexation%nkinmrsrfcplxrxn = &
   277)         surface_complexation%nkinmrsrfcplxrxn + 1
   278)       ! if site fractions is not specified, we can calculate this
   279)       ! based on a uniform distribution and the number of rates
   280)       if (.not.associated(srfcplx_rxn%site_fractions) .and. &
   281)           associated(srfcplx_rxn%rates)) then
   282)         allocate(srfcplx_rxn%site_fractions(size(srfcplx_rxn%rates)))
   283)         ! it is possible to specify the number of site fractions
   284)         srfcplx_rxn%site_fractions = 1.d0 / dble(size(srfcplx_rxn%rates))
   285)       endif
   286)       ! check to ensure that rates for multirate surface complexation 
   287)       ! are aligned with surface fractions
   288)       if (size(srfcplx_rxn%rates) /= size(srfcplx_rxn%site_fractions)) then
   289)         write(word,*) size(srfcplx_rxn%rates)
   290)         write(string,*) size(srfcplx_rxn%site_fractions)
   291)         option%io_buffer = 'Number of kinetic rates (' // &
   292)           trim(adjustl(word)) // &
   293)           ') does not match the number of surface fractions (' // &
   294)           trim(adjustl(string)) // ').'
   295)         call printErrMsg(option)
   296)       endif
   297)       tempreal = 0.d0
   298)       do i = 1, size(srfcplx_rxn%site_fractions)
   299)         tempreal = tempreal + srfcplx_rxn%site_fractions(i)
   300)         srfcplx_rxn%rates(i) = srfcplx_rxn%rates(i) * &
   301)           srfcplx_rxn%kinmr_scale_factor
   302)       enddo
   303)     
   304)       if (dabs(1.d0 - tempreal) > 1.d-6) then
   305)         write(string,*) tempreal
   306)         option%io_buffer = 'The sum of the surface site fractions for ' // &
   307)           'multirate kinetic sorption does not add up to 1.d0 (' // &
   308)           trim(adjustl(string)) // '.'
   309)         call printErrMsg(option)
   310)       endif
   311)     case(SRFCMPLX_RXN_KINETIC)
   312)       ! match up rates with their corresponding surface complex
   313)       cur_srfcplx => srfcplx_rxn%complex_list
   314)       do
   315)         if (.not.associated(cur_srfcplx)) exit
   316)         found = PETSC_FALSE
   317)         nullify(prev_srfcplx_rate)
   318)         cur_srfcplx_rate => rate_list
   319)         do
   320)           if (.not.associated(cur_srfcplx_rate)) exit
   321)           ! check for same name
   322)           if (StringCompare(cur_srfcplx_rate%name, &
   323)                             cur_srfcplx%name, &
   324)                             MAXWORDLENGTH)) then
   325)             ! set rates
   326)             cur_srfcplx%forward_rate = cur_srfcplx_rate%forward_rate
   327)             cur_srfcplx%backward_rate = cur_srfcplx_rate%backward_rate
   328)             ! remove srfcplx_rate from list of rates
   329)             if (associated(prev_srfcplx_rate)) then
   330)               prev_srfcplx_rate%next => cur_srfcplx_rate%next
   331)             else
   332)               rate_list => cur_srfcplx_rate%next
   333)             endif
   334)             ! destroy the object
   335)             call SurfaceComplexDestroy(cur_srfcplx_rate)
   336)             found = PETSC_TRUE
   337)             exit
   338)           endif
   339)           prev_srfcplx_rate => cur_srfcplx_rate
   340)           cur_srfcplx_rate => cur_srfcplx_rate%next
   341)         enddo
   342)         if (.not.found) then
   343)           option%io_buffer = 'Rates for surface complex ' // &
   344)             trim(cur_srfcplx%name) // ' not found in kinetic rate list'
   345)           call printErrMsg(option)
   346)         endif
   347)         cur_srfcplx => cur_srfcplx%next
   348)       enddo
   349)       ! check to ensure that rates are matched
   350)       if (associated(rate_list)) then
   351)         option%io_buffer = '# of rates is greater than # of surface complexes'
   352)         call printErrMsg(option)
   353)       endif
   354)       nullify(cur_srfcplx)
   355)       nullify(prev_srfcplx)
   356)       nullify(rate_list)
   357)       nullify(cur_srfcplx_rate)
   358)       nullify(prev_srfcplx_rate)                  
   359)       surface_complexation%nkinsrfcplx = &
   360)         surface_complexation%nkinsrfcplx + temp_srfcplx_count
   361)       surface_complexation%nkinsrfcplxrxn = &
   362)         surface_complexation%nkinsrfcplxrxn + 1
   363)   end select
   364)   srfcplx_rxn%free_site_id = srfcplx_rxn%id
   365)               
   366)   nullify(srfcplx_rxn)
   367) 
   368) end subroutine SurfaceComplexationRead
   369) 
   370) ! ************************************************************************** !
   371) 
   372) subroutine SrfCplxProcessConstraint(surface_complexation,constraint_name, &
   373)                                     constraint,option)
   374)   ! 
   375)   ! Initializes constraints based on surface complex
   376)   ! species in system
   377)   ! 
   378)   ! Author: Glenn Hammond
   379)   ! Date: 01/07/13
   380)   ! 
   381) 
   382)   use Option_module
   383)   use Input_Aux_module
   384)   use String_module
   385)   use Utility_module  
   386)   
   387)   implicit none
   388)   
   389)   type(surface_complexation_type), pointer :: surface_complexation
   390)   character(len=MAXWORDLENGTH) :: constraint_name
   391)   type(srfcplx_constraint_type), pointer :: constraint
   392)   type(option_type) :: option
   393)   
   394)   PetscBool :: found
   395)   PetscInt :: isrfcplx, jsrfcplx
   396)   
   397)   character(len=MAXWORDLENGTH) :: srfcplx_name(surface_complexation%nkinsrfcplx)
   398)   PetscReal :: constraint_conc(surface_complexation%nkinsrfcplx)
   399)   
   400)   if (.not.associated(constraint)) return
   401) 
   402)   if (surface_complexation%nkinsrfcplx == 0) then
   403)     option%io_buffer = 'Surface complexation specified in constraint "' // &
   404)       trim(constraint_name) // '" requires that kinetic surface ' // &
   405)       'complexation be defined in the CHEMISTRY section.'
   406)     call printErrMsg(option)
   407)   endif
   408)   
   409)   srfcplx_name = ''
   410)   do isrfcplx = 1, surface_complexation%nkinsrfcplx
   411)     found = PETSC_FALSE
   412)     do jsrfcplx = 1, surface_complexation%nkinsrfcplx
   413)       if (StringCompare(constraint%names(isrfcplx), &
   414)                         surface_complexation%srfcplx_names(&
   415)                           !TODO(geh): fix 0 index
   416)                         surface_complexation%kinsrfcplx_to_name(jsrfcplx,0)), &
   417)                         MAXWORDLENGTH)) then
   418)         found = PETSC_TRUE
   419)         exit
   420)       endif
   421)     enddo
   422)     if (.not.found) then
   423)       option%io_buffer = &
   424)                 'Surface complex ' // trim(constraint%names(isrfcplx)) // &
   425)                 'from CONSTRAINT ' // trim(constraint_name) // &
   426)                 ' not found among kinetic surface complexes.'
   427)       call printErrMsg(option)
   428)     else
   429)       constraint_conc(jsrfcplx) = &
   430)         constraint%constraint_conc(isrfcplx)
   431)       srfcplx_name(jsrfcplx) = constraint%names(isrfcplx)
   432)     endif  
   433)   enddo
   434)   constraint%names = srfcplx_name
   435)   constraint%constraint_conc = constraint_conc
   436)   
   437) end subroutine SrfCplxProcessConstraint
   438) 
   439) ! ************************************************************************** !
   440) 
   441) subroutine RTotalSorbEqSurfCplx(rt_auxvar,global_auxvar,material_auxvar, &
   442)                                 reaction,option)
   443)   ! 
   444)   ! Computes the total sorbed component concentrations and
   445)   ! derivative with respect to free-ion for equilibrium
   446)   ! surface complexation
   447)   ! 
   448)   ! Author: Glenn Hammond
   449)   ! Date: 10/22/08; 05/26/09
   450)   ! 
   451) 
   452)   use Option_module
   453)   use Matrix_Block_Aux_module
   454)   
   455)   implicit none
   456)   
   457)   type(reactive_transport_auxvar_type) :: rt_auxvar
   458)   type(global_auxvar_type) :: global_auxvar
   459)   class(material_auxvar_type) :: material_auxvar
   460)   type(reaction_type) :: reaction
   461)   type(option_type) :: option
   462)   
   463)   PetscInt :: irxn, ieqrxn
   464)   PetscReal, pointer :: colloid_array_ptr(:)
   465)   PetscInt, parameter :: iphase = 1
   466)   type(matrix_block_auxvar_type), pointer :: colloid_matrix_block_ptr
   467)   type(surface_complexation_type), pointer :: surface_complexation
   468) 
   469)   surface_complexation => reaction%surface_complexation
   470)   
   471)   if (reaction%ncollcomp > 0) then  
   472)     rt_auxvar%colloid%total_eq_mob = 0.d0
   473)     rt_auxvar%colloid%dRj_dCj%dtotal = 0.d0
   474)     colloid_array_ptr => rt_auxvar%colloid%total_eq_mob
   475)     colloid_matrix_block_ptr => rt_auxvar%colloid%dRj_dCj
   476)   else
   477)     nullify(colloid_matrix_block_ptr)
   478)     nullify(colloid_array_ptr)
   479)   endif
   480) 
   481)   ! Surface Complexation
   482)   do ieqrxn = 1, surface_complexation%neqsrfcplxrxn
   483)   
   484)     irxn = surface_complexation%eqsrfcplxrxn_to_srfcplxrxn(ieqrxn)
   485)     
   486)     !TODO(geh): clean up colloidpointers
   487)     call RTotalSorbEqSurfCplx1(rt_auxvar,global_auxvar,material_auxvar, &
   488)                                reaction,option, &
   489)                                irxn, &
   490)                                rt_auxvar%srfcplxrxn_free_site_conc(irxn), &
   491)                                rt_auxvar%eqsrfcplx_conc, &
   492)                                rt_auxvar%total_sorb_eq, &
   493)                                rt_auxvar%dtotal_sorb_eq, &
   494)                                colloid_array_ptr, &
   495)                                colloid_matrix_block_ptr)    
   496)   
   497)   enddo ! irxn
   498)   
   499)   ! units of total_sorb = mol/m^3
   500)   ! units of dtotal_sorb = kg water/m^3 bulk
   501)   
   502) end subroutine RTotalSorbEqSurfCplx
   503) 
   504) ! ************************************************************************** !
   505) 
   506) subroutine RTotalSorbMultiRateAsEQ(rt_auxvar,global_auxvar,material_auxvar, &
   507)                                    reaction,option)
   508)   ! 
   509)   ! Calculates the multirate surface complexation
   510)   ! reaction as if it were equilibrium.
   511)   ! 
   512)   ! Author: Glenn Hammond
   513)   ! Date: 03/19/12
   514)   ! 
   515) 
   516)   use Option_module
   517)   use Matrix_Block_Aux_module
   518)   
   519)   implicit none
   520)   
   521)   type(reactive_transport_auxvar_type) :: rt_auxvar
   522)   type(global_auxvar_type) :: global_auxvar
   523)   class(material_auxvar_type) :: material_auxvar
   524)   type(reaction_type) :: reaction
   525)   type(option_type) :: option
   526)   
   527)   PetscInt :: irxn, ikinmrrxn
   528)   PetscInt, parameter :: iphase = 1
   529)   type(surface_complexation_type), pointer :: surface_complexation
   530)   PetscReal :: total_sorb_eq(reaction%naqcomp)
   531)   PetscReal :: dtotal_sorb_eq(reaction%naqcomp,reaction%naqcomp)
   532)   PetscReal, pointer :: null_array_ptr(:)
   533)   type(matrix_block_auxvar_type), pointer :: null_matrix_block
   534) 
   535)   surface_complexation => reaction%surface_complexation
   536) 
   537)   nullify(null_array_ptr)
   538)   nullify(null_matrix_block)
   539) 
   540)   ! Surface Complexation
   541)   do ikinmrrxn = 1, surface_complexation%nkinmrsrfcplxrxn
   542)   
   543)     irxn = surface_complexation%kinmrsrfcplxrxn_to_srfcplxrxn(ikinmrrxn)
   544) 
   545)     total_sorb_eq = 0.d0
   546)     dtotal_sorb_eq = 0.d0  
   547) 
   548)     call RTotalSorbEqSurfCplx1(rt_auxvar,global_auxvar,material_auxvar, &
   549)                                reaction,option, &
   550)                                irxn, &
   551)                                rt_auxvar%srfcplxrxn_free_site_conc(irxn), &
   552)                                null_array_ptr, &
   553)                                total_sorb_eq, &
   554)                                dtotal_sorb_eq, &
   555)                                null_array_ptr, &
   556)                                null_matrix_block)    
   557) 
   558)     rt_auxvar%kinmr_total_sorb(:,0,ikinmrrxn) = total_sorb_eq(:)
   559)   
   560)   enddo ! irxn
   561)   
   562) end subroutine RTotalSorbMultiRateAsEQ
   563) 
   564) ! ************************************************************************** !
   565) 
   566) subroutine RMultiRateSorption(Res,Jac,compute_derivative,rt_auxvar, &
   567)                               global_auxvar,material_auxvar,reaction,option)
   568)   ! 
   569)   ! Computes contribution to the accumualtion term due
   570)   ! due to multirate sorption
   571)   ! 
   572)   ! Author: Glenn Hammond
   573)   ! Date: 05/20/09; 03/16/12
   574)   ! 
   575) 
   576)   use Option_module
   577)   use Matrix_Block_Aux_module
   578)   
   579)   implicit none
   580) 
   581)   PetscBool :: compute_derivative
   582)   type(reactive_transport_auxvar_type) :: rt_auxvar
   583)   type(global_auxvar_type) :: global_auxvar
   584)   class(material_auxvar_type) :: material_auxvar
   585)   type(reaction_type) :: reaction
   586)   PetscReal :: Res(reaction%ncomp)
   587)   PetscReal :: Jac(reaction%ncomp,reaction%ncomp)
   588)   type(option_type) :: option
   589)   
   590)   PetscInt :: irxn, ikinmrrxn
   591)   type(surface_complexation_type), pointer :: surface_complexation
   592)   
   593)   PetscInt :: irate
   594)   PetscInt, parameter :: iphase = 1
   595)   PetscReal :: kdt, one_plus_kdt, k_over_one_plus_kdt
   596)   PetscReal :: total_sorb_eq(reaction%naqcomp)
   597)   PetscReal :: dtotal_sorb_eq(reaction%naqcomp,reaction%naqcomp)
   598)   PetscReal, pointer :: null_array_ptr(:)
   599)   type(matrix_block_auxvar_type), pointer :: null_matrix_block
   600) 
   601)   surface_complexation => reaction%surface_complexation
   602) 
   603)   nullify(null_array_ptr)
   604)   nullify(null_matrix_block)
   605) 
   606)   ! only zero out the zero index.  The other indices hold values from 
   607)   ! the previous time step
   608)   rt_auxvar%kinmr_total_sorb(:,0,:) = 0.d0
   609)   
   610)   ! Surface Complexation
   611)   do ikinmrrxn = 1, surface_complexation%nkinmrsrfcplxrxn
   612)   
   613)     irxn = surface_complexation%kinmrsrfcplxrxn_to_srfcplxrxn(ikinmrrxn)
   614) 
   615)     total_sorb_eq = 0.d0
   616)     dtotal_sorb_eq = 0.d0  
   617) 
   618)     call RTotalSorbEqSurfCplx1(rt_auxvar,global_auxvar,material_auxvar, &
   619)                                reaction,option, &
   620)                                irxn, &
   621)                                rt_auxvar%srfcplxrxn_free_site_conc(irxn), &
   622)                                null_array_ptr, &
   623)                                total_sorb_eq, &
   624)                                dtotal_sorb_eq, &
   625)                                null_array_ptr, &
   626)                                null_matrix_block)
   627) 
   628)     ! WARNING: this assumes site fraction multiplicative factor 
   629)     do irate = 1, surface_complexation%kinmr_nrate(ikinmrrxn)
   630)       kdt = surface_complexation%kinmr_rate(irate,ikinmrrxn) * &
   631)             option%tran_dt
   632)       one_plus_kdt = 1.d0 + kdt
   633)       k_over_one_plus_kdt = &
   634)         surface_complexation%kinmr_rate(irate,ikinmrrxn)/one_plus_kdt
   635)         
   636)       ! this is the constribution to the accumulation term in the residual
   637)       Res(:) = Res(:) + material_auxvar%volume * k_over_one_plus_kdt * &
   638)         (surface_complexation%kinmr_frac(irate,ikinmrrxn)*total_sorb_eq(:) - &
   639)          rt_auxvar%kinmr_total_sorb(:,irate,ikinmrrxn))
   640)       
   641)       if (compute_derivative) then
   642)         Jac = Jac + material_auxvar%volume * k_over_one_plus_kdt * &
   643)           surface_complexation%kinmr_frac(irate,ikinmrrxn) * dtotal_sorb_eq
   644)       endif
   645) 
   646)     enddo
   647) 
   648)     ! store the target equilibrium concentration to update the sorbed 
   649)     ! concentration at the end of the time step.
   650)     rt_auxvar%kinmr_total_sorb(:,0,ikinmrrxn) = total_sorb_eq
   651)   
   652)   enddo ! ikinmrrxn
   653)   
   654) end subroutine RMultiRateSorption
   655) 
   656) ! ************************************************************************** !
   657) 
   658) subroutine RTotalSorbEqSurfCplx1(rt_auxvar,global_auxvar,material_auxvar, &
   659)                                  reaction,option, &
   660)                                  irxn,external_free_site_conc, &
   661)                                  external_srfcplx_conc, &
   662)                                  external_total_sorb, &
   663)                                  external_dtotal_sorb, &
   664)                                  external_total_mob, &
   665)                                  external_dRj_dCj)
   666)   ! 
   667)   ! Computes the total sorbed component concentrations and
   668)   ! derivative with respect to free-ion for equilibrium
   669)   ! surface complexation
   670)   ! 
   671)   ! Author: Glenn Hammond
   672)   ! Date: 10/22/08; 05/26/09; 03/16/12
   673)   ! 
   674) 
   675)   use Option_module
   676)   use Matrix_Block_Aux_module
   677)   
   678)   implicit none
   679)   
   680)   type(reactive_transport_auxvar_type) :: rt_auxvar
   681)   type(global_auxvar_type) :: global_auxvar
   682)   class(material_auxvar_type) :: material_auxvar
   683)   type(reaction_type) :: reaction
   684)   type(option_type) :: option
   685)   PetscInt :: irxn
   686)   PetscReal :: external_free_site_conc
   687)   PetscReal, pointer :: external_srfcplx_conc(:)
   688)   PetscReal :: external_total_sorb(reaction%naqcomp)
   689)   PetscReal :: external_dtotal_sorb(reaction%naqcomp,reaction%naqcomp)
   690)   PetscReal, pointer :: external_total_mob(:)
   691)   type(matrix_block_auxvar_type), pointer :: external_dRj_dCj
   692)   
   693)   PetscInt :: i, j, k, icplx, icomp, jcomp, ncomp, ncplx
   694)   PetscReal :: srfcplx_conc(reaction%surface_complexation%nsrfcplx)
   695)   PetscReal :: dSx_dmi(reaction%naqcomp)
   696)   PetscReal :: nui_Si_over_Sx
   697)   PetscReal :: ln_free_site
   698)   PetscReal :: lnQK, tempreal, tempreal1, tempreal2, total
   699)   PetscInt, parameter :: iphase = 1
   700)   PetscReal, parameter :: tol = 1.d-12
   701)   PetscReal :: ln_conc(reaction%naqcomp)
   702)   PetscReal :: ln_act(reaction%naqcomp)
   703)   PetscBool :: one_more
   704)   PetscReal :: res, dres_dfree_site, dfree_site_conc
   705)   PetscReal :: free_site_conc, rel_change_in_free_site_conc
   706)   PetscReal :: site_density(2)
   707)   PetscReal :: mobile_fraction
   708)   PetscInt :: num_types_of_sites
   709)   PetscInt :: isite
   710)   PetscInt :: num_iterations
   711)   PetscReal :: damping_factor
   712)   type(surface_complexation_type), pointer :: surface_complexation
   713) 
   714)   surface_complexation => reaction%surface_complexation
   715)   
   716)   ln_conc = log(rt_auxvar%pri_molal)
   717)   ln_act = ln_conc+log(rt_auxvar%pri_act_coef)
   718)   
   719)   ncplx = surface_complexation%srfcplxrxn_to_complex(0,irxn)
   720)     
   721)   free_site_conc = external_free_site_conc
   722)   srfcplx_conc = 0.d0
   723) 
   724)   select case(surface_complexation%srfcplxrxn_surf_type(irxn))
   725)     case(MINERAL_SURFACE)
   726)       site_density(1) = surface_complexation%srfcplxrxn_site_density(irxn)* &
   727)                 rt_auxvar%mnrl_volfrac(surface_complexation% &
   728)                                          srfcplxrxn_to_surf(irxn))
   729)       num_types_of_sites = 1
   730)     case(ROCK_SURFACE)
   731)       site_density(1) = surface_complexation%srfcplxrxn_site_density(irxn)* &
   732)                         material_auxvar%soil_particle_density * &
   733)                         (1.d0-material_auxvar%porosity)
   734)       num_types_of_sites = 1
   735)     case(COLLOID_SURFACE)
   736)       mobile_fraction = reaction%colloid_mobile_fraction( &
   737)                           surface_complexation%srfcplxrxn_to_surf(irxn))
   738)       site_density(1) = (1.d0-mobile_fraction)* &
   739)                         surface_complexation%srfcplxrxn_site_density(irxn)
   740)       site_density(2) = mobile_fraction* &
   741)                         surface_complexation%srfcplxrxn_site_density(irxn)
   742)       num_types_of_sites = 2 ! two types of sites (mobile and immobile) with
   743)                              ! separate site densities
   744)     case(NULL_SURFACE)
   745)       site_density(1) = surface_complexation%srfcplxrxn_site_density(irxn)
   746)       num_types_of_sites = 1
   747)   end select
   748)     
   749)   do isite=1, num_types_of_sites
   750)     ! isite == 1 - immobile (colloids, minerals, etc.)
   751)     ! isite == 2 - mobile (colloids)
   752)     
   753)     if (site_density(isite) < 1.d-40) cycle
   754)     
   755)     ! get a pointer to the first complex (there will always be at least 1)
   756)     ! in order to grab free site conc
   757)     one_more = PETSC_FALSE
   758)     num_iterations = 0
   759)     damping_factor = 1.d0
   760)     do
   761) 
   762)       num_iterations = num_iterations + 1
   763)       total = free_site_conc
   764)       ln_free_site = log(free_site_conc)
   765)       do j = 1, ncplx
   766)         icplx = surface_complexation%srfcplxrxn_to_complex(j,irxn)
   767)         ! compute secondary species concentration
   768)         lnQK = -surface_complexation%srfcplx_logK(icplx)*LOG_TO_LN
   769) 
   770)         ! activity of water
   771)         if (surface_complexation%srfcplxh2oid(icplx) > 0) then
   772)           lnQK = lnQK + surface_complexation%srfcplxh2ostoich(icplx)* &
   773)                         rt_auxvar%ln_act_h2o
   774)         endif
   775) 
   776)         lnQK = lnQK + surface_complexation%srfcplx_free_site_stoich(icplx)* &
   777)                       ln_free_site
   778)         
   779)         ncomp = surface_complexation%srfcplxspecid(0,icplx)
   780)         do i = 1, ncomp
   781)           icomp = surface_complexation%srfcplxspecid(i,icplx)
   782)           lnQK = lnQK + surface_complexation%srfcplxstoich(i,icplx)* &
   783)                         ln_act(icomp)
   784)         enddo
   785)         srfcplx_conc(icplx) = exp(lnQK)
   786)         total = total + surface_complexation%srfcplx_free_site_stoich(icplx)* &
   787)                         srfcplx_conc(icplx) 
   788)           
   789)       enddo
   790)         
   791)       if (one_more) exit
   792)         
   793)       if (surface_complexation%srfcplxrxn_stoich_flag(irxn)) then 
   794)         ! stoichiometry for free sites in one of reactions is not 1, thus must
   795)         ! use nonlinear iteration to solve
   796)         res = site_density(isite)-total
   797)           
   798)         dres_dfree_site = 1.d0
   799) 
   800)         do j = 1, ncplx
   801)           icplx = surface_complexation%srfcplxrxn_to_complex(j,irxn)
   802)           dres_dfree_site = dres_dfree_site + &
   803)             surface_complexation%srfcplx_free_site_stoich(icplx)* &
   804)             srfcplx_conc(icplx)/free_site_conc
   805)         enddo
   806) 
   807)         dfree_site_conc = res / dres_dfree_site
   808)         
   809)         ! if number of Newton iterations is excessive, try damping the update
   810)         if (num_iterations > 1000) then
   811)           damping_factor = 0.5d0
   812)         endif
   813)         free_site_conc = free_site_conc + damping_factor*dfree_site_conc
   814)         rel_change_in_free_site_conc = dabs(dfree_site_conc/free_site_conc)
   815)         
   816)         if (rel_change_in_free_site_conc < tol) then
   817)           one_more = PETSC_TRUE
   818)         endif
   819)         
   820)       else
   821)         
   822)         total = total / free_site_conc
   823)         free_site_conc = site_density(isite) / total  
   824)           
   825)         one_more = PETSC_TRUE 
   826)         
   827)       endif
   828) 
   829)     enddo ! generic do
   830)       
   831)     external_free_site_conc = free_site_conc
   832)    
   833) !!!!!!!!!!!!
   834)     ! 2.3-46
   835) 
   836)     ! Sx = free site
   837)     ! mi = molality of component i
   838)     dSx_dmi = 0.d0
   839)     tempreal = 0.d0
   840)     do j = 1, ncplx
   841)       icplx = surface_complexation%srfcplxrxn_to_complex(j,irxn)
   842)       ncomp = surface_complexation%srfcplxspecid(0,icplx)
   843)       do i = 1, ncomp
   844)         icomp = surface_complexation%srfcplxspecid(i,icplx)
   845)         ! sum of nu_li * nu_i * S_i
   846)         dSx_dmi(icomp) = dSx_dmi(icomp) + &
   847)           surface_complexation%srfcplxstoich(i,icplx)* &
   848)           surface_complexation%srfcplx_free_site_stoich(icplx)* &
   849)           srfcplx_conc(icplx)
   850)       enddo
   851)       ! sum of nu_i^2 * S_i
   852)       tempreal = tempreal + &
   853)         surface_complexation%srfcplx_free_site_stoich(icplx)* & 
   854)         surface_complexation%srfcplx_free_site_stoich(icplx)* &
   855)         srfcplx_conc(icplx)
   856)     enddo 
   857)     ! divide denominator by Sx
   858)     tempreal = tempreal / free_site_conc
   859)     ! add 1.d0 to denominator
   860)     tempreal = tempreal + 1.d0
   861)     ! divide numerator by denominator
   862)     dSx_dmi = -dSx_dmi / tempreal
   863)     ! convert from dlogm to dm
   864)     dSx_dmi = dSx_dmi / rt_auxvar%pri_molal
   865) !!!!!!!!!!!!
   866)     
   867)     if (isite == 1 .and. associated(external_srfcplx_conc)) then
   868)       external_srfcplx_conc(:) = external_srfcplx_conc(:) + srfcplx_conc(:)
   869)     endif
   870) #if 0
   871) !geh: for use if we decide to map surface complexes to master
   872)     if (isite == 1 .and. associated(ext_srfcplx_to_rxn_srfcplx_map)) then
   873)       do i = 1, surface_complexation%nsrfcplx
   874)         icplx = ext_srfcplx_to_rxn_srfcplx_map(i)
   875)         if (icplx > 0) then
   876)           external_srfcplx_conc(icplx) = external_srfcplx_conc(icplx) + &
   877)                                          srfcplx_conc(i)
   878)         endif
   879)       enddo
   880)     endif    
   881) #endif
   882)    
   883)     do k = 1, ncplx
   884)       icplx = surface_complexation%srfcplxrxn_to_complex(k,irxn)
   885) 
   886)       ncomp = surface_complexation%srfcplxspecid(0,icplx)
   887)       if (isite == 1) then ! immobile sites  
   888)         do i = 1, ncomp
   889)           icomp = surface_complexation%srfcplxspecid(i,icplx)
   890)           external_total_sorb(icomp) = external_total_sorb(icomp) + &
   891)             surface_complexation%srfcplxstoich(i,icplx)*srfcplx_conc(icplx)
   892)         enddo
   893)       else ! mobile sites
   894)         do i = 1, ncomp
   895)           icomp = reaction%pri_spec_to_coll_spec(surface_complexation% &
   896)                                                  srfcplxspecid(i,icplx))
   897)           external_total_mob(icomp) = external_total_mob(icomp) + &
   898)             surface_complexation%srfcplxstoich(i,icplx)*srfcplx_conc(icplx)
   899)         enddo
   900)       endif
   901)         
   902)       ! for 2.3-47 which feeds into 2.3-50
   903)       nui_Si_over_Sx = surface_complexation%srfcplx_free_site_stoich(icplx)* &
   904)                         srfcplx_conc(icplx)/ &
   905)                         free_site_conc
   906) 
   907)       do j = 1, ncomp
   908)         jcomp = surface_complexation%srfcplxspecid(j,icplx)
   909)         tempreal = surface_complexation%srfcplxstoich(j,icplx)* &
   910)                    srfcplx_conc(icplx) / &
   911)                    rt_auxvar%pri_molal(jcomp)+ &
   912)                    nui_Si_over_Sx*dSx_dmi(jcomp)
   913)         if (isite == 1) then ! immobile sites                  
   914)           do i = 1, ncomp
   915)             icomp = surface_complexation%srfcplxspecid(i,icplx)
   916)             external_dtotal_sorb(icomp,jcomp) = &
   917)                                 external_dtotal_sorb(icomp,jcomp) + &
   918)                                 surface_complexation%srfcplxstoich(i,icplx)* &
   919)                                 tempreal
   920)           enddo ! i
   921)         else ! mobile sites
   922)           do i = 1, ncomp
   923)             icomp = surface_complexation%srfcplxspecid(i,icplx)
   924)             external_dRj_dCj%dtotal(icomp,jcomp,1) = &
   925)                                  external_dRj_dCj%dtotal(icomp,jcomp,1) + &
   926)                                  surface_complexation%srfcplxstoich(i,icplx)* &
   927)                                  tempreal
   928)           enddo ! i
   929)         endif
   930)       enddo ! j
   931)     enddo ! k
   932)   enddo ! isite
   933)   
   934) end subroutine RTotalSorbEqSurfCplx1
   935) 
   936) ! ************************************************************************** !
   937) 
   938) subroutine RKineticSurfCplx(Res,Jac,compute_derivative,rt_auxvar, &
   939)                             global_auxvar,material_auxvar,reaction,option)
   940)   ! 
   941)   ! Computes contribution to residual and jacobian for
   942)   ! kinetic surface complexation reactions
   943)   ! 
   944)   ! Author: Glenn Hammond
   945)   ! Date: 12/07/09
   946)   ! 
   947) 
   948)   use Option_module
   949)   use Material_Aux_class
   950)   
   951)   implicit none
   952)   
   953)   PetscBool :: compute_derivative
   954)   type(reactive_transport_auxvar_type) :: rt_auxvar
   955)   type(global_auxvar_type) :: global_auxvar
   956)   class(material_auxvar_type) :: material_auxvar
   957)   type(reaction_type) :: reaction
   958)   PetscReal :: Res(reaction%ncomp)
   959)   PetscReal :: Jac(reaction%ncomp,reaction%ncomp)
   960)   type(option_type) :: option
   961) 
   962)   PetscInt :: i, j, k, l, icplx, icomp, jcomp, lcomp, ncomp, ncplx
   963)   PetscReal :: ln_conc(reaction%naqcomp)
   964)   PetscReal :: ln_act(reaction%naqcomp)
   965)   PetscInt :: irxn, isite, ikinrxn
   966)   PetscReal :: dt
   967)   type(surface_complexation_type), pointer :: surface_complexation
   968) 
   969)   PetscReal :: numerator_sum(reaction%surface_complexation%nkinsrfcplxrxn)
   970)   PetscReal :: denominator_sum(reaction%surface_complexation%nkinsrfcplxrxn)
   971) 
   972)   PetscReal :: denominator
   973)   PetscReal :: fac
   974)   PetscReal :: fac_sum(reaction%naqcomp)
   975)   PetscReal :: lnQ(reaction%surface_complexation%nkinsrfcplx)
   976)   PetscReal :: Q(reaction%surface_complexation%nkinsrfcplx)
   977)   PetscReal :: srfcplx_conc_k(maxval(reaction%surface_complexation%srfcplxrxn_to_complex(0,:)),1)
   978)   PetscReal :: srfcplx_conc_kp1(maxval(reaction%surface_complexation%srfcplxrxn_to_complex(0,:)),1)
   979)   
   980)   surface_complexation => reaction%surface_complexation
   981)   
   982)   ln_conc = log(rt_auxvar%pri_molal)
   983)   ln_act = ln_conc+log(rt_auxvar%pri_act_coef)
   984) 
   985) ! Members of the rt aux var object: mol/m^3
   986) ! PetscReal, pointer :: kinsrfcplx_conc(:)          ! S_{i\alpha}^k
   987) ! PetscReal, pointer :: kinsrfcplx_conc_kp1(:)      ! S_{i\alpha}^k+1
   988) ! PetscReal, pointer :: kinsrfcplx_free_site_conc(:) ! S_\alpha
   989)   
   990) ! units
   991) ! k_f: dm^3/mol/sec
   992) ! k_b: 1/sec
   993) ! Res: mol/sec
   994) 
   995)   dt = option%tran_dt
   996)   
   997) ! compute ion activity product and store: units mol/L
   998)   lnQ = 0.d0
   999)   do ikinrxn = 1, surface_complexation%nkinsrfcplxrxn
  1000)     irxn = surface_complexation%kinsrfcplxrxn_to_srfcplxrxn(ikinrxn)
  1001)     ncplx = surface_complexation%srfcplxrxn_to_complex(0,irxn)
  1002)     do k = 1, ncplx ! ncplx in rxn
  1003)       icplx = surface_complexation%srfcplxrxn_to_complex(k,irxn)
  1004)       if (surface_complexation%srfcplxh2oid(icplx) > 0) then
  1005)         lnQ(icplx) = lnQ(icplx) + surface_complexation%srfcplxh2ostoich(icplx)* &
  1006)           rt_auxvar%ln_act_h2o
  1007)       endif
  1008)     
  1009)       ncomp = surface_complexation%srfcplxspecid(0,icplx)
  1010)       do i = 1, ncomp
  1011)         icomp = surface_complexation%srfcplxspecid(i,icplx)
  1012)         lnQ(icplx) = lnQ(icplx) + surface_complexation%srfcplxstoich(i,icplx)* &
  1013)           ln_act(icomp)
  1014)       enddo
  1015)       Q(icplx) = exp(lnQ(icplx))
  1016)     enddo
  1017)   enddo
  1018)     
  1019)   ! compute summation in numerator of 5.1-29: units mol/m^3
  1020)   numerator_sum = 0.d0
  1021)   do ikinrxn = 1, surface_complexation%nkinsrfcplxrxn
  1022)     irxn = surface_complexation%kinsrfcplxrxn_to_srfcplxrxn(ikinrxn)
  1023)     isite = surface_complexation%srfcplxrxn_to_surf(irxn)
  1024)     ncplx = surface_complexation%srfcplxrxn_to_complex(0,irxn)
  1025)     do k = 1, ncplx ! ncplx in rxn
  1026)       icplx = surface_complexation%srfcplxrxn_to_complex(k,irxn)
  1027)       numerator_sum(isite) = numerator_sum(isite) + &
  1028)                   rt_auxvar%kinsrfcplx_conc(icplx,ikinrxn)/ &
  1029)                   (1.d0+surface_complexation%kinsrfcplx_backward_rate(icplx,ikinrxn)*dt)
  1030)     enddo
  1031)   enddo
  1032) 
  1033)   do ikinrxn = 1, surface_complexation%nkinsrfcplxrxn
  1034)     irxn = surface_complexation%kinsrfcplxrxn_to_srfcplxrxn(ikinrxn)
  1035)     isite = surface_complexation%srfcplxrxn_to_surf(irxn)
  1036)     numerator_sum(isite) = surface_complexation%srfcplxrxn_site_density(isite) - &
  1037)                            numerator_sum(isite)
  1038)   enddo
  1039)   
  1040)   ! compute summation in denominator of 5.1-29
  1041)   denominator_sum = 1.d0
  1042)   do ikinrxn = 1, surface_complexation%nkinsrfcplxrxn
  1043)     irxn = surface_complexation%kinsrfcplxrxn_to_srfcplxrxn(ikinrxn)
  1044)     isite = surface_complexation%srfcplxrxn_to_surf(irxn)
  1045)     ncplx = surface_complexation%srfcplxrxn_to_complex(0,irxn)
  1046)     do k = 1, ncplx ! ncplx in rxn
  1047)       icplx = surface_complexation%srfcplxrxn_to_complex(k,irxn)
  1048)       denominator_sum(isite) = denominator_sum(isite) + &
  1049)                                (surface_complexation%kinsrfcplx_forward_rate(icplx,ikinrxn)*dt)/ &
  1050)                                (1.d0+surface_complexation%kinsrfcplx_backward_rate(icplx,ikinrxn)*dt)* &
  1051)                                Q(icplx)
  1052)     enddo
  1053)   enddo
  1054) 
  1055) ! compute surface complex conc. at new time step (5.1-30)  
  1056)   do ikinrxn = 1, surface_complexation%nkinsrfcplxrxn
  1057)     irxn = surface_complexation%kinsrfcplxrxn_to_srfcplxrxn(ikinrxn)
  1058)     isite = surface_complexation%srfcplxrxn_to_surf(irxn)
  1059)     ncplx = surface_complexation%srfcplxrxn_to_complex(0,irxn)
  1060)     do k = 1, ncplx ! ncplx in rxn
  1061)       icplx = surface_complexation%srfcplxrxn_to_complex(k,irxn)
  1062)       srfcplx_conc_k(icplx,ikinrxn) = rt_auxvar%kinsrfcplx_conc(icplx,ikinrxn)
  1063)       denominator = 1.d0 + surface_complexation%kinsrfcplx_backward_rate(icplx,ikinrxn)*dt
  1064)       srfcplx_conc_kp1(icplx,ikinrxn) = (srfcplx_conc_k(icplx,ikinrxn) + &
  1065)                                 surface_complexation%kinsrfcplx_forward_rate(icplx,ikinrxn)*dt * &
  1066)                                 numerator_sum(isite)/denominator_sum(isite)* &
  1067)                                 Q(icplx))/denominator
  1068)       rt_auxvar%kinsrfcplx_conc_kp1(icplx,ikinrxn) = srfcplx_conc_kp1(icplx,ikinrxn)
  1069)     enddo
  1070)     rt_auxvar%kinsrfcplx_free_site_conc(isite) = numerator_sum(isite)/ &
  1071)                                                denominator_sum(isite)
  1072)   enddo
  1073) 
  1074) ! compute residual (5.1-34)
  1075)   
  1076)   do ikinrxn = 1, surface_complexation%nkinsrfcplxrxn
  1077)     irxn = surface_complexation%kinsrfcplxrxn_to_srfcplxrxn(ikinrxn)
  1078)     ncplx = surface_complexation%srfcplxrxn_to_complex(0,irxn)
  1079)     do k = 1, ncplx ! ncplx in rxn
  1080)       icplx = surface_complexation%srfcplxrxn_to_complex(k,irxn)
  1081)       ncomp = surface_complexation%srfcplxspecid(0,icplx)
  1082)       do i = 1, ncomp
  1083)         icomp = surface_complexation%srfcplxspecid(i,icplx)
  1084)         Res(icomp) = Res(icomp) + surface_complexation%srfcplxstoich(i,icplx)* &
  1085)                      (srfcplx_conc_kp1(icplx,ikinrxn)-rt_auxvar%kinsrfcplx_conc(icplx,ikinrxn))/ &
  1086)                      dt * material_auxvar%volume
  1087)       enddo
  1088)     enddo
  1089)   enddo
  1090) 
  1091)   if (compute_derivative) then
  1092)   ! compute jacobian (5.1-39)
  1093)     fac_sum = 0.d0
  1094)     do ikinrxn = 1, surface_complexation%nkinsrfcplxrxn
  1095)       irxn = surface_complexation%kinsrfcplxrxn_to_srfcplxrxn(ikinrxn)
  1096)       ncplx = surface_complexation%srfcplxrxn_to_complex(0,irxn)
  1097)       do k = 1, ncplx ! ncplx in rxn
  1098)         icplx = surface_complexation%srfcplxrxn_to_complex(k,irxn)
  1099)         denominator = 1.d0 + surface_complexation%kinsrfcplx_backward_rate(icplx,ikinrxn)*dt
  1100)         fac = surface_complexation%kinsrfcplx_forward_rate(icplx,ikinrxn)/denominator
  1101)         ncomp = surface_complexation%srfcplxspecid(0,icplx)
  1102)         do j = 1, ncomp
  1103)           jcomp = surface_complexation%srfcplxspecid(j,icplx)
  1104)           fac_sum(jcomp) = fac_sum(jcomp) + surface_complexation%srfcplxstoich(j,icplx)* &
  1105)             fac * Q(icplx)
  1106)         enddo
  1107)       enddo
  1108)     enddo
  1109) 
  1110)     do ikinrxn = 1, surface_complexation%nkinsrfcplxrxn
  1111)       irxn = surface_complexation%kinsrfcplxrxn_to_srfcplxrxn(ikinrxn)
  1112)       isite = surface_complexation%srfcplxrxn_to_surf(irxn)
  1113)       ncplx = surface_complexation%srfcplxrxn_to_complex(0,irxn)
  1114)       do k = 1, ncplx ! ncplx in rxn
  1115)         icplx = surface_complexation%srfcplxrxn_to_complex(k,irxn)
  1116)         denominator = 1.d0 + surface_complexation%kinsrfcplx_backward_rate(icplx,ikinrxn)*dt
  1117)         fac = surface_complexation%kinsrfcplx_forward_rate(icplx,ikinrxn)/denominator
  1118)         ncomp = surface_complexation%srfcplxspecid(0,icplx)
  1119)         do j = 1, ncomp
  1120)           jcomp = surface_complexation%srfcplxspecid(j,icplx)
  1121)           do l = 1, ncomp
  1122)             lcomp = surface_complexation%srfcplxspecid(l,icplx)
  1123)             Jac(jcomp,lcomp) = Jac(jcomp,lcomp) + &
  1124)               (surface_complexation%srfcplxstoich(j,icplx) * fac * numerator_sum(isite) * &
  1125)               Q(icplx) * (surface_complexation%srfcplxstoich(l,icplx) - &
  1126)               dt * fac_sum(lcomp)/denominator_sum(isite)))/denominator_sum(isite) * &
  1127)               exp(-ln_conc(lcomp)) * material_auxvar%volume
  1128)           enddo
  1129)         enddo
  1130)       enddo
  1131)     enddo
  1132)   endif
  1133)   
  1134)   ! units of total_sorb = mol/m^3
  1135)   ! units of dtotal_sorb = kg water/m^3 bulk
  1136)   
  1137) end subroutine RKineticSurfCplx
  1138) 
  1139) end module Reaction_Surface_Complexation_module

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