pm_miscible.F90       coverage:  0.00 %func     0.00 %block


     1) module PM_Miscible_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_miscible_type
    21)   contains
    22)     procedure, public :: Read => PMMiscibleRead
    23)     procedure, public :: InitializeTimestep => PMMiscibleInitializeTimestep
    24)     procedure, public :: Residual => PMMiscibleResidual
    25)     procedure, public :: Jacobian => PMMiscibleJacobian
    26)     procedure, public :: UpdateTimestep => PMMiscibleUpdateTimestep
    27)     procedure, public :: PreSolve => PMMisciblePreSolve
    28)     procedure, public :: PostSolve => PMMisciblePostSolve
    29) #if 0
    30)     procedure, public :: CheckUpdatePre => PMMiscibleCheckUpdatePre
    31)     procedure, public :: CheckUpdatePost => PMMiscibleCheckUpdatePost
    32) #endif
    33)     procedure, public :: TimeCut => PMMiscibleTimeCut
    34)     procedure, public :: UpdateSolution => PMMiscibleUpdateSolution
    35)     procedure, public :: UpdateAuxVars => PMMiscibleUpdateAuxVars
    36)     procedure, public :: MaxChange => PMMiscibleMaxChange
    37)     procedure, public :: ComputeMassBalance => PMMiscibleComputeMassBalance
    38)     procedure, public :: InputRecord => PMMiscibleInputRecord
    39)     procedure, public :: Destroy => PMMiscibleDestroy
    40)   end type pm_miscible_type
    41)   
    42)   public :: PMMiscibleCreate
    43)   
    44) contains
    45) 
    46) ! ************************************************************************** !
    47) 
    48) function PMMiscibleCreate()
    49)   ! 
    50)   ! Creates Miscible process models shell
    51)   ! 
    52)   ! Author: Gautam Bisht
    53)   ! Date: 11/27/13
    54)   ! 
    55) 
    56)   implicit none
    57)   
    58)   class(pm_miscible_type), pointer :: PMMiscibleCreate
    59) 
    60)   class(pm_miscible_type), pointer :: miscible_pm
    61)   
    62) #ifdef PM__DEBUG  
    63)   print *, 'PMMiscibleCreate()'
    64) #endif  
    65) 
    66)   allocate(miscible_pm)
    67) 
    68)   call PMSubsurfaceFlowCreate(miscible_pm)
    69)   miscible_pm%name = 'PMMiscible'
    70) 
    71)   PMMiscibleCreate => miscible_pm
    72)   
    73) end function PMMiscibleCreate
    74) 
    75) ! ************************************************************************** !
    76) 
    77) subroutine PMMiscibleRead(this,input)
    78)   ! 
    79)   ! Reads input file parameters associated with the Miscible process model
    80)   ! 
    81)   ! Author: Glenn Hammond
    82)   ! Date: 01/29/15
    83)   use Input_Aux_module
    84)   use String_module
    85)   use Utility_module
    86)   use EOS_Water_module  
    87)   use Option_module
    88)   use Miscible_Aux_module
    89)  
    90)   implicit none
    91)   
    92)   class(pm_miscible_type) :: this
    93)   type(input_type), pointer :: input
    94)   
    95)   character(len=MAXWORDLENGTH) :: word
    96)   character(len=MAXSTRINGLENGTH) :: error_string
    97)   type(option_type), pointer :: option
    98)   PetscBool :: found
    99) 
   100)   option => this%option
   101)   
   102)   error_string = 'Miscible Options'
   103)   
   104)   input%ierr = 0
   105)   do
   106)   
   107)     call InputReadPflotranString(input,option)
   108)     if (InputError(input)) exit
   109)     if (InputCheckExit(input,option)) exit
   110)     
   111)     call InputReadWord(input,option,word,PETSC_TRUE)
   112)     call InputErrorMsg(input,option,'keyword',error_string)
   113)     call StringToUpper(word)
   114) 
   115)     found = PETSC_FALSE
   116)     call PMSubsurfaceFlowReadSelectCase(this,input,word,found,option)
   117)     if (found) cycle
   118)     
   119)     select case(trim(word))
   120)       case default
   121)         call InputKeywordUnrecognized(word,error_string,option)
   122)     end select
   123)   enddo
   124)   
   125) end subroutine PMMiscibleRead
   126) 
   127) ! ************************************************************************** !
   128) 
   129) subroutine PMMiscibleInitializeTimestep(this)
   130)   ! 
   131)   ! Should not need this as it is called in PreSolve.
   132)   ! 
   133)   ! Author: Gautam Bisht
   134)   ! Date: 11/27/13
   135)   ! 
   136) 
   137)   use Miscible_module, only : MiscibleInitializeTimestep
   138)   
   139)   implicit none
   140)   
   141)   class(pm_miscible_type) :: this
   142) 
   143)   call PMSubsurfaceFlowInitializeTimestepA(this)         
   144) 
   145)   if (this%option%print_screen_flag) then
   146)     write(*,'(/,2("=")," MISCIBLE FLOW ",63("="))')
   147)   endif
   148)   
   149)   call MiscibleInitializeTimestep(this%realization)
   150)   call PMSubsurfaceFlowInitializeTimestepB(this)         
   151)   
   152) end subroutine PMMiscibleInitializeTimestep
   153) 
   154) ! ************************************************************************** !
   155) 
   156) subroutine PMMisciblePreSolve(this)
   157)   ! 
   158)   ! Author: Gautam Bisht
   159)   ! Date: 11/27/13
   160)   ! 
   161) 
   162)   implicit none
   163)   
   164)   class(pm_miscible_type) :: this
   165)   
   166) end subroutine PMMisciblePreSolve
   167) 
   168) ! ************************************************************************** !
   169) 
   170) subroutine PMMisciblePostSolve(this)
   171)   ! 
   172)   ! PMMiscibleUpdatePostSolve:
   173)   ! 
   174)   ! Author: Gautam Bisht
   175)   ! Date: 11/27/13
   176)   ! 
   177) 
   178)   implicit none
   179)   
   180)   class(pm_miscible_type) :: this
   181)   
   182) end subroutine PMMisciblePostSolve
   183) 
   184) ! ************************************************************************** !
   185) 
   186) subroutine PMMiscibleUpdateTimestep(this,dt,dt_min,dt_max,iacceleration, &
   187)                                     num_newton_iterations,tfac)
   188)   ! 
   189)   ! Author: Gautam Bisht
   190)   ! Date: 11/27/13
   191)   ! 
   192) 
   193)   implicit none
   194)   
   195)   class(pm_miscible_type) :: this
   196)   PetscReal :: dt
   197)   PetscReal :: dt_min,dt_max
   198)   PetscInt :: iacceleration
   199)   PetscInt :: num_newton_iterations
   200)   PetscReal :: tfac(:)
   201)   
   202)   PetscReal :: fac
   203)   PetscReal :: ut
   204)   PetscReal :: up
   205)   PetscReal :: utmp
   206)   PetscReal :: uc
   207)   PetscReal :: uus
   208)   PetscReal :: dtt
   209)   PetscReal :: dt_p
   210)   PetscReal :: dt_tfac
   211)   PetscInt :: ifac
   212)   
   213) #ifdef PM_MISCIBLE_DEBUG  
   214)   call printMsg(this%option,'PMMiscible%UpdateTimestep()')
   215) #endif
   216)   
   217)   if (iacceleration > 0) then
   218)     fac = 0.5d0
   219)     if (num_newton_iterations >= iacceleration) then
   220)       fac = 0.33d0
   221)       ut = 0.d0
   222)     else
   223)       up = this%pressure_change_governor/(this%max_pressure_change+0.1)
   224)       utmp = this%temperature_change_governor/(this%max_temperature_change+1.d-5)
   225)       uc = this%xmol_change_governor/(this%max_xmol_change+1.d-6)
   226)       uus= this%saturation_change_governor/(this%max_saturation_change+1.d-6)      
   227)       ut = min(up,utmp,uc,uus)
   228)     endif
   229)     dtt = fac * dt * (1.d0 + ut)
   230)   else
   231)     ifac = max(min(num_newton_iterations,size(tfac)),1)
   232)     dt_tfac = tfac(ifac) * dt
   233) 
   234)     fac = 0.5d0
   235)     up = this%pressure_change_governor/(this%max_pressure_change+0.1)
   236)     dt_p = fac * dt * (1.d0 + up)
   237) 
   238)     dtt = min(dt_tfac,dt_p)
   239)   endif
   240)   
   241)   if (dtt > 2.d0 * dt) dtt = 2.d0 * dt
   242)   if (dtt > dt_max) dtt = dt_max
   243)   ! geh: There used to be code here that cut the time step if it is too
   244)   !      large relative to the simulation time.  This has been removed.
   245)   dtt = max(dtt,dt_min)
   246)   dt = dtt
   247) 
   248)   call PMSubsurfaceFlowLimitDTByCFL(this,dt)
   249)   
   250) end subroutine PMMiscibleUpdateTimestep
   251) 
   252) ! ************************************************************************** !
   253) 
   254) subroutine PMMiscibleResidual(this,snes,xx,r,ierr)
   255)   ! 
   256)   ! Author: Gautam Bisht
   257)   ! Date: 11/27/13
   258)   ! 
   259) 
   260)   use Miscible_module, only : MiscibleResidual
   261) 
   262)   implicit none
   263)   
   264)   class(pm_miscible_type) :: this
   265)   SNES :: snes
   266)   Vec :: xx
   267)   Vec :: r
   268)   PetscErrorCode :: ierr
   269)   
   270)   call MiscibleResidual(snes,xx,r,this%realization,ierr)
   271) 
   272) end subroutine PMMiscibleResidual
   273) 
   274) ! ************************************************************************** !
   275) 
   276) subroutine PMMiscibleJacobian(this,snes,xx,A,B,ierr)
   277)   ! 
   278)   ! Author: Gautam Bisht
   279)   ! Date: 11/27/13
   280)   ! 
   281) 
   282)   use Miscible_module, only : MiscibleJacobian
   283) 
   284)   implicit none
   285)   
   286)   class(pm_miscible_type) :: this
   287)   SNES :: snes
   288)   Vec :: xx
   289)   Mat :: A, B
   290)   PetscErrorCode :: ierr
   291)   
   292)   call MiscibleJacobian(snes,xx,A,B,this%realization,ierr)
   293) 
   294) end subroutine PMMiscibleJacobian
   295)     
   296) #if 0
   297) 
   298) ! ************************************************************************** !
   299) 
   300) subroutine PMMiscibleCheckUpdatePre(this,line_search,X,dX,changed,ierr)
   301)   ! 
   302)   ! Author: Gautam Bisht
   303)   ! Date: 11/27/13
   304)   ! 
   305) 
   306)   use Miscible_module, only : MiscibleCheckUpdatePre
   307) 
   308)   implicit none
   309)   
   310)   class(pm_miscible_type) :: this
   311)   SNESLineSearch :: line_search
   312)   Vec :: X
   313)   Vec :: dX
   314)   PetscBool :: changed
   315)   PetscErrorCode :: ierr
   316)   
   317)   call MiscibleCheckUpdatePre(line_search,X,dX,changed,this%realization,ierr)
   318) 
   319) end subroutine PMMiscibleCheckUpdatePre
   320) 
   321) ! ************************************************************************** !
   322) 
   323) subroutine PMMiscibleCheckUpdatePost(this,line_search,X0,dX,X1,dX_changed, &
   324)                                   X1_changed,ierr)
   325)   ! 
   326)   ! Author: Gautam Bisht
   327)   ! Date: 11/27/13
   328)   ! 
   329) 
   330)   use Miscible_module, only : MiscibleCheckUpdatePost
   331) 
   332)   implicit none
   333)   
   334)   class(pm_miscible_type) :: this
   335)   SNESLineSearch :: line_search
   336)   Vec :: X0
   337)   Vec :: dX
   338)   Vec :: X1
   339)   PetscBool :: dX_changed
   340)   PetscBool :: X1_changed
   341)   PetscErrorCode :: ierr
   342)   
   343)   call MiscibleCheckUpdatePost(line_search,X0,dX,X1,dX_changed, &
   344)                                X1_changed,this%realization,ierr)
   345) 
   346) end subroutine PMMiscibleCheckUpdatePost
   347) #endif
   348) 
   349) ! ************************************************************************** !
   350) 
   351) subroutine PMMiscibleTimeCut(this)
   352)   ! 
   353)   ! Author: Gautam Bisht
   354)   ! Date: 11/27/13
   355)   ! 
   356) 
   357)   use Miscible_module, only : MiscibleTimeCut
   358) 
   359)   implicit none
   360)   
   361)   class(pm_miscible_type) :: this
   362)   
   363)   call PMSubsurfaceFlowTimeCut(this)
   364)   call MiscibleTimeCut(this%realization)
   365) 
   366) end subroutine PMMiscibleTimeCut
   367) 
   368) ! ************************************************************************** !
   369) 
   370) subroutine PMMiscibleUpdateSolution(this)
   371)   ! 
   372)   ! Author: Gautam Bisht
   373)   ! Date: 11/27/13
   374)   ! 
   375) 
   376)   use Miscible_module, only : MiscibleUpdateSolution
   377) 
   378)   implicit none
   379)   
   380)   class(pm_miscible_type) :: this
   381)   
   382)   call PMSubsurfaceFlowUpdateSolution(this)
   383)   call MiscibleUpdateSolution(this%realization)
   384) 
   385) end subroutine PMMiscibleUpdateSolution     
   386) 
   387) ! ************************************************************************** !
   388) 
   389) subroutine PMMiscibleUpdateAuxVars(this)
   390)   ! 
   391)   ! Author: Glenn Hammond
   392)   ! Date: 04/21/14
   393) 
   394)   use Miscible_module, only : MiscibleUpdateAuxVars
   395)   
   396)   implicit none
   397)   
   398)   class(pm_miscible_type) :: this
   399) 
   400)   call MiscibleUpdateAuxVars(this%realization)
   401) 
   402) end subroutine PMMiscibleUpdateAuxVars   
   403) 
   404) ! ************************************************************************** !
   405) 
   406) subroutine PMMiscibleMaxChange(this)
   407)   ! 
   408)   ! Not needed given PMMiscibleMaxChange is called in PostSolve
   409)   ! 
   410)   ! Author: Gautam Bisht
   411)   ! Date: 11/27/13
   412)   ! 
   413) 
   414)   use Mphase_module, only : MphaseMaxChange
   415) 
   416)   implicit none
   417)   
   418)   class(pm_miscible_type) :: this
   419)   
   420)   !geh: yes, call Mphase.  No need to replicate code  
   421)   call MphaseMaxChange(this%realization,this%max_pressure_change, &
   422)                        this%max_temperature_change, &
   423)                        this%max_saturation_change, &
   424)                        this%max_xmol_change)
   425)   if (this%option%print_screen_flag) then
   426)     write(*,'("  --> max chng: dpmx= ",1pe12.4, &
   427)       & " dtmpmx= ",1pe12.4," dcmx= ",1pe12.4," dsmx= ",1pe12.4)') &
   428)           this%max_pressure_change,this%max_temperature_change, &
   429)           this%max_xmol_change,this%max_saturation_change
   430)   endif
   431)   if (this%option%print_file_flag) then
   432)     write(this%option%fid_out,'("  --> max chng: dpmx= ",1pe12.4, &
   433)       & " dtmpmx= ",1pe12.4," dcmx= ",1pe12.4," dsmx= ",1pe12.4)') &
   434)           this%max_pressure_change,this%max_temperature_change, &
   435)           this%max_xmol_change,this%max_saturation_change
   436)   endif     
   437) 
   438) end subroutine PMMiscibleMaxChange
   439) 
   440) ! ************************************************************************** !
   441) 
   442) subroutine PMMiscibleComputeMassBalance(this,mass_balance_array)
   443)   ! 
   444)   ! Author: Gautam Bisht
   445)   ! Date: 11/27/13
   446)   ! 
   447) 
   448)   use Miscible_module, only : MiscibleComputeMassBalance
   449) 
   450)   implicit none
   451)   
   452)   class(pm_miscible_type) :: this
   453)   PetscReal :: mass_balance_array(:)
   454)   
   455)   !geh: currently does not include "trapped" mass
   456)   !call MiscibleComputeMassBalance(this%realization,mass_balance_array)
   457) 
   458) end subroutine PMMiscibleComputeMassBalance
   459) 
   460) ! ************************************************************************** !
   461) 
   462) subroutine PMMiscibleInputRecord(this)
   463)   ! 
   464)   ! Writes ingested information to the input record file.
   465)   ! 
   466)   ! Author: Jenn Frederick, SNL
   467)   ! Date: 03/21/2016
   468)   ! 
   469)   
   470)   implicit none
   471)   
   472)   class(pm_miscible_type) :: this
   473) 
   474)   character(len=MAXWORDLENGTH) :: word
   475)   PetscInt :: id
   476) 
   477)   id = INPUT_RECORD_UNIT
   478) 
   479)   write(id,'(a29)',advance='no') 'pm: '
   480)   write(id,'(a)') this%name
   481)   write(id,'(a29)',advance='no') 'mode: '
   482)   write(id,'(a)') 'miscible'
   483) 
   484) end subroutine PMMiscibleInputRecord
   485) 
   486) ! ************************************************************************** !
   487) 
   488) subroutine PMMiscibleDestroy(this)
   489)   ! 
   490)   ! Destroys Miscible process model
   491)   ! 
   492)   ! Author: Gautam Bisht
   493)   ! Date: 11/27/13
   494)   ! 
   495) 
   496)   use Miscible_module, only : MiscibleDestroy
   497) 
   498)   implicit none
   499)   
   500)   class(pm_miscible_type) :: this
   501)   
   502)   if (associated(this%next)) then
   503)     call this%next%Destroy()
   504)   endif
   505) 
   506)   ! preserve this ordering
   507)   call MiscibleDestroy(this%realization)
   508)   call PMSubsurfaceFlowDestroy(this)
   509)   
   510) end subroutine PMMiscibleDestroy
   511)   
   512) end module PM_Miscible_class

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