pm_ufd_decay.F90       coverage:  52.17 %func     79.72 %block


     1) module PM_UFD_Decay_class
     2) 
     3)   use PM_Base_class
     4)   use Realization_Subsurface_class
     5)   use Option_module
     6)   
     7)   use PFLOTRAN_Constants_module
     8) 
     9)   implicit none
    10) 
    11)   private
    12) 
    13) #include "petsc/finclude/petscsys.h"
    14) 
    15)   
    16)   type, public, extends(pm_base_type) :: pm_ufd_decay_type
    17)     class(realization_subsurface_type), pointer :: realization
    18)     PetscInt, pointer :: element(:)
    19)     PetscInt, pointer :: element_isotopes(:,:)
    20)     PetscInt, pointer :: isotope_to_primary_species(:)
    21)     PetscInt, pointer :: isotope_to_mineral(:)
    22)     PetscReal, pointer :: isotope_decay_rate(:)
    23)     PetscInt, pointer :: isotope_daughters(:,:)
    24)     PetscInt, pointer :: isotope_parents(:,:)
    25)     PetscReal, pointer :: element_solubility(:)
    26)     PetscReal, pointer :: element_Kd(:,:)
    27)     PetscInt :: num_elements
    28)     PetscInt :: num_isotopes
    29)     PetscInt, pointer :: num_isotopes_per_element(:)
    30)     type(isotope_type), pointer :: isotope_list
    31)     type(element_type), pointer :: element_list
    32)   contains
    33) !geh: commented out subroutines can only be called externally
    34)     procedure, public :: Setup => PMUFDDecayInit
    35)     procedure, public :: Read => PMUFDDecayRead
    36) !    procedure, public :: SetupSolvers => PMUFDDecaySetupSolvers
    37)     procedure, public :: PMUFDDecaySetRealization
    38)     procedure, public :: InitializeRun => PMUFDDecayInitializeRun
    39) !!    procedure, public :: FinalizeRun => PMUFDDecayFinalizeRun
    40)     procedure, public :: InitializeTimestep => PMUFDDecayInitializeTimestep
    41)     procedure, public :: FinalizeTimestep => PMUFDDecayFinalizeTimestep
    42) !    procedure, public :: PreSolve => PMUFDDecayPreSolve
    43)     procedure, public :: Solve => PMUFDDecaySolve
    44) !    procedure, public :: PostSolve => PMUFDDecayPostSolve
    45) !    procedure, public :: AcceptSolution => PMUFDDecayAcceptSolution
    46) !    procedure, public :: TimeCut => PMUFDDecayTimeCut
    47) !    procedure, public :: UpdateSolution => PMUFDDecayUpdateSolution
    48) !    procedure, public :: UpdateAuxvars => PMUFDDecayUpdateAuxvars
    49) !    procedure, public :: Checkpoint => PMUFDDecayCheckpoint    
    50) !    procedure, public :: Restart => PMUFDDecayRestart  
    51)     procedure, public :: InputRecord => PMUFDDecayInputRecord
    52)     procedure, public :: Destroy => PMUFDDecayDestroy
    53)   end type pm_ufd_decay_type
    54)   
    55)   type :: isotope_type
    56)     character(len=MAXWORDLENGTH) :: name
    57)     character(len=MAXWORDLENGTH) :: element
    58)     PetscInt :: iisotope
    59)     PetscInt :: ielement
    60)     PetscReal :: decay_rate
    61)     type(daughter_type), pointer :: daughter_list
    62)     type(isotope_type), pointer :: next
    63)   end type isotope_type
    64)   
    65)   type :: daughter_type
    66)     character(len=MAXWORDLENGTH) :: name
    67)     PetscReal :: stoichiometry
    68)     type(daughter_type), pointer :: next
    69)   end type
    70)   
    71)   type :: element_type
    72)     character(len=MAXWORDLENGTH) :: name
    73)     PetscInt :: ielement
    74)     PetscReal :: solubility
    75)     PetscReal, pointer :: Kd(:)
    76)     character(len=MAXWORDLENGTH), pointer :: Kd_material_name(:)
    77)     type(element_type), pointer :: next
    78)   end type
    79)   
    80)   public :: PMUFDDecayCreate, &
    81)             PMUFDDecayInit !, &
    82) !            PMUFDDecaySetupSolvers, &
    83) !            PMUFDDecayInitializeTimestepA, &
    84) !            PMUFDDecayInitializeTimestepB, &
    85) !            PMUFDDecayInitializeRun, &
    86) !            PMUFDDecayUpdateSolution, &
    87) !            PMUFDDecayUpdatePropertiesNI, &
    88) !            PMUFDDecayTimeCut, &
    89) !            PMUFDDecayDestroy
    90)   
    91) contains
    92) 
    93) ! ************************************************************************** !
    94) 
    95) function PMUFDDecayCreate()
    96)   ! 
    97)   ! Creates the UFD decay process model
    98)   ! 
    99)   ! Author: Glenn Hammond
   100)   ! Date: 06/24/15
   101) 
   102)   implicit none
   103)   
   104)   class(pm_ufd_decay_type), pointer :: PMUFDDecayCreate
   105)   
   106)   allocate(PMUFDDecayCreate)
   107)   call PMBaseInit(PMUFDDecayCreate)
   108) 
   109)   nullify(PMUFDDecayCreate%realization)
   110)   nullify(PMUFDDecayCreate%element)
   111)   nullify(PMUFDDecayCreate%element_isotopes)
   112)   nullify(PMUFDDecayCreate%isotope_to_primary_species)
   113)   nullify(PMUFDDecayCreate%isotope_to_mineral)
   114)   nullify(PMUFDDecayCreate%isotope_decay_rate)
   115)   nullify(PMUFDDecayCreate%isotope_daughters)
   116)   nullify(PMUFDDecayCreate%isotope_parents)
   117)   nullify(PMUFDDecayCreate%element_solubility)
   118)   nullify(PMUFDDecayCreate%element_Kd)
   119)   nullify(PMUFDDecayCreate%num_isotopes_per_element)
   120)   nullify(PMUFDDecayCreate%isotope_list)
   121)   nullify(PMUFDDecayCreate%element_list)
   122) 
   123)   call PMBaseInit(PMUFDDecayCreate)
   124) 
   125) end function PMUFDDecayCreate
   126) 
   127) ! ************************************************************************** !
   128) 
   129) subroutine PMUFDDecayRead(this,input)
   130)   ! 
   131)   ! Reads input file parameters associated with the ufd decay process model
   132)   ! 
   133)   ! Author: Glenn Hammond
   134)   ! Date: 06/24/15
   135) 
   136)   use Input_Aux_module
   137)   use String_module
   138)   use Utility_module
   139)   use Option_module
   140)   
   141)   implicit none
   142)   
   143)   class(pm_ufd_decay_type) :: this
   144)   type(input_type), pointer :: input
   145)   
   146)   type(option_type), pointer :: option
   147)   character(len=MAXWORDLENGTH) :: word
   148)   character(len=MAXSTRINGLENGTH) :: error_string
   149)   type(isotope_type), pointer :: isotope, prev_isotope
   150)   type(daughter_type), pointer :: daughter, prev_daughter
   151)   type(element_type), pointer :: element, prev_element
   152)   PetscInt :: i
   153)   character(len=MAXWORDLENGTH) :: Kd_material_name(20)
   154)   PetscReal :: Kd(20)
   155) 
   156)   option => this%option
   157)   
   158)   option%io_buffer = 'pflotran card:: UFD_Decay'
   159)   call printMsg(option)
   160)   
   161)   input%ierr = 0
   162)   nullify(prev_isotope)
   163)   nullify(prev_element)
   164)   do
   165)   
   166)     call InputReadPflotranString(input,option)
   167)     if (InputError(input)) exit
   168)     if (InputCheckExit(input,option)) exit
   169)     
   170)     call InputReadWord(input,option,word,PETSC_TRUE)
   171)     call InputErrorMsg(input,option,'keyword',error_string)
   172)     call StringToUpper(word)
   173)     
   174)     select case(trim(word))
   175)       case('ELEMENT')
   176)         error_string = 'UFD Decay, Element'
   177)         element => ElementCreate()
   178)         call InputReadWord(input,option,element%name,PETSC_TRUE)
   179)         call InputErrorMsg(input,option,'name',error_string)
   180)         do
   181)           call InputReadPflotranString(input,option)
   182)           if (InputError(input)) exit
   183)           if (InputCheckExit(input,option)) exit
   184)           call InputReadWord(input,option,word,PETSC_TRUE)
   185)           call InputErrorMsg(input,option,'keyword',error_string)
   186)           call StringToUpper(word)
   187)           select case(trim(word))
   188)             case('SOLUBILITY')
   189)               call InputReadDouble(input,option,element%solubility)
   190)               call InputErrorMsg(input,option,'solubility',error_string)
   191)             case('KD')
   192)               i = 0
   193)               Kd(:) = UNINITIALIZED_DOUBLE
   194)               Kd_material_name(:) = ''
   195)               do
   196)                 call InputReadPflotranString(input,option)
   197)                 if (InputError(input)) exit
   198)                 if (InputCheckExit(input,option)) exit
   199)                 i = i + 1
   200)                 if (i > 20) then
   201)                   write(word,*) i-1
   202)                   option%io_buffer = 'Kd array in PMUFDDecayRead() must be &
   203)                     &allocated larger than ' // trim(adjustl(word)) // '.'
   204)                   call printErrMsg(option)
   205)                 endif
   206)                 call InputReadWord(input,option,word,PETSC_TRUE)
   207)                 call InputErrorMsg(input,option,'Kd material name', &
   208)                                    error_string)
   209)                 Kd_material_name(i) = word
   210)                 call InputReadDouble(input,option,Kd(i))
   211)                 call InputErrorMsg(input,option,'Kd',error_string)
   212)               enddo
   213)               if (i == 0) then
   214)                 option%io_buffer = 'No KD/material combinations specified &
   215)                   &under ' // trim(error_string) // '.'
   216)                 call printErrMsg(option)
   217)               endif
   218)               allocate(element%Kd(i))
   219)               element%Kd = Kd(1:i)
   220)               allocate(element%Kd_material_name(i))
   221)               element%Kd_material_name = Kd_material_name(1:i)
   222)             case default
   223)               call InputKeywordUnrecognized(word,error_string,option)
   224)           end select
   225)         enddo
   226)         if (associated(prev_element)) then
   227)           prev_element%next => element
   228)         else
   229)           this%element_list => element
   230)         endif
   231)         prev_element => element
   232)         nullify(element)
   233)       case('ISOTOPE')
   234)         nullify(prev_daughter)
   235)         error_string = 'UFD Decay, Isotope'
   236)         isotope => IsotopeCreate()
   237)         call InputReadWord(input,option,isotope%name,PETSC_TRUE)
   238)         call InputErrorMsg(input,option,'name',error_string)
   239)         error_string = 'UFD Decay, Isotope, ' // trim(isotope%name)
   240)         do
   241)           call InputReadPflotranString(input,option)
   242)           if (InputError(input)) exit
   243)           if (InputCheckExit(input,option)) exit
   244)           call InputReadWord(input,option,word,PETSC_TRUE)
   245)           call InputErrorMsg(input,option,'keyword',error_string)
   246)           call StringToUpper(word)
   247)           select case(trim(word))
   248)             case('ELEMENT')
   249)               call InputReadWord(input,option,isotope%element,PETSC_TRUE)
   250)               call InputErrorMsg(input,option,'element name',error_string)
   251)             case('DECAY_RATE')
   252)               call InputReadDouble(input,option,isotope%decay_rate)
   253)               call InputErrorMsg(input,option,'decay rate',error_string)
   254)             case('DAUGHTER')
   255)               daughter => IsotopeDaughterCreate()
   256)               call InputReadWord(input,option,daughter%name,PETSC_TRUE)
   257)               call InputErrorMsg(input,option,'daughter name',error_string)
   258)               call InputReadDouble(input,option,daughter%stoichiometry)
   259)               call InputErrorMsg(input,option,'daughter stoichiometry', &
   260)                                  error_string)
   261)               if (associated(prev_daughter)) then
   262)                 prev_daughter%next => daughter
   263)               else
   264)                 isotope%daughter_list => daughter
   265)               endif
   266)               prev_daughter => daughter
   267)               nullify(daughter)
   268)             case default
   269)               call InputKeywordUnrecognized(word,error_string,option)
   270)           end select
   271)         enddo
   272)         if (associated(prev_isotope)) then
   273)           prev_isotope%next => isotope
   274)         else
   275)           this%isotope_list => isotope
   276)         endif
   277)         prev_isotope => isotope
   278)         nullify(isotope)
   279)       case default
   280)         call InputKeywordUnrecognized(word,error_string,option)
   281)     end select
   282)   enddo
   283)   
   284) end subroutine PMUFDDecayRead
   285) 
   286) 
   287) ! ************************************************************************** !
   288) 
   289) function ElementCreate()
   290)   ! 
   291)   ! Creates isotope object
   292)   ! 
   293)   ! Author: Glenn Hammond
   294)   ! Date: 11/20/15
   295)   ! 
   296)   implicit none
   297) 
   298)   type(element_type), pointer :: ElementCreate
   299)   
   300)   allocate(ElementCreate)
   301)   ElementCreate%name = ''
   302)   ElementCreate%ielement = UNINITIALIZED_INTEGER
   303)   ElementCreate%solubility = UNINITIALIZED_DOUBLE
   304)   nullify(ElementCreate%Kd)
   305)   nullify(ElementCreate%Kd_material_name)
   306)   nullify(ElementCreate%next)
   307)   
   308) end function ElementCreate
   309) 
   310) ! ************************************************************************** !
   311) 
   312) function IsotopeCreate()
   313)   ! 
   314)   ! Creates isotope object
   315)   ! 
   316)   ! Author: Glenn Hammond
   317)   ! Date: 11/20/15
   318)   ! 
   319)   implicit none
   320) 
   321)   type(isotope_type), pointer :: IsotopeCreate
   322)   
   323)   allocate(IsotopeCreate)
   324)   IsotopeCreate%name = ''
   325)   IsotopeCreate%element = ''
   326)   IsotopeCreate%iisotope = UNINITIALIZED_INTEGER
   327)   IsotopeCreate%ielement = UNINITIALIZED_INTEGER
   328)   IsotopeCreate%decay_rate = UNINITIALIZED_DOUBLE
   329)   nullify(IsotopeCreate%daughter_list)
   330)   nullify(IsotopeCreate%next)
   331)   
   332) end function IsotopeCreate
   333) 
   334) ! ************************************************************************** !
   335) 
   336) function IsotopeDaughterCreate()
   337)   ! 
   338)   ! Creates daughter object
   339)   ! 
   340)   ! Author: Glenn Hammond
   341)   ! Date: 11/20/15
   342)   ! 
   343)   implicit none
   344) 
   345)   type(daughter_type), pointer :: IsotopeDaughterCreate
   346)   
   347)   allocate(IsotopeDaughterCreate)
   348)   IsotopeDaughterCreate%name = ''
   349)   IsotopeDaughterCreate%stoichiometry = UNINITIALIZED_DOUBLE
   350)   nullify(IsotopeDaughterCreate%next)
   351) 
   352) end function IsotopeDaughterCreate
   353) 
   354) ! ************************************************************************** !
   355) 
   356) subroutine PMUFDDecayInit(this)
   357)   ! 
   358)   ! Initializes variables associated with the UFD decay process model
   359)   ! 
   360)   ! Author: Glenn Hammond
   361)   ! Date: 06/24/15
   362) 
   363)   use Option_module
   364)   use String_module
   365)   use Grid_module
   366)   use Reaction_Aux_module
   367)   use Reaction_Mineral_Aux_module
   368)   use Reactive_Transport_Aux_module
   369)   use Material_module
   370)   
   371)   implicit none
   372)   
   373)   class(pm_ufd_decay_type) :: this
   374)   
   375)   type(option_type), pointer :: option
   376)   type(reaction_type), pointer :: reaction
   377)   type(reactive_transport_auxvar_type), pointer :: rt_auxvars(:)
   378)   type(grid_type), pointer :: grid
   379)   type(isotope_type), pointer :: isotope, isotope2
   380)   type(daughter_type), pointer :: daughter
   381)   type(element_type), pointer :: element
   382)   PetscInt, allocatable :: num_isotopes_per_element(:)
   383)   character(len=MAXWORDLENGTH) :: word
   384)   type(material_property_ptr_type), pointer :: material_property_array(:)
   385)   type(material_property_type), pointer :: material_property
   386)   
   387)   PetscInt :: icount
   388)   PetscInt :: ghosted_id
   389)   PetscInt :: max_daughters_per_isotope
   390)   PetscInt :: max_parents_per_isotope
   391)   PetscBool :: found
   392)   PetscInt :: iisotope
   393)   
   394)   option => this%realization%option
   395)   grid => this%realization%patch%grid
   396)   reaction => this%realization%reaction
   397)   rt_auxvars => this%realization%patch%aux%RT%auxvars
   398)   material_property_array => this%realization%patch%material_property_array
   399)   
   400)   do ghosted_id = 1, grid%ngmax
   401)     allocate(rt_auxvars(ghosted_id)%total_sorb_eq(reaction%naqcomp))
   402)     rt_auxvars(ghosted_id)%total_sorb_eq = 0.d0
   403)   enddo
   404)   
   405)   max_daughters_per_isotope = 0
   406)   max_parents_per_isotope = 0
   407) 
   408)   ! sum the number of isotopes, elements, max number of isotopes per element
   409)   
   410)   element => this%element_list
   411)   icount = 0
   412)   do
   413)     if (.not.associated(element)) exit
   414)     icount = icount + 1
   415)     element%ielement = icount
   416)     element => element%next
   417)   enddo
   418)   this%num_elements = icount
   419)   allocate(num_isotopes_per_element(this%num_elements))
   420)   num_isotopes_per_element = 0
   421)   allocate(this%element_solubility(this%num_elements))
   422)   this%element_solubility = 0.d0
   423)   allocate(this%element_Kd(this%num_elements,size(material_property_array)))
   424)   this%element_Kd = UNINITIALIZED_DOUBLE
   425)   element => this%element_list
   426)   do
   427)     if (.not.associated(element)) exit
   428)     this%element_solubility(element%ielement) = element%solubility
   429)     if (size(element%Kd) /= size(material_property_array)) then
   430)       write(word,*) size(element%Kd)
   431)       option%io_buffer = 'Incorrect number of Kds (' // &
   432)         trim(adjustl(word)) // ') specified for number of materials ('
   433)       write(word,*) size(material_property_array)
   434)       option%io_buffer = trim(option%io_buffer) // &
   435)         trim(adjustl(word)) // ') for UFD Decay element "' // &
   436)         trim(element%name) // '".'
   437)     endif
   438)     do icount = 1, size(element%Kd_material_name)
   439)       material_property => &
   440)         MaterialPropGetPtrFromArray(element%Kd_material_name(icount), &
   441)                                     material_property_array)
   442)       this%element_Kd(element%ielement,material_property%internal_id) = &
   443)         element%Kd(icount)
   444)     enddo
   445)     do icount = 1, size(material_property_array)
   446)       if (UnInitialized(this%element_Kd(element%ielement,icount))) then
   447)         option%io_buffer = 'Uninitialized KD in UFD Decay element "' // &
   448)           trim(element%name) // '" for material "' // &
   449)           trim(material_property_array(icount)%ptr%name) // '".'
   450)         call printErrMsg(option)
   451)       endif
   452)     enddo
   453)     element => element%next
   454)   enddo  
   455)   
   456)   this%num_isotopes = 0
   457)   isotope => this%isotope_list
   458)   do
   459)     if (.not.associated(isotope)) exit
   460)     this%num_isotopes = this%num_isotopes + 1
   461)     isotope%iisotope = this%num_isotopes
   462)     found = PETSC_FALSE
   463)     element => this%element_list
   464)     do
   465)       if (.not.associated(element)) exit
   466)       if (StringCompare(isotope%element,element%name)) then
   467)         found = PETSC_TRUE
   468)         isotope%ielement = element%ielement
   469)         num_isotopes_per_element(element%ielement) = &
   470)           num_isotopes_per_element(element%ielement) + 1
   471)         exit
   472)       endif
   473)       element => element%next
   474)     enddo
   475)     if (.not.found) then
   476)       option%io_buffer = 'Element "' // trim(isotope%element) // &
   477)         '" of isotope "' // trim(isotope%name) // &
   478)         '" not found among list of elements.'
   479)       call printErrMsg(option)
   480)     endif
   481)     daughter => isotope%daughter_list
   482)     icount = 0
   483)     do
   484)       if (.not.associated(daughter)) exit
   485)       icount = icount + 1
   486)       daughter => daughter%next
   487)     enddo
   488)     max_daughters_per_isotope = max(max_daughters_per_isotope,icount)
   489)     isotope => isotope%next
   490)   enddo
   491)   
   492)   allocate(this%isotope_to_primary_species(this%num_isotopes))
   493)   this%isotope_to_primary_species = UNINITIALIZED_INTEGER
   494)   allocate(this%isotope_to_mineral(this%num_isotopes))
   495)   this%isotope_to_mineral = UNINITIALIZED_INTEGER
   496)   allocate(this%element_isotopes(0:maxval(num_isotopes_per_element), &
   497)            this%num_elements))
   498)   this%element_isotopes = UNINITIALIZED_INTEGER
   499)   this%element_isotopes(0,:) = 0
   500)   allocate(this%isotope_decay_rate(this%num_isotopes))
   501)   this%isotope_decay_rate = UNINITIALIZED_DOUBLE
   502)   allocate(this%isotope_daughters(0:max_daughters_per_isotope,this%num_isotopes))
   503)   this%isotope_daughters = UNINITIALIZED_INTEGER
   504)   this%isotope_daughters(0,:) = 0
   505)   
   506)   isotope => this%isotope_list
   507)   do
   508)     if (.not.associated(isotope)) exit
   509)     found = PETSC_FALSE
   510)     this%isotope_to_primary_species(isotope%iisotope) = &
   511)       GetPrimarySpeciesIDFromName(isotope%name,reaction,option)
   512)     word = isotope%name
   513)     word = trim(word) // '(s)'
   514)     this%isotope_to_mineral(isotope%iisotope) = &
   515)       GetKineticMineralIDFromName(word,reaction%mineral,option)
   516)     this%element_isotopes(0,isotope%ielement) = &
   517)       this%element_isotopes(0,isotope%ielement) + 1
   518)     this%element_isotopes(this%element_isotopes(0,isotope%ielement), &
   519)                           isotope%ielement) = isotope%iisotope
   520)     this%isotope_decay_rate(isotope%iisotope) = isotope%decay_rate
   521)     daughter => isotope%daughter_list
   522)     icount = 0
   523)     do
   524)       if (.not.associated(daughter)) exit
   525)       icount = icount + 1
   526)       isotope2 => this%isotope_list
   527)       found = PETSC_FALSE
   528)       do
   529)         if (.not.associated(isotope2)) exit
   530)         if (StringCompare(daughter%name,isotope2%name)) then
   531)           found = PETSC_TRUE
   532)           this%isotope_daughters(icount,isotope%iisotope) = isotope2%iisotope
   533)           this%isotope_daughters(0,isotope%iisotope) = icount
   534)           exit
   535)         endif
   536)         isotope2 => isotope2%next
   537)       enddo
   538)       if (.not.found) then
   539)         option%io_buffer = 'Daughter "' // trim(daughter%name) // &
   540)                            '" not found among isotope list.'
   541)         call printErrMsg(option)
   542)       endif
   543)       daughter => daughter%next
   544)     enddo
   545)     isotope => isotope%next
   546)   enddo
   547)   
   548)   allocate(this%isotope_parents(1,this%num_isotopes))
   549)   this%isotope_parents = 0
   550)   do iisotope = 1, this%num_isotopes
   551)     do icount = 1, this%isotope_daughters(0,iisotope)
   552)       this%isotope_parents(1,this%isotope_daughters(icount,iisotope)) = &
   553)         this%isotope_parents(1,this%isotope_daughters(icount,iisotope)) + 1
   554)     enddo
   555)   enddo
   556)   max_parents_per_isotope = maxval(this%isotope_parents)
   557)   deallocate(this%isotope_parents)
   558)   allocate(this%isotope_parents(0:max_parents_per_isotope,this%num_isotopes))
   559)   this%isotope_parents = 0
   560)   do iisotope = 1, this%num_isotopes
   561)     do icount = 1, this%isotope_daughters(0,iisotope)
   562)       this%isotope_parents(0,this%isotope_daughters(icount,iisotope)) = &
   563)         this%isotope_parents(0,this%isotope_daughters(icount,iisotope)) + 1
   564)       this%isotope_parents(this%isotope_parents(0, &
   565)                                     this%isotope_daughters(icount,iisotope)), &
   566)                            this%isotope_daughters(icount,iisotope)) = iisotope
   567)     enddo
   568)   enddo
   569)   
   570) #if 0  
   571)   this%num_elements = 3
   572)   max_num_isotopes_per_element = 3
   573)   this%num_isotopes = 6
   574)   this%isotope_to_primary_species = [1, 2, 3, 4, 5, 6]
   575)   this%element_isotopes = 0
   576)   this%element_isotopes(0:1,1) = [1,1]
   577)   this%element_isotopes(0:3,2) = [3,2,3,4]
   578)   this%element_isotopes(0:2,3) = [2,5,6]
   579)   this%isotope_decay_rate = 0.d0
   580)   this%isotope_decay_rate(:) = [1.29d-15, 5.08d-11, 1.03d-14, 1.38d-13, 2.78d-12,2.78d-13]
   581)   allocate(this%isotope_daughters(0:max_daughters,this%num_isotopes))
   582)   this%isotope_daughters = 0
   583)   this%isotope_daughters(:,:) = 0
   584)   this%isotope_daughters(0,1) = 1  ! 1 -> 3
   585)   this%isotope_daughters(1,1) = 3
   586)   this%isotope_daughters(0,3) = 1   ! 3 -> 5
   587)   this%isotope_daughters(1,3) = 5
   588)   allocate(this%element_solubility(this%num_elements))
   589)   this%element_solubility = 0.d0
   590) #endif
   591) 
   592) end subroutine PMUFDDecayInit
   593) 
   594) ! ************************************************************************** !
   595) 
   596) subroutine PMUFDDecaySetRealization(this,realization)
   597)   ! 
   598)   ! Author: Glenn Hammond
   599)   ! Date: 06/24/15
   600) 
   601)   use Realization_Subsurface_class
   602) 
   603)   implicit none
   604)   
   605)   class(pm_ufd_decay_type) :: this
   606)   class(realization_subsurface_type), pointer :: realization
   607)   
   608)   this%realization => realization
   609)   this%realization_base => realization
   610) 
   611) end subroutine PMUFDDecaySetRealization
   612) 
   613) ! ************************************************************************** !
   614) 
   615) recursive subroutine PMUFDDecayInitializeRun(this)
   616)   ! 
   617)   ! Initializes the time stepping
   618)   ! 
   619)   ! Author: Glenn Hammond
   620)   ! Date: 06/24/15
   621) 
   622)   use Patch_module
   623)   use Grid_module
   624)   use Reactive_Transport_Aux_module
   625)   
   626)   implicit none
   627) 
   628)   class(pm_ufd_decay_type) :: this
   629)   
   630)   type(patch_type), pointer :: patch
   631)   type(grid_type), pointer :: grid
   632)   type(reactive_transport_auxvar_type), pointer :: rt_auxvars(:)
   633)   
   634)   PetscReal :: kd_kgw_m3b
   635)   PetscInt :: local_id, ghosted_id
   636)   PetscInt :: iele, iiso, ipri, i, imat
   637)   
   638)   patch => this%realization%patch
   639)   grid => patch%grid
   640)   rt_auxvars => patch%aux%RT%auxvars
   641)   
   642)   ! set initial sorbed concentration in equilibrium with aqueous phase
   643)   do local_id = 1, grid%nlmax
   644)     ghosted_id = grid%nL2G(local_id)
   645)     imat = patch%imat(ghosted_id) 
   646)     if (imat <= 0) cycle
   647)     do iele = 1, this%num_elements
   648)       kd_kgw_m3b = this%element_Kd(iele,imat)
   649)       do i = 1, this%element_isotopes(0,iele)
   650)         iiso = this%element_isotopes(i,iele)
   651)         ipri = this%isotope_to_primary_species(iiso)
   652)         rt_auxvars(ghosted_id)%total_sorb_eq(ipri) = &
   653)           rt_auxvars(ghosted_id)%pri_molal(ipri) * kd_kgw_m3b
   654)       enddo
   655)     enddo      
   656)   enddo
   657) 
   658) end subroutine PMUFDDecayInitializeRun
   659) 
   660) ! ************************************************************************** !
   661) 
   662) subroutine PMUFDDecayInitializeTimestep(this)
   663)   ! 
   664)   ! Author: Glenn Hammond
   665)   ! Date: 06/24/15
   666) 
   667)   use Global_module
   668)   
   669)   implicit none
   670)   
   671)   class(pm_ufd_decay_type) :: this
   672) 
   673)   if (this%option%print_screen_flag) then
   674)     write(*,'(/,2("=")," USED FUEL DISPOSITION DECAY MODEL ",43("="))')
   675)   endif
   676) 
   677) end subroutine PMUFDDecayInitializeTimestep
   678) 
   679) ! ************************************************************************** !
   680) 
   681) subroutine PMUFDDecayPreSolve(this)
   682)   ! 
   683)   ! Author: Glenn Hammond
   684)   ! Date: 06/24/15
   685) 
   686)   use Grid_module
   687)   use Global_Aux_module
   688)   use Reactive_Transport_Aux_module
   689) 
   690)   implicit none
   691)   
   692)   class(pm_ufd_decay_type) :: this
   693)   
   694)   type(grid_type), pointer :: grid
   695)   type(global_auxvar_type), pointer :: global_auxvars(:)
   696)   type(reactive_transport_auxvar_type), pointer :: rt_auxvars(:)
   697)   PetscInt :: local_id
   698)   PetscInt :: ghosted_id
   699)   
   700)   grid => this%realization%patch%grid
   701)   global_auxvars => this%realization%patch%aux%Global%auxvars
   702)   rt_auxvars => this%realization%patch%aux%RT%auxvars
   703)   
   704) end subroutine PMUFDDecayPreSolve
   705) 
   706) ! ************************************************************************** !
   707) 
   708) subroutine PMUFDDecaySolve(this,time,ierr)
   709)   ! 
   710)   ! Author: Glenn Hammond
   711)   ! Date: 06/24/15
   712)   !
   713)   use Option_module
   714)   use Reaction_Aux_module
   715)   use Patch_module
   716)   use Grid_module
   717)   use Field_module
   718)   use Reactive_Transport_Aux_module
   719)   use Global_Aux_module
   720)   use Material_Aux_class
   721)   
   722)   implicit none
   723) 
   724) #include "petsc/finclude/petscvec.h"
   725) #include "petsc/finclude/petscvec.h90"
   726) 
   727)   class(pm_ufd_decay_type) :: this
   728)   
   729)   PetscReal :: time
   730)   PetscErrorCode :: ierr
   731)   
   732)   type(option_type), pointer :: option
   733)   type(reaction_type), pointer :: reaction
   734)   type(patch_type), pointer :: patch
   735)   type(grid_type), pointer :: grid
   736)   type(field_type), pointer :: field
   737)   type(reactive_transport_auxvar_type), pointer :: rt_auxvars(:)
   738)   type(global_auxvar_type), pointer :: global_auxvars(:)
   739)   class(material_auxvar_type), pointer :: material_auxvars(:)
   740)   PetscInt :: local_id
   741)   PetscInt :: ghosted_id
   742)   PetscInt :: iele, i, p, g, ip, ig, iiso, ipri, imnrl, imat
   743)   PetscReal :: dt
   744)   PetscReal :: vol, por, sat, den_w_kg, vps
   745)   PetscReal :: conc_iso_aq0, conc_iso_sorb0, conc_iso_ppt0
   746)   PetscReal :: conc_ele_aq1, conc_ele_sorb1, conc_ele_ppt1
   747)   PetscReal :: mass_iso_aq0, mass_iso_sorb0, mass_iso_ppt0
   748)   PetscReal :: mass_ele_aq1, mass_ele_sorb1, mass_ele_ppt1, mass_cmp_tot1
   749)   PetscReal :: mass_iso_tot0(this%num_isotopes), mass_iso_tot_star(this%num_isotopes)
   750)   PetscReal :: mass_iso_tot1(this%num_isotopes), delta_mass_iso_tot(this%num_isotopes)
   751)   PetscReal :: mass_ele_tot1
   752)   PetscReal :: coeff(this%num_isotopes)
   753)   PetscReal :: mass_old(this%num_isotopes)
   754)   PetscReal :: mol_fraction_iso(this%num_isotopes)
   755)   PetscReal :: kd_kgw_m3b
   756)   PetscBool :: above_solubility
   757)   PetscReal, pointer :: xx_p(:)
   758) 
   759)   ierr = 0
   760)   
   761)   option => this%realization%option
   762)   reaction => this%realization%reaction
   763)   patch => this%realization%patch
   764)   field => this%realization%field
   765)   grid => patch%grid
   766)   rt_auxvars => patch%aux%RT%auxvars
   767)   global_auxvars => patch%aux%Global%auxvars
   768)   material_auxvars => patch%aux%Material%auxvars
   769)   
   770)   dt = option%tran_dt
   771)   call VecGetArrayF90(field%tran_xx,xx_p,ierr);CHKERRQ(ierr)
   772) 
   773)   do local_id = 1, grid%nlmax
   774)     ghosted_id = grid%nL2G(local_id)
   775)     imat = patch%imat(ghosted_id)
   776)     if (imat <= 0) cycle
   777)     vol = material_auxvars(ghosted_id)%volume
   778)     den_w_kg = global_auxvars(ghosted_id)%den_kg(1)
   779)     por = material_auxvars(ghosted_id)%porosity
   780)     sat = global_auxvars(ghosted_id)%sat(1)
   781)     vps = vol * por * sat ! m^3 water
   782)     
   783)     ! sum up mass of each isotope across phases and decay
   784)     do iele = 1, this%num_elements
   785)       do i = 1, this%element_isotopes(0,iele)
   786)         iiso = this%element_isotopes(i,iele)
   787)         ipri = this%isotope_to_primary_species(iiso)
   788)         imnrl = this%isotope_to_mineral(iiso)
   789)         ! # indicated time level (0 = prev time level, 1 = new time level) 
   790)         conc_iso_aq0 = xx_p((local_id-1)*reaction%ncomp+ipri) * &
   791)                        den_w_kg / 1000.d0  ! mol/L
   792)         !conc_iso_aq0 = rt_auxvars(ghosted_id)%total(ipri,1) ! mol/L
   793)         conc_iso_sorb0 = rt_auxvars(ghosted_id)%total_sorb_eq(ipri) ! mol/m^3 bulk
   794)         conc_iso_ppt0 = rt_auxvars(ghosted_id)%mnrl_volfrac(imnrl) ! m^3 mnrl/m^3 bulk
   795)         mass_iso_aq0 = conc_iso_aq0*vps*1.d3 ! mol/L * m^3 water * 1000 L /m^3 = mol
   796)         mass_iso_sorb0 = conc_iso_sorb0 * vol ! mol/m^3 bulk * m^3 bulk = mol
   797)         mass_iso_ppt0 = conc_iso_ppt0 * vol / &  ! m^3 mnrl/m^3 bulk * m^3 bulk / (m^3 mnrl/mol mnrl) = mol
   798)                         reaction%mineral%kinmnrl_molar_vol(imnrl)
   799)         mass_iso_tot0(iiso) = mass_iso_aq0 + mass_iso_sorb0 + mass_iso_ppt0
   800)       enddo
   801)     enddo 
   802)     
   803)     ! save the mass from the previous time step:
   804)     mass_old(:) = mass_iso_tot0(:)
   805) 
   806)     ! FIRST PASS decay ==============================================
   807)     do i = 1,this%num_isotopes
   808)       ! update the initial value of the isotope coefficient:
   809)       coeff(i) = mass_old(i)
   810)       ! loop through the isotope's parents:
   811)       do p = 1,this%isotope_parents(0,i)
   812)         ip = this%isotope_parents(p,i)
   813)         coeff(i) = coeff(i) - (this%isotope_decay_rate(ip) * mass_old(ip)) / &
   814)           (this%isotope_decay_rate(i) - this%isotope_decay_rate(ip))
   815)         ! loop through the isotope's parent's parents:
   816)         do g = 1,this%isotope_parents(0,p)
   817)           ig = this%isotope_parents(g,p)
   818)           coeff(i) = coeff(i) - &
   819)             ((this%isotope_decay_rate(ip) * this%isotope_decay_rate(ig) * &        
   820)             mass_old(ig)) / ((this%isotope_decay_rate(ip) - &
   821)             this%isotope_decay_rate(ig)) * (this%isotope_decay_rate(i) - &
   822)             this%isotope_decay_rate(ig)))) + ((this%isotope_decay_rate(ip) * &
   823)             this%isotope_decay_rate(ig) * mass_old(ig)) / &
   824)             ((this%isotope_decay_rate(ip) - this%isotope_decay_rate(ig)) * &
   825)             (this%isotope_decay_rate(i) - this%isotope_decay_rate(ip))))
   826)         enddo ! grandparent loop
   827)       enddo ! parent loop
   828)     enddo ! isotope loop
   829)     ! SECOND PASS decay =============================================
   830)     do i = 1,this%num_isotopes
   831)       ! decay the isotope species:
   832)       mass_iso_tot1(i) = coeff(i)*exp(-1.d0*this%isotope_decay_rate(i)*dt)
   833)       ! loop through the isotope's parents:
   834)       do p = 1,this%isotope_parents(0,i)
   835)         ip = this%isotope_parents(p,i)
   836)         mass_iso_tot1(i) = mass_iso_tot1(i) + &
   837)               (((this%isotope_decay_rate(ip) * mass_old(ip)) / &
   838)               (this%isotope_decay_rate(i) - this%isotope_decay_rate(ip))) * &
   839)                exp(-1.d0 * this%isotope_decay_rate(ip) * dt))
   840)         ! loop through the isotope's parent's parents:
   841)         do g = 1,this%isotope_parents(0,p)
   842)           ig = this%isotope_parents(g,p)
   843)           mass_iso_tot1(i) = mass_iso_tot1(i) - &
   844)             ((this%isotope_decay_rate(ip) * this%isotope_decay_rate(ig) * &
   845)             mass_old(ig) * exp(-1.d0 * this%isotope_decay_rate(ip) * dt)) / &
   846)             ((this%isotope_decay_rate(ip) - this%isotope_decay_rate(ig)) * &
   847)             (this%isotope_decay_rate(i) - this%isotope_decay_rate(ip)))) + &
   848)             ((this%isotope_decay_rate(ip) * this%isotope_decay_rate(ig) * &
   849)             mass_old(ig) * exp(-1.d0 * this%isotope_decay_rate(ig) * dt)) / &
   850)             ((this%isotope_decay_rate(ip) - this%isotope_decay_rate(ig)) * &
   851)             (this%isotope_decay_rate(i) - this%isotope_decay_rate(ig))))
   852)         enddo ! grandparent loop
   853)       enddo ! parent loop
   854)     enddo ! isotope loop
   855)     
   856)     mass_iso_tot1 = max(mass_iso_tot1,1.d-90)
   857) 
   858)     do iele = 1, this%num_elements
   859)       ! calculate mole fractions
   860)       mass_ele_tot1 = 0.d0
   861)       mol_fraction_iso = 0.d0
   862)       do i = 1, this%element_isotopes(0,iele)
   863)         iiso = this%element_isotopes(i,iele)
   864)         mass_ele_tot1 = mass_ele_tot1 + mass_iso_tot1(iiso)
   865)       enddo
   866)       do i = 1, this%element_isotopes(0,iele)
   867)         iiso = this%element_isotopes(i,iele)
   868)         mol_fraction_iso(i) = mass_iso_tot1(iiso) / mass_ele_tot1
   869)       enddo
   870)       ! split mass between phases
   871)       kd_kgw_m3b = this%element_Kd(iele,imat)
   872)       conc_ele_aq1 = mass_ele_tot1 / (1.d0+kd_kgw_m3b/(den_w_kg*por*sat)) / &
   873)                          (vps*1.d3)
   874)       above_solubility = conc_ele_aq1 > this%element_solubility(iele)
   875)       if (above_solubility) then
   876)         conc_ele_aq1 = this%element_solubility(iele)
   877)       endif
   878)       ! assume identical sorption for all isotopes in element
   879)       conc_ele_sorb1 = conc_ele_aq1 / den_w_kg * 1.d3 * kd_kgw_m3b
   880)       mass_ele_aq1 = conc_ele_aq1*vps*1.d3
   881)       mass_ele_sorb1 = conc_ele_sorb1 * vol
   882)       ! roundoff can erroneously result in precipitate. this conditional avoids
   883)       ! such an issue
   884)       if (above_solubility) then
   885)         mass_ele_ppt1 = max(mass_ele_tot1 - mass_ele_aq1 - mass_ele_sorb1,0.d0)
   886)       else
   887)         mass_ele_ppt1 = 0.d0
   888)       endif
   889)       conc_ele_ppt1 = mass_ele_ppt1 * &
   890)                       reaction%mineral%kinmnrl_molar_vol(imnrl) / vol
   891)       ! store mass in data structures
   892)       do i = 1, this%element_isotopes(0,iele)
   893)         iiso = this%element_isotopes(i,iele)
   894)         ipri = this%isotope_to_primary_species(iiso)
   895)         imnrl = this%isotope_to_mineral(iiso)
   896)         rt_auxvars(ghosted_id)%total(ipri,1) = conc_ele_aq1 * mol_fraction_iso(i)
   897)         rt_auxvars(ghosted_id)%pri_molal(ipri) = conc_ele_aq1 / den_w_kg * 1.d3 * mol_fraction_iso(i)
   898)         rt_auxvars(ghosted_id)%total_sorb_eq(ipri) = conc_ele_sorb1 * mol_fraction_iso(i)
   899)         rt_auxvars(ghosted_id)%mnrl_volfrac(imnrl) = conc_ele_ppt1 * mol_fraction_iso(i)
   900)         ! need to copy primary molalities back into transport solution Vec
   901)         xx_p((local_id-1)*reaction%ncomp+ipri) = rt_auxvars(ghosted_id)%pri_molal(ipri)
   902)       enddo
   903)     enddo      
   904)   enddo
   905) 
   906)   call VecRestoreArrayF90(field%tran_xx,xx_p,ierr);CHKERRQ(ierr)
   907)   if (reaction%use_log_formulation) then
   908)     call VecCopy(field%tran_xx,field%tran_log_xx,ierr);CHKERRQ(ierr)
   909)     call VecLog(field%tran_log_xx,ierr);CHKERRQ(ierr)
   910)   endif  
   911) !  call DiscretizationGlobalToLocal(this%realization%discretization, &
   912) !                                   field%tran_xx,field%tran_xx_loc,NTRANDOF)
   913)   
   914) end subroutine PMUFDDecaySolve
   915) 
   916) ! ************************************************************************** !
   917) 
   918) subroutine PMUFDDecayPostSolve(this)
   919)   ! 
   920)   ! PMUFDDecayUpdatePostSolve:
   921)   ! 
   922)   ! Author: Glenn Hammond
   923)   ! Date: 06/24/15
   924)   ! 
   925)   implicit none
   926)   
   927)   class(pm_ufd_decay_type) :: this
   928)   
   929) end subroutine PMUFDDecayPostSolve
   930) 
   931) ! ************************************************************************** !
   932) 
   933) function PMUFDDecayAcceptSolution(this)
   934)   ! 
   935)   ! Author: Glenn Hammond
   936)   ! Date: 06/24/15
   937) 
   938)   implicit none
   939)   
   940)   class(pm_ufd_decay_type) :: this
   941)   
   942)   PetscBool :: PMUFDDecayAcceptSolution
   943)   
   944)   ! do nothing
   945)   PMUFDDecayAcceptSolution = PETSC_TRUE
   946)   
   947) end function PMUFDDecayAcceptSolution
   948) 
   949) ! ************************************************************************** !
   950) 
   951) subroutine PMUFDDecayUpdatePropertiesTS(this)
   952)   ! 
   953)   ! Updates parameters/properties at each Newton iteration
   954)   !
   955)   ! Author: Glenn Hammond
   956)   ! Date: 06/24/15
   957) 
   958)   implicit none
   959)   
   960)   class(pm_ufd_decay_type) :: this
   961)   
   962) end subroutine PMUFDDecayUpdatePropertiesTS
   963) 
   964) ! ************************************************************************** !
   965) 
   966) subroutine PMUFDDecayTimeCut(this)
   967)   ! 
   968)   ! Author: Glenn Hammond
   969)   ! Date: 06/24/15
   970)   
   971)   implicit none
   972)   
   973)   class(pm_ufd_decay_type) :: this
   974)   
   975)   PetscErrorCode :: ierr
   976)   
   977) end subroutine PMUFDDecayTimeCut
   978) 
   979) ! ************************************************************************** !
   980) 
   981) subroutine PMUFDDecayFinalizeTimestep(this)
   982)   ! 
   983)   ! Author: Glenn Hammond
   984)   ! Date: 06/24/15
   985) 
   986)   implicit none
   987)   
   988)   class(pm_ufd_decay_type) :: this
   989)   
   990) end subroutine PMUFDDecayFinalizeTimestep
   991) 
   992) ! ************************************************************************** !
   993) 
   994) subroutine PMUFDDecayUpdateSolution(this)
   995)   ! 
   996)   ! Author: Glenn Hammond
   997)   ! Date: 06/24/15
   998) 
   999)   implicit none
  1000)   
  1001)   class(pm_ufd_decay_type) :: this
  1002)   
  1003)   PetscErrorCode :: ierr
  1004) 
  1005) end subroutine PMUFDDecayUpdateSolution  
  1006) 
  1007) ! ************************************************************************** !
  1008) 
  1009) subroutine PMUFDDecayUpdateAuxvars(this)
  1010)   ! 
  1011)   ! Author: Glenn Hammond
  1012)   ! Date: 06/24/15
  1013) 
  1014)   implicit none
  1015)   
  1016)   class(pm_ufd_decay_type) :: this
  1017) 
  1018)   this%option%io_buffer = 'PMUFDDecayUpdateAuxvars() must be extended.'
  1019)   call printErrMsg(this%option)
  1020) 
  1021) end subroutine PMUFDDecayUpdateAuxvars   
  1022) 
  1023) ! ************************************************************************** !
  1024) 
  1025) subroutine PMUFDDecayCheckpoint(this,viewer)
  1026)   ! 
  1027)   ! Checkpoints data associated with Subsurface PM
  1028)   ! 
  1029)   ! Author: Glenn Hammond
  1030)   ! Date: 06/24/15
  1031) 
  1032)   implicit none
  1033) #include "petsc/finclude/petscviewer.h"      
  1034) 
  1035)   class(pm_ufd_decay_type) :: this
  1036)   PetscViewer :: viewer
  1037)   
  1038) end subroutine PMUFDDecayCheckpoint
  1039) 
  1040) ! ************************************************************************** !
  1041) 
  1042) subroutine PMUFDDecayRestart(this,viewer)
  1043)   ! 
  1044)   ! Restarts data associated with Subsurface PM
  1045)   ! 
  1046)   ! Author: Glenn Hammond
  1047)   ! Date: 06/24/15
  1048) 
  1049)   implicit none
  1050) #include "petsc/finclude/petscviewer.h"      
  1051) 
  1052)   class(pm_ufd_decay_type) :: this
  1053)   PetscViewer :: viewer
  1054)   
  1055) end subroutine PMUFDDecayRestart
  1056) 
  1057) ! ************************************************************************** !
  1058) 
  1059) recursive subroutine PMUFDDecayFinalizeRun(this)
  1060)   ! 
  1061)   ! Finalizes the time stepping
  1062)   ! 
  1063)   ! Author: Glenn Hammond
  1064)   ! Date: 06/24/15
  1065)   
  1066)   implicit none
  1067)   
  1068)   class(pm_ufd_decay_type) :: this
  1069)   
  1070)   ! do something here
  1071)   
  1072)   if (associated(this%next)) then
  1073)     call this%next%FinalizeRun()
  1074)   endif  
  1075)   
  1076) end subroutine PMUFDDecayFinalizeRun
  1077) 
  1078) ! ************************************************************************** !
  1079) 
  1080) subroutine PMUFDDecayInputRecord(this)
  1081)   ! 
  1082)   ! Writes ingested information to the input record file.
  1083)   ! 
  1084)   ! Author: Jenn Frederick, SNL
  1085)   ! Date: 03/21/2016
  1086)   ! 
  1087)   
  1088)   implicit none
  1089)   
  1090)   class(pm_ufd_decay_type) :: this
  1091) 
  1092)   character(len=MAXWORDLENGTH) :: word
  1093)   PetscInt :: id
  1094) 
  1095)   id = INPUT_RECORD_UNIT
  1096) 
  1097)   write(id,'(a29)',advance='no') 'pm: '
  1098)   write(id,'(a)') this%name
  1099) 
  1100) end subroutine PMUFDDecayInputRecord
  1101) 
  1102) ! ************************************************************************** !
  1103) 
  1104) subroutine PMUFDDecayDestroy(this)
  1105)   ! 
  1106)   ! Destroys Subsurface process model
  1107)   ! 
  1108)   ! Author: Glenn Hammond
  1109)   ! Date: 06/24/15
  1110)   use Utility_module, only : DeallocateArray
  1111)   use Option_module
  1112) 
  1113)   implicit none
  1114)   
  1115)   class(pm_ufd_decay_type) :: this
  1116)   
  1117)   type(element_type), pointer :: cur_element, prev_element
  1118)   type(isotope_type), pointer :: cur_isotope, prev_isotope
  1119)   type(daughter_type), pointer :: cur_daughter, prev_daughter
  1120)     
  1121)   call DeallocateArray(this%element)
  1122)   call DeallocateArray(this%element_isotopes)
  1123)   call DeallocateArray(this%isotope_to_primary_species)
  1124)   call DeallocateArray(this%isotope_to_mineral)
  1125)   call DeallocateArray(this%isotope_decay_rate)
  1126)   call DeallocateArray(this%isotope_daughters)
  1127)   call DeallocateArray(this%isotope_parents)
  1128)   call DeallocateArray(this%element_solubility)
  1129)   call DeallocateArray(this%element_Kd)
  1130)   call DeallocateArray(this%num_isotopes_per_element)
  1131)   
  1132)   cur_isotope => this%isotope_list
  1133)   do
  1134)     if (.not.associated(cur_isotope)) exit
  1135)     cur_daughter => cur_isotope%daughter_list
  1136)     do
  1137)       if (.not.associated(cur_daughter)) exit
  1138)       prev_daughter => cur_daughter
  1139)       cur_daughter => cur_daughter%next
  1140)       nullify(prev_daughter%next)
  1141)       deallocate(prev_daughter)
  1142)       nullify(prev_daughter)
  1143)     enddo
  1144)     prev_isotope => cur_isotope
  1145)     cur_isotope => cur_isotope%next
  1146)     nullify(prev_isotope%next)
  1147)     deallocate(prev_isotope)
  1148)     nullify(prev_isotope)
  1149)   enddo
  1150)   
  1151)   cur_element => this%element_list
  1152)   do
  1153)     if (.not.associated(cur_element)) exit
  1154)     prev_element => cur_element
  1155)     cur_element => cur_element%next
  1156)     call DeallocateArray(prev_element%Kd)
  1157)     call DeallocateArray(prev_element%Kd_material_name)
  1158)     nullify(prev_element%next)
  1159)     deallocate(prev_element)
  1160)     nullify(prev_element)
  1161)   enddo
  1162)   
  1163) end subroutine PMUFDDecayDestroy
  1164)   
  1165) end module PM_UFD_Decay_class

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