pm_mphase.F90       coverage:  73.33 %func     51.76 %block


     1) module PM_Mphase_class
     2) 
     3)   use PM_Base_class
     4)   use PM_Subsurface_Flow_class
     5)   
     6)   use PFLOTRAN_Constants_module
     7) 
     8)   implicit none
     9) 
    10)   private
    11) 
    12) #include "petsc/finclude/petscsys.h"
    13) 
    14) #include "petsc/finclude/petscvec.h"
    15) #include "petsc/finclude/petscvec.h90"
    16) #include "petsc/finclude/petscmat.h"
    17) #include "petsc/finclude/petscmat.h90"
    18) #include "petsc/finclude/petscsnes.h"
    19) 
    20)   type, public, extends(pm_subsurface_flow_type) :: pm_mphase_type
    21)   contains
    22)     procedure, public :: Read => PMMphaseRead
    23)     procedure, public :: InitializeTimestep => PMMphaseInitializeTimestep
    24)     procedure, public :: Residual => PMMphaseResidual
    25)     procedure, public :: Jacobian => PMMphaseJacobian
    26)     procedure, public :: UpdateTimestep => PMMphaseUpdateTimestep
    27)     procedure, public :: PreSolve => PMMphasePreSolve
    28)     procedure, public :: PostSolve => PMMphasePostSolve
    29) #if 0
    30)     procedure, public :: CheckUpdatePre => PMMphaseCheckUpdatePre
    31)     procedure, public :: CheckUpdatePost => PMMphaseCheckUpdatePost
    32) #endif
    33)     procedure, public :: TimeCut => PMMphaseTimeCut
    34)     procedure, public :: UpdateSolution => PMMphaseUpdateSolution
    35)     procedure, public :: UpdateAuxVars => PMMphaseUpdateAuxVars
    36)     procedure, public :: MaxChange => PMMphaseMaxChange
    37)     procedure, public :: ComputeMassBalance => PMMphaseComputeMassBalance
    38)     procedure, public :: InputRecord => PMMphaseInputRecord
    39)     procedure, public :: Destroy => PMMphaseDestroy
    40)   end type pm_mphase_type
    41)   
    42)   public :: PMMphaseCreate
    43)   
    44) contains
    45) 
    46) ! ************************************************************************** !
    47) 
    48) function PMMphaseCreate()
    49)   ! 
    50)   ! Creates Mphase process models shell
    51)   ! 
    52)   ! Author: Glenn Hammond
    53)   ! Date: 03/14/13
    54)   ! 
    55) 
    56)   implicit none
    57)   
    58)   class(pm_mphase_type), pointer :: PMMphaseCreate
    59) 
    60)   class(pm_mphase_type), pointer :: mphase_pm
    61) 
    62)   allocate(mphase_pm)
    63) 
    64)   call PMSubsurfaceFlowCreate(mphase_pm)
    65)   mphase_pm%name = 'PMMphase'
    66) 
    67)   PMMphaseCreate => mphase_pm
    68)   
    69) end function PMMphaseCreate
    70) 
    71) ! ************************************************************************** !
    72) 
    73) subroutine PMMphaseRead(this,input)
    74)   ! 
    75)   ! Reads input file parameters associated with the Mphase process model
    76)   ! 
    77)   ! Author: Glenn Hammond
    78)   ! Date: 01/29/15
    79)   use Input_Aux_module
    80)   use String_module
    81)   use Utility_module
    82)   use EOS_Water_module  
    83)   use Option_module
    84)   use Mphase_Aux_module
    85)  
    86)   implicit none
    87)   
    88)   class(pm_mphase_type) :: this
    89)   type(input_type), pointer :: input
    90)   
    91)   character(len=MAXWORDLENGTH) :: word
    92)   character(len=MAXSTRINGLENGTH) :: error_string
    93)   type(option_type), pointer :: option
    94)   PetscBool :: found
    95) 
    96)   option => this%option
    97)   
    98)   error_string = 'Mphase Options'
    99)   
   100)   input%ierr = 0
   101)   do
   102)   
   103)     call InputReadPflotranString(input,option)
   104)     if (InputError(input)) exit
   105)     if (InputCheckExit(input,option)) exit
   106)     
   107)     call InputReadWord(input,option,word,PETSC_TRUE)
   108)     call InputErrorMsg(input,option,'keyword',error_string)
   109)     call StringToUpper(word)
   110) 
   111)     found = PETSC_FALSE
   112)     call PMSubsurfaceFlowReadSelectCase(this,input,word,found,option)
   113)     if (found) cycle
   114)     
   115)     select case(trim(word))
   116)       case default
   117)         call InputKeywordUnrecognized(word,error_string,option)
   118)     end select
   119)   enddo
   120)   
   121) end subroutine PMMphaseRead
   122) 
   123) ! ************************************************************************** !
   124) 
   125) subroutine PMMphaseInitializeTimestep(this)
   126)   ! 
   127)   ! Should not need this as it is called in PreSolve.
   128)   ! 
   129)   ! Author: Glenn Hammond
   130)   ! Date: 03/14/13
   131)   ! 
   132) 
   133)   use Mphase_module, only : MphaseInitializeTimestep
   134)   
   135)   implicit none
   136)   
   137)   class(pm_mphase_type) :: this
   138) 
   139)   call PMSubsurfaceFlowInitializeTimestepA(this)         
   140) 
   141)   if (this%option%print_screen_flag) then
   142)     write(*,'(/,2("=")," MPHASE FLOW ",65("="))')
   143)   endif
   144)   
   145)   call MphaseInitializeTimestep(this%realization)
   146)   call PMSubsurfaceFlowInitializeTimestepB(this)         
   147)   
   148) end subroutine PMMphaseInitializeTimestep
   149) 
   150) ! ************************************************************************** !
   151) 
   152) subroutine PMMphasePreSolve(this)
   153)   ! 
   154)   ! Author: Glenn Hammond
   155)   ! Date: 03/14/13
   156)   use Option_module
   157)   use Reaction_Aux_module
   158)   use Reactive_Transport_Aux_module
   159)   use Grid_module
   160)   use Patch_module
   161)   use Global_Aux_module
   162)   use Coupler_module
   163)   use Connection_module  
   164) 
   165)   implicit none
   166)   
   167)   class(pm_mphase_type) :: this
   168) 
   169)   type(reaction_type), pointer :: reaction
   170)   type(patch_type), pointer :: patch
   171)   type(option_type), pointer :: option
   172)   type(grid_type), pointer :: grid
   173)   type(reactive_transport_auxvar_type), pointer :: rt_auxvars(:)
   174)   type(global_auxvar_type), pointer :: global_auxvars(:)
   175)   type(coupler_type), pointer :: boundary_condition
   176)   type(connection_set_type), pointer :: cur_connection_set
   177)   PetscInt :: local_id
   178)   PetscInt :: ghosted_id
   179)   PetscInt :: sum_connection
   180)   PetscInt :: iconn
   181)   PetscInt :: na_id, cl_id
   182)   
   183)   reaction => this%realization%reaction
   184)   option => this%realization%option
   185) #if 1
   186)   if (associated(reaction)) then
   187)     if (associated(reaction%species_idx)) then
   188)       patch => this%realization%patch
   189)       global_auxvars => patch%aux%Global%auxvars
   190)       if (associated(global_auxvars(1)%m_nacl)) then
   191)         na_id = reaction%species_idx%na_ion_id
   192)         cl_id = reaction%species_idx%cl_ion_id
   193)  
   194)         grid => patch%grid
   195)         rt_auxvars => patch%aux%RT%auxvars
   196)         global_auxvars => patch%aux%Global%auxvars
   197) 
   198)         if (na_id > 0 .and. cl_id > 0) then
   199)           do ghosted_id = 1, grid%ngmax
   200)             if (grid%nG2L(ghosted_id) < 0) cycle ! bypass ghosted corner cells
   201)             !geh - Ignore inactive cells with inactive materials
   202)             if (patch%imat(ghosted_id) <= 0) cycle
   203)             global_auxvars(ghosted_id)%m_nacl(1) = &
   204)               rt_auxvars(ghosted_id)%pri_molal(na_id)
   205)             global_auxvars(ghosted_id)%m_nacl(2) = &
   206)               rt_auxvars(ghosted_id)%pri_molal(cl_id)
   207)           enddo
   208)         else    
   209)           do ghosted_id = 1, grid%ngmax
   210)             if (grid%nG2L(ghosted_id) < 0) cycle ! bypass ghosted corner cells
   211)             !geh - Ignore inactive cells with inactive materials
   212)             if (patch%imat(ghosted_id) <= 0) cycle
   213)             global_auxvars(ghosted_id)%m_nacl = option%m_nacl
   214)           enddo
   215)         endif
   216)       
   217)         rt_auxvars => this%realization%patch%aux%RT%auxvars_bc
   218)         global_auxvars => this%realization%patch%aux%Global%auxvars_bc
   219)     
   220)         boundary_condition => patch%boundary_condition_list%first
   221)         sum_connection = 0 
   222)         do 
   223)           if (.not.associated(boundary_condition)) exit
   224)           cur_connection_set => boundary_condition%connection_set
   225)           if (na_id > 0 .and. cl_id > 0) then
   226)             do iconn = 1, cur_connection_set%num_connections
   227)               sum_connection = sum_connection + 1
   228)               local_id = cur_connection_set%id_dn(iconn)
   229)               ghosted_id = grid%nL2G(local_id)
   230)               if (patch%imat(ghosted_id) <= 0) cycle
   231)               global_auxvars(sum_connection)%m_nacl(1) = &
   232)                 rt_auxvars(sum_connection)%pri_molal(na_id)
   233)               global_auxvars(sum_connection)%m_nacl(2) = &
   234)                 rt_auxvars(sum_connection)%pri_molal(cl_id)
   235)             enddo
   236)           else    
   237)             do iconn = 1, cur_connection_set%num_connections
   238)               sum_connection = sum_connection + 1
   239)               local_id = cur_connection_set%id_dn(iconn)
   240)               ghosted_id = grid%nL2G(local_id)
   241)               if (patch%imat(ghosted_id) <= 0) cycle
   242)               global_auxvars(sum_connection)%m_nacl = option%m_nacl
   243)             enddo
   244)           endif        
   245)           boundary_condition => boundary_condition%next
   246)         enddo
   247)       endif
   248)     endif
   249)   endif   
   250) #endif
   251) 
   252) end subroutine PMMphasePreSolve
   253) 
   254) ! ************************************************************************** !
   255) 
   256) subroutine PMMphasePostSolve(this)
   257)   ! 
   258)   ! PMMphaseUpdatePostSolve:
   259)   ! 
   260)   ! Author: Glenn Hammond
   261)   ! Date: 03/14/13
   262) 
   263)   implicit none
   264)   
   265)   class(pm_mphase_type) :: this
   266)   
   267) end subroutine PMMphasePostSolve
   268) 
   269) ! ************************************************************************** !
   270) 
   271) subroutine PMMphaseUpdateTimestep(this,dt,dt_min,dt_max,iacceleration, &
   272)                                     num_newton_iterations,tfac)
   273)   ! 
   274)   ! Author: Glenn Hammond
   275)   ! Date: 03/14/13
   276)   ! 
   277) 
   278)   implicit none
   279)   
   280)   class(pm_mphase_type) :: this
   281)   PetscReal :: dt
   282)   PetscReal :: dt_min,dt_max
   283)   PetscInt :: iacceleration
   284)   PetscInt :: num_newton_iterations
   285)   PetscReal :: tfac(:)
   286)   
   287)   PetscReal :: fac
   288)   PetscReal :: ut
   289)   PetscReal :: up
   290)   PetscReal :: utmp
   291)   PetscReal :: uc
   292)   PetscReal :: uus
   293)   PetscReal :: dtt
   294)   PetscReal :: dt_p
   295)   PetscReal :: dt_tfac
   296)   PetscInt :: ifac
   297)   
   298)   if (iacceleration > 0) then
   299)     fac = 0.5d0
   300)     if (num_newton_iterations >= iacceleration) then
   301)       fac = 0.33d0
   302)       ut = 0.d0
   303)     else
   304)       up = this%pressure_change_governor/(this%max_pressure_change+0.1)
   305)       utmp = this%temperature_change_governor/(this%max_temperature_change+1.d-5)
   306)       uc = this%xmol_change_governor/(this%max_xmol_change+1.d-6)
   307)       uus= this%saturation_change_governor/(this%max_saturation_change+1.d-6)
   308)       ut = min(up,utmp,uc,uus)
   309)     endif
   310)     dtt = fac * dt * (1.d0 + ut)
   311)   else
   312)     ifac = max(min(num_newton_iterations,size(tfac)),1)
   313)     dt_tfac = tfac(ifac) * dt
   314) 
   315)     fac = 0.5d0
   316)     up = this%pressure_change_governor/(this%max_pressure_change+0.1)
   317)     dt_p = fac * dt * (1.d0 + up)
   318) 
   319)     dtt = min(dt_tfac,dt_p)
   320)   endif
   321)   
   322)   if (dtt > 2.d0 * dt) dtt = 2.d0 * dt
   323)   if (dtt > dt_max) dtt = dt_max
   324)   ! geh: There used to be code here that cut the time step if it is too
   325)   !      large relative to the simulation time.  This has been removed.
   326)   dtt = max(dtt,dt_min)
   327)   dt = dtt
   328) 
   329)   call PMSubsurfaceFlowLimitDTByCFL(this,dt)
   330)   
   331) end subroutine PMMphaseUpdateTimestep
   332) 
   333) ! ************************************************************************** !
   334) 
   335) subroutine PMMphaseResidual(this,snes,xx,r,ierr)
   336)   ! 
   337)   ! Author: Glenn Hammond
   338)   ! Date: 03/14/13
   339)   ! 
   340) 
   341)   use Mphase_module, only : MphaseResidual
   342) 
   343)   implicit none
   344)   
   345)   class(pm_mphase_type) :: this
   346)   SNES :: snes
   347)   Vec :: xx
   348)   Vec :: r
   349)   PetscErrorCode :: ierr
   350)   
   351)   call MphaseResidual(snes,xx,r,this%realization,ierr)
   352) 
   353) end subroutine PMMphaseResidual
   354) 
   355) ! ************************************************************************** !
   356) 
   357) subroutine PMMphaseJacobian(this,snes,xx,A,B,ierr)
   358)   ! 
   359)   ! Author: Glenn Hammond
   360)   ! Date: 03/14/13
   361)   ! 
   362) 
   363)   use Mphase_module, only : MphaseJacobian
   364) 
   365)   implicit none
   366)   
   367)   class(pm_mphase_type) :: this
   368)   SNES :: snes
   369)   Vec :: xx
   370)   Mat :: A, B
   371)   PetscErrorCode :: ierr
   372)   
   373)   call MphaseJacobian(snes,xx,A,B,this%realization,ierr)
   374) 
   375) end subroutine PMMphaseJacobian
   376)     
   377) #if 0
   378) 
   379) ! ************************************************************************** !
   380) 
   381) subroutine PMMphaseCheckUpdatePre(this,line_search,X,dX,changed,ierr)
   382)   ! 
   383)   ! Author: Glenn Hammond
   384)   ! Date: 03/14/13
   385)   ! 
   386) 
   387)   use Mphase_module, only : MphaseCheckUpdatePre
   388) 
   389)   implicit none
   390)   
   391)   class(pm_mphase_type) :: this
   392)   SNESLineSearch :: line_search
   393)   Vec :: X
   394)   Vec :: dX
   395)   PetscBool :: changed
   396)   PetscErrorCode :: ierr
   397)   
   398)   call MphaseCheckUpdatePre(line_search,X,dX,changed,this%realization,ierr)
   399) 
   400) end subroutine PMMphaseCheckUpdatePre
   401) 
   402) ! ************************************************************************** !
   403) 
   404) subroutine PMMphaseCheckUpdatePost(this,line_search,X0,dX,X1,dX_changed, &
   405)                                    X1_changed,ierr)
   406)   ! 
   407)   ! Author: Glenn Hammond
   408)   ! Date: 03/14/13
   409)   ! 
   410) 
   411)   use Mphase_module, only : MphaseCheckUpdatePost
   412) 
   413)   implicit none
   414)   
   415)   class(pm_mphase_type) :: this
   416)   SNESLineSearch :: line_search
   417)   Vec :: X0
   418)   Vec :: dX
   419)   Vec :: X1
   420)   PetscBool :: dX_changed
   421)   PetscBool :: X1_changed
   422)   PetscErrorCode :: ierr
   423)   
   424)   call MphaseCheckUpdatePost(line_search,X0,dX,X1,dX_changed, &
   425)                                X1_changed,this%realization,ierr)
   426) 
   427) end subroutine PMMphaseCheckUpdatePost
   428) #endif
   429) 
   430) ! ************************************************************************** !
   431) 
   432) subroutine PMMphaseTimeCut(this)
   433)   ! 
   434)   ! Author: Glenn Hammond
   435)   ! Date: 03/14/13
   436)   ! 
   437) 
   438)   use Mphase_module, only : MphaseTimeCut
   439) 
   440)   implicit none
   441)   
   442)   class(pm_mphase_type) :: this
   443)   
   444)   call PMSubsurfaceFlowTimeCut(this)
   445)   call MphaseTimeCut(this%realization)
   446) 
   447) end subroutine PMMphaseTimeCut
   448) 
   449) ! ************************************************************************** !
   450) 
   451) subroutine PMMphaseUpdateSolution(this)
   452)   ! 
   453)   ! Author: Glenn Hammond
   454)   ! Date: 03/14/13
   455)   ! 
   456) 
   457)   use Mphase_module, only : MphaseUpdateSolution
   458) 
   459)   implicit none
   460)   
   461)   class(pm_mphase_type) :: this
   462)   
   463)   call PMSubsurfaceFlowUpdateSolution(this)
   464)   call MphaseUpdateSolution(this%realization)
   465) 
   466) end subroutine PMMphaseUpdateSolution     
   467) 
   468) ! ************************************************************************** !
   469) 
   470) subroutine PMMphaseUpdateAuxVars(this)
   471)   ! 
   472)   ! Author: Glenn Hammond
   473)   ! Date: 04/21/14
   474) 
   475)   use Mphase_module, only : MphaseUpdateAuxVars
   476)   
   477)   implicit none
   478)   
   479)   class(pm_mphase_type) :: this
   480) 
   481)   call MphaseUpdateAuxVars(this%realization)
   482) 
   483) end subroutine PMMphaseUpdateAuxVars   
   484) 
   485) ! ************************************************************************** !
   486) 
   487) subroutine PMMphaseMaxChange(this)
   488)   ! 
   489)   ! Not needed given MphaseMaxChange is called in PostSolve
   490)   ! 
   491)   ! Author: Glenn Hammond
   492)   ! Date: 03/14/13
   493)   ! 
   494) 
   495)   use Mphase_module, only : MphaseMaxChange
   496) 
   497)   implicit none
   498)   
   499)   class(pm_mphase_type) :: this
   500)   
   501)   call MphaseMaxChange(this%realization,this%max_pressure_change, &
   502)                        this%max_temperature_change, &
   503)                        this%max_saturation_change, &
   504)                        this%max_xmol_change)
   505)   if (this%option%print_screen_flag) then
   506)     write(*,'("  --> max chng: dpmx= ",1pe12.4, &
   507)       & " dtmpmx= ",1pe12.4," dcmx= ",1pe12.4," dsmx= ",1pe12.4)') &
   508)           this%max_pressure_change,this%max_temperature_change, &
   509)           this%max_xmol_change,this%max_saturation_change
   510)   endif
   511)   if (this%option%print_file_flag) then
   512)     write(this%option%fid_out,'("  --> max chng: dpmx= ",1pe12.4, &
   513)       & " dtmpmx= ",1pe12.4," dcmx= ",1pe12.4," dsmx= ",1pe12.4)') &
   514)           this%max_pressure_change,this%max_temperature_change, &
   515)           this%max_xmol_change,this%max_saturation_change
   516)   endif   
   517) 
   518) end subroutine PMMphaseMaxChange
   519) 
   520) ! ************************************************************************** !
   521) 
   522) subroutine PMMphaseComputeMassBalance(this,mass_balance_array)
   523)   ! 
   524)   ! Author: Glenn Hammond
   525)   ! Date: 03/14/13
   526)   ! 
   527) 
   528)   use Mphase_module, only : MphaseComputeMassBalance
   529) 
   530)   implicit none
   531)   
   532)   class(pm_mphase_type) :: this
   533)   PetscReal :: mass_balance_array(:)
   534)   
   535)   !geh: currently does not include "trapped" mass
   536)   !call MphaseComputeMassBalance(this%realization,mass_balance_array)
   537) 
   538) end subroutine PMMphaseComputeMassBalance
   539) 
   540) ! ************************************************************************** !
   541) 
   542) subroutine PMMphaseInputRecord(this)
   543)   ! 
   544)   ! Writes ingested information to the input record file.
   545)   ! 
   546)   ! Author: Jenn Frederick, SNL
   547)   ! Date: 03/21/2016
   548)   ! 
   549)   
   550)   implicit none
   551)   
   552)   class(pm_mphase_type) :: this
   553) 
   554)   character(len=MAXWORDLENGTH) :: word
   555)   PetscInt :: id
   556) 
   557)   id = INPUT_RECORD_UNIT
   558) 
   559)   write(id,'(a29)',advance='no') 'pm: '
   560)   write(id,'(a)') this%name
   561)   write(id,'(a29)',advance='no') 'mode: '
   562)   write(id,'(a)') 'mphase'
   563) 
   564) end subroutine PMMphaseInputRecord
   565) 
   566) ! ************************************************************************** !
   567) 
   568) subroutine PMMphaseDestroy(this)
   569)   ! 
   570)   ! Destroys Mphase process model
   571)   ! 
   572)   ! Author: Glenn Hammond
   573)   ! Date: 03/14/13
   574)   ! 
   575) 
   576)   use Mphase_module, only : MphaseDestroy
   577) 
   578)   implicit none
   579)   
   580)   class(pm_mphase_type) :: this
   581)   
   582)   if (associated(this%next)) then
   583)     call this%next%Destroy()
   584)   endif
   585) 
   586)   ! preserve this ordering
   587)   call MphaseDestroy(this%realization)
   588)   call PMSubsurfaceFlowDestroy(this)
   589)   
   590) end subroutine PMMphaseDestroy
   591)   
   592) end module PM_Mphase_class

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