pm_th.F90       coverage:  72.22 %func     45.45 %block


     1) module PM_TH_class
     2) 
     3)   use PM_Base_class
     4)   use PM_Subsurface_Flow_class
     5) !geh: using TH_module here fails with gfortran (internal compiler error)
     6) !  use TH_module
     7)   use Realization_Subsurface_class
     8)   use Communicator_Base_module
     9)   use Option_module
    10)   
    11)   use PFLOTRAN_Constants_module
    12) 
    13)   implicit none
    14) 
    15)   private
    16) 
    17) #include "petsc/finclude/petscsys.h"
    18) 
    19) #include "petsc/finclude/petscvec.h"
    20) #include "petsc/finclude/petscvec.h90"
    21) #include "petsc/finclude/petscmat.h"
    22) #include "petsc/finclude/petscmat.h90"
    23) #include "petsc/finclude/petscsnes.h"
    24) 
    25)   type, public, extends(pm_subsurface_flow_type) :: pm_th_type
    26)     class(communicator_type), pointer :: commN
    27)   contains
    28)     procedure, public :: Setup => PMTHSetup
    29)     procedure, public :: Read => PMTHRead
    30)     procedure, public :: InitializeTimestep => PMTHInitializeTimestep
    31)     procedure, public :: Residual => PMTHResidual
    32)     procedure, public :: Jacobian => PMTHJacobian
    33)     procedure, public :: UpdateTimestep => PMTHUpdateTimestep
    34)     procedure, public :: PreSolve => PMTHPreSolve
    35)     procedure, public :: PostSolve => PMTHPostSolve
    36)     procedure, public :: CheckUpdatePre => PMTHCheckUpdatePre
    37)     procedure, public :: CheckUpdatePost => PMTHCheckUpdatePost
    38)     procedure, public :: TimeCut => PMTHTimeCut
    39)     procedure, public :: UpdateSolution => PMTHUpdateSolution
    40)     procedure, public :: UpdateAuxVars => PMTHUpdateAuxVars
    41)     procedure, public :: MaxChange => PMTHMaxChange
    42)     procedure, public :: ComputeMassBalance => PMTHComputeMassBalance
    43)     procedure, public :: InputRecord => PMTHInputRecord
    44)     procedure, public :: Destroy => PMTHDestroy
    45)   end type pm_th_type
    46)   
    47)   public :: PMTHCreate
    48)   
    49) contains
    50) 
    51) ! ************************************************************************** !
    52) 
    53) function PMTHCreate()
    54)   ! 
    55)   ! This routine
    56)   ! 
    57)   ! Author: Gautam Bisht, LBNL
    58)   ! Date: 03/90/13
    59)   ! 
    60) 
    61)   implicit none
    62)   
    63)   class(pm_th_type), pointer :: PMTHCreate
    64) 
    65)   class(pm_th_type), pointer :: th_pm
    66)   
    67) #ifdef PM_TH_DEBUG
    68)   print *, 'PMTHCreate()'
    69) #endif  
    70) 
    71)   allocate(th_pm)
    72) 
    73)   nullify(th_pm%commN)
    74) 
    75)   call PMSubsurfaceFlowCreate(th_pm)
    76)   th_pm%name = 'PMTH'
    77) 
    78)   PMTHCreate => th_pm
    79)   
    80) end function PMTHCreate
    81) 
    82) ! ************************************************************************** !
    83) 
    84) subroutine PMTHRead(this,input)
    85)   ! 
    86)   ! Reads input file parameters associated with the TH process model
    87)   ! 
    88)   ! Author: Glenn Hammond
    89)   ! Date: 01/29/15
    90)   use Input_Aux_module
    91)   use String_module
    92)   use Utility_module
    93)   use EOS_Water_module  
    94)   use Option_module
    95)   use TH_Aux_module
    96)  
    97)   implicit none
    98)   
    99)   class(pm_th_type) :: this
   100)   type(input_type), pointer :: input
   101)   
   102)   character(len=MAXWORDLENGTH) :: word
   103)   character(len=MAXSTRINGLENGTH) :: error_string
   104)   type(option_type), pointer :: option
   105)   PetscBool :: found
   106) 
   107)   option => this%option
   108)   
   109)   error_string = 'TH Options'
   110)   
   111)   input%ierr = 0
   112)   do
   113)   
   114)     call InputReadPflotranString(input,option)
   115)     if (InputError(input)) exit
   116)     if (InputCheckExit(input,option)) exit
   117)     
   118)     call InputReadWord(input,option,word,PETSC_TRUE)
   119)     call InputErrorMsg(input,option,'keyword',error_string)
   120)     call StringToUpper(word)
   121) 
   122)     found = PETSC_FALSE
   123)     call PMSubsurfaceFlowReadSelectCase(this,input,word,found,option)
   124)     if (found) cycle
   125)     
   126)     select case(trim(word))
   127)       case('ITOL_SCALED_RESIDUAL')
   128)         call InputReadDouble(input,option,th_itol_scaled_res)
   129)         call InputDefaultMsg(input,option,'th_itol_scaled_res')
   130)         this%check_post_convergence = PETSC_TRUE
   131)       case('ITOL_RELATIVE_UPDATE')
   132)         call InputReadDouble(input,option,th_itol_rel_update)
   133)         call InputDefaultMsg(input,option,'th_itol_rel_update')
   134)         this%check_post_convergence = PETSC_TRUE        
   135)       case('FREEZING')
   136)         option%use_th_freezing = PETSC_TRUE
   137)         option%io_buffer = ' TH: using FREEZING submode!'
   138)         call printMsg(option)
   139)         ! Override the default setting for TH-mode with freezing
   140)         call EOSWaterSetDensity('PAINTER')
   141)         call EOSWaterSetEnthalpy('PAINTER')
   142)       case('ICE_MODEL')
   143)         call InputReadWord(input,option,word,PETSC_FALSE)
   144)         call StringToUpper(word)
   145)         select case (trim(word))
   146)           case ('PAINTER_EXPLICIT')
   147)             option%ice_model = PAINTER_EXPLICIT
   148)           case ('PAINTER_KARRA_IMPLICIT')
   149)             option%ice_model = PAINTER_KARRA_IMPLICIT
   150)           case ('PAINTER_KARRA_EXPLICIT')
   151)             option%ice_model = PAINTER_KARRA_EXPLICIT
   152)           case ('PAINTER_KARRA_EXPLICIT_NOCRYO')
   153)             option%ice_model = PAINTER_KARRA_EXPLICIT_NOCRYO
   154)           case ('DALL_AMICO')
   155)             option%ice_model = DALL_AMICO
   156)           case default
   157)             option%io_buffer = 'Cannot identify the specificed ice model.' // &
   158)              'Specify PAINTER_EXPLICIT or PAINTER_KARRA_IMPLICIT' // &
   159)              ' or PAINTER_KARRA_EXPLICIT or PAINTER_KARRA_EXPLICIT_NOCRYO ' // &
   160)              ' or DALL_AMICO.'
   161)             call printErrMsg(option)
   162)           end select
   163)       case default
   164)         call InputKeywordUnrecognized(word,error_string,option)
   165)     end select
   166)   enddo
   167)   
   168) end subroutine PMTHRead
   169) 
   170) ! ************************************************************************** !
   171) 
   172) subroutine PMTHSetup(this)
   173)   ! 
   174)   ! This routine
   175)   ! 
   176)   ! Author: Gautam Bisht, LBNL
   177)   ! Date: 03/90/13
   178)   ! 
   179) 
   180)   use Discretization_module
   181)   use Communicator_Structured_class
   182)   use Communicator_Unstructured_class
   183)   use Grid_module 
   184) 
   185)   implicit none
   186)   
   187)   class(pm_th_type) :: this
   188) 
   189)   call PMSubsurfaceFlowSetup(this)
   190)   
   191)   ! set up communicator
   192)   select case(this%realization%discretization%itype)
   193)     case(STRUCTURED_GRID)
   194)       this%commN => StructuredCommunicatorCreate()
   195)     case(UNSTRUCTURED_GRID)
   196)       this%commN => UnstructuredCommunicatorCreate()
   197)   end select
   198)   call this%commN%SetDM(this%realization%discretization%dm_nflowdof)
   199) 
   200) end subroutine PMTHSetup
   201) 
   202) ! ************************************************************************** !
   203) 
   204) subroutine PMTHInitializeTimestep(this)
   205)   ! 
   206)   ! This routine
   207)   ! 
   208)   ! Author: Gautam Bisht, LBNL
   209)   ! Date: 03/90/13
   210)   ! 
   211) 
   212)   use TH_module, only : THInitializeTimestep
   213) 
   214)   implicit none
   215)   
   216)   class(pm_th_type) :: this
   217) 
   218)   call PMSubsurfaceFlowInitializeTimestepA(this)
   219) 
   220)   ! update porosity
   221)   call this%comm1%LocalToLocal(this%realization%field%icap_loc, &
   222)                                this%realization%field%icap_loc)
   223)   call this%comm1%LocalToLocal(this%realization%field%ithrm_loc, &
   224)                                this%realization%field%ithrm_loc)
   225)   call this%comm1%LocalToLocal(this%realization%field%iphas_loc, &
   226)                                this%realization%field%iphas_loc)
   227) 
   228)   if (this%option%print_screen_flag) then
   229)     write(*,'(/,2("=")," TH FLOW ",69("="))')
   230)   endif
   231)   
   232)   call THInitializeTimestep(this%realization)
   233)   call PMSubsurfaceFlowInitializeTimestepB(this)
   234)   
   235) end subroutine PMTHInitializeTimestep
   236) 
   237) ! ************************************************************************** !
   238) 
   239) subroutine PMTHPreSolve(this)
   240)   ! 
   241)   ! This routine
   242)   ! 
   243)   ! Author: Gautam Bisht, LBNL
   244)   ! Date: 03/90/13
   245)   ! 
   246) 
   247)   use Global_module
   248) 
   249)   implicit none
   250)   
   251)   class(pm_th_type) :: this
   252) 
   253) end subroutine PMTHPreSolve
   254) 
   255) ! ************************************************************************** !
   256) 
   257) subroutine PMTHPostSolve(this)
   258)   ! 
   259)   ! Date: 03/14/13
   260)   ! 
   261) 
   262)   use Global_module
   263) 
   264)   implicit none
   265)   
   266)   class(pm_th_type) :: this
   267)   
   268) end subroutine PMTHPostSolve
   269) 
   270) ! ************************************************************************** !
   271) 
   272) subroutine PMTHUpdateTimestep(this,dt,dt_min,dt_max,iacceleration, &
   273)                               num_newton_iterations,tfac)
   274)   ! 
   275)   ! This routine
   276)   ! 
   277)   ! Author: Gautam Bisht, LBNL
   278)   ! Date: 03/90/13
   279)   ! 
   280) 
   281)   implicit none
   282)   
   283)   class(pm_th_type) :: this
   284)   PetscReal :: dt
   285)   PetscReal :: dt_min,dt_max
   286)   PetscInt :: iacceleration
   287)   PetscInt :: num_newton_iterations
   288)   PetscReal :: tfac(:)
   289)   
   290)   PetscReal :: fac
   291)   PetscReal :: ut
   292)   PetscReal :: up
   293)   PetscReal :: utmp
   294)   PetscReal :: dtt
   295)   PetscReal :: dt_u
   296)   PetscReal :: dt_tfac
   297)   PetscInt :: ifac
   298)   
   299) #ifdef PM_TH_DEBUG
   300)   call printMsg(this%option,'PMTH%UpdateTimestep()')
   301) #endif
   302)   
   303)   if (iacceleration > 0) then
   304)     fac = 0.5d0
   305)     if (num_newton_iterations >= iacceleration) then
   306)       fac = 0.33d0
   307)       ut = 0.d0
   308)     else
   309)       up = this%pressure_change_governor/(this%max_pressure_change+0.1)
   310)       utmp = this%temperature_change_governor/ &
   311)              (this%max_temperature_change+1.d-5)
   312)       ut = min(up,utmp)
   313)     endif
   314)     dtt = fac * dt * (1.d0 + ut)
   315)   else
   316)     ifac = max(min(num_newton_iterations,size(tfac)),1)
   317)     dt_tfac = tfac(ifac) * dt
   318) 
   319)     fac = 0.5d0
   320)     up = this%pressure_change_governor/(this%max_pressure_change+0.1)
   321)     utmp = this%temperature_change_governor/ &
   322)            (this%max_temperature_change+1.d-5)
   323)     ut = min(up,utmp)
   324)     dt_u = fac * dt * (1.d0 + ut)
   325) 
   326)     dtt = min(dt_tfac,dt_u)
   327)   endif
   328) 
   329)   if (dtt > 2.d0 * dt) dtt = 2.d0 * dt
   330)   if (dtt > dt_max) dtt = dt_max
   331)   ! geh: There used to be code here that cut the time step if it is too
   332)   !      large relative to the simulation time.  This has been removed.
   333)   dtt = max(dtt,dt_min)
   334)   dt = dtt
   335) 
   336)   call PMSubsurfaceFlowLimitDTByCFL(this,dt)
   337)   
   338) end subroutine PMTHUpdateTimestep
   339) 
   340) ! ************************************************************************** !
   341) 
   342) subroutine PMTHResidual(this,snes,xx,r,ierr)
   343)   ! 
   344)   ! This routine
   345)   ! 
   346)   ! Author: Gautam Bisht, LBNL
   347)   ! Date: 03/90/13
   348)   ! 
   349) 
   350)   use TH_module, only : THResidual
   351) 
   352)   implicit none
   353)   
   354)   class(pm_th_type) :: this
   355)   SNES :: snes
   356)   Vec :: xx
   357)   Vec :: r
   358)   PetscErrorCode :: ierr
   359)   
   360)   call THResidual(snes,xx,r,this%realization,ierr)
   361) 
   362) end subroutine PMTHResidual
   363) 
   364) ! ************************************************************************** !
   365) 
   366) subroutine PMTHJacobian(this,snes,xx,A,B,ierr)
   367)   ! 
   368)   ! This routine
   369)   ! 
   370)   ! Author: Gautam Bisht, LBNL
   371)   ! Date: 03/90/13
   372)   ! 
   373) 
   374)   use TH_module, only : THJacobian
   375) 
   376)   implicit none
   377)   
   378)   class(pm_th_type) :: this
   379)   SNES :: snes
   380)   Vec :: xx
   381)   Mat :: A, B
   382)   PetscErrorCode :: ierr
   383)   
   384)   call THJacobian(snes,xx,A,B,this%realization,ierr)
   385) 
   386) end subroutine PMTHJacobian
   387) 
   388) ! ************************************************************************** !
   389) 
   390) subroutine PMTHCheckUpdatePre(this,line_search,X,dX,changed,ierr)
   391)   ! 
   392)   ! This routine
   393)   ! 
   394)   ! Author: Gautam Bisht, LBNL
   395)   ! Date: 03/90/13
   396)   ! 
   397) 
   398)   use Realization_Subsurface_class
   399)   use Grid_module
   400)   use Field_module
   401)   use Option_module
   402)   use Saturation_Function_module
   403)   use Patch_module
   404)   use TH_Aux_module
   405)   use Global_Aux_module
   406) 
   407)   implicit none
   408)   
   409)   class(pm_th_type) :: this
   410)   SNESLineSearch :: line_search
   411)   Vec :: X
   412)   Vec :: dX
   413)   PetscBool :: changed
   414)   PetscErrorCode :: ierr
   415)   
   416)   PetscReal, pointer :: X_p(:)
   417)   PetscReal, pointer :: dX_p(:)
   418)   PetscReal, pointer :: r_p(:)
   419)   type(grid_type), pointer :: grid
   420)   type(option_type), pointer :: option
   421)   type(patch_type), pointer :: patch
   422)   type(field_type), pointer :: field
   423)   type(TH_auxvar_type), pointer :: TH_auxvars(:)
   424)   type(global_auxvar_type), pointer :: global_auxvars(:)  
   425)   PetscInt :: local_id, ghosted_id
   426)   PetscReal :: P0, P1, P_R, delP, delP_old
   427)   PetscReal :: scale, press_limit, temp_limit
   428)   PetscInt :: iend, istart
   429)   
   430)   patch => this%realization%patch
   431)   grid => patch%grid
   432)   option => this%realization%option
   433)   field => this%realization%field
   434)   TH_auxvars => patch%aux%TH%auxvars
   435)   global_auxvars => patch%aux%Global%auxvars
   436) 
   437)   if (Initialized(this%pressure_change_limit)) then
   438) 
   439)     call VecGetArrayF90(dX,dX_p,ierr);CHKERRQ(ierr)
   440)     call VecGetArrayF90(X,X_p,ierr);CHKERRQ(ierr)
   441) 
   442)     press_limit = dabs(this%pressure_change_limit)
   443)     do local_id = 1, grid%nlmax
   444)       ghosted_id = grid%nL2G(local_id)
   445)       if (patch%imat(ghosted_id) <= 0) cycle
   446)       iend = local_id*option%nflowdof
   447)       istart = iend-option%nflowdof+1
   448)       P0 = X_p(istart)
   449)       delP = dX_p(istart)
   450)       if (press_limit < dabs(delP)) then
   451)         write(option%io_buffer,'("dP_trunc:",1i7,2es15.7)') &         
   452)           grid%nG2A(grid%nL2G(local_id)),press_limit,dabs(delP)
   453)         call printMsgAnyRank(option)
   454)       endif
   455)       delP = sign(min(dabs(delP),press_limit),delP)
   456)       dX_p(istart) = delP
   457)     enddo
   458)     
   459)     call VecRestoreArrayF90(dX,dX_p,ierr);CHKERRQ(ierr)
   460)     call VecRestoreArrayF90(X,X_p,ierr);CHKERRQ(ierr)
   461) 
   462)   endif
   463)   
   464)   if (dabs(this%temperature_change_limit) > 0.d0) then
   465)       
   466)     call VecGetArrayF90(dX,dX_p,ierr);CHKERRQ(ierr)
   467)     call VecGetArrayF90(X,X_p,ierr);CHKERRQ(ierr)
   468) 
   469)     temp_limit = dabs(this%temperature_change_limit)
   470)     do local_id = 1, grid%nlmax
   471)       ghosted_id = grid%nL2G(local_id)
   472)       iend = local_id*option%nflowdof
   473)       istart = iend-option%nflowdof+1
   474)       P0 = X_p(iend)
   475)       delP = dX_p(iend)
   476)       if (abs(delP) > abs(temp_limit)) then
   477)         write(option%io_buffer,'("dT_trunc:",1i7,2es15.7)') &
   478)           grid%nG2A(grid%nL2G(local_id)),temp_limit,dabs(delP)
   479)         call printMsgAnyRank(option)
   480)       endif
   481)       delP = sign(min(dabs(delP),temp_limit),delP)
   482)       dX_p(iend) = delP
   483)     enddo
   484)     
   485)     call VecRestoreArrayF90(dX,dX_p,ierr);CHKERRQ(ierr)
   486)     call VecRestoreArrayF90(X,X_p,ierr);CHKERRQ(ierr)
   487)     
   488)   endif
   489) 
   490) 
   491)   if (Initialized(this%pressure_dampening_factor)) then
   492)     ! P^p+1 = P^p - dP^p
   493)     P_R = option%reference_pressure
   494)     scale = this%pressure_dampening_factor
   495) 
   496)     call VecGetArrayF90(dX,dX_p,ierr);CHKERRQ(ierr)
   497)     call VecGetArrayF90(X,X_p,ierr);CHKERRQ(ierr)
   498)     call VecGetArrayF90(field%flow_r,r_p,ierr);CHKERRQ(ierr)
   499)     do local_id = 1, grid%nlmax
   500)       iend = local_id*option%nflowdof
   501)       istart = iend-option%nflowdof+1
   502)       delP = dX_p(istart)
   503)       P0 = X_p(istart)
   504)       P1 = P0 - delP
   505)       if (P0 < P_R .and. P1 > P_R) then
   506)         write(option%io_buffer,'("U -> S:",1i7,2f12.1)') &
   507)           grid%nG2A(grid%nL2G(local_id)),P0,P1 
   508)         call printMsgAnyRank(option)
   509) #if 0
   510)         ghosted_id = grid%nL2G(local_id)
   511)         call RichardsPrintAuxVars(rich_auxvars(ghosted_id), &
   512)                                   global_auxvars(ghosted_id),ghosted_id)
   513)         write(option%io_buffer,'("Residual:",es15.7)') r_p(istart)
   514)         call printMsgAnyRank(option)
   515) #endif
   516)       else if (P1 < P_R .and. P0 > P_R) then
   517)         write(option%io_buffer,'("S -> U:",1i7,2f12.1)') &
   518)           grid%nG2A(grid%nL2G(local_id)),P0,P1
   519)         call printMsgAnyRank(option)
   520) #if 0
   521)         ghosted_id = grid%nL2G(local_id)
   522)         call RichardsPrintAuxVars(rich_auxvars(ghosted_id), &
   523)                                   global_auxvars(ghosted_id),ghosted_id)
   524)         write(option%io_buffer,'("Residual:",es15.7)') r_p(istart)
   525)         call printMsgAnyRank(option)
   526) #endif
   527)       endif
   528)       ! transition from unsaturated to saturated
   529)       if (P0 < P_R .and. P1 > P_R) then
   530)         dX_p(istart) = scale*delP
   531)       endif
   532)     enddo
   533)     call VecRestoreArrayF90(dX,dX_p,ierr);CHKERRQ(ierr)
   534)     call VecRestoreArrayF90(X,X_p,ierr);CHKERRQ(ierr)
   535)     call VecGetArrayF90(field%flow_r,r_p,ierr);CHKERRQ(ierr)
   536)   endif
   537) 
   538) end subroutine PMTHCheckUpdatePre
   539) 
   540) ! ************************************************************************** !
   541) 
   542) subroutine PMTHCheckUpdatePost(this,line_search,X0,dX,X1,dX_changed, &
   543)                                   X1_changed,ierr)
   544)   ! 
   545)   ! This routine
   546)   ! 
   547)   ! Author: Gautam Bisht, LBNL
   548)   ! Date: 03/90/13
   549)   ! 
   550) 
   551)   use Realization_Subsurface_class
   552)   use Grid_module
   553)   use Field_module
   554)   use Option_module
   555)   use Secondary_Continuum_Aux_module
   556)   use TH_module
   557)   use TH_Aux_module
   558)   use Global_Aux_module
   559)   use Material_Aux_class
   560)   use Patch_module
   561) 
   562)   implicit none
   563)   
   564)   class(pm_th_type) :: this
   565)   SNESLineSearch :: line_search
   566)   Vec :: X0
   567)   Vec :: dX
   568)   Vec :: X1
   569)   PetscBool :: dX_changed
   570)   PetscBool :: X1_changed
   571)   PetscErrorCode :: ierr
   572)   
   573)   PetscReal, pointer :: X1_p(:)
   574)   PetscReal, pointer :: dX_p(:)
   575)   PetscReal, pointer :: ithrm_loc_p(:)
   576)   PetscReal, pointer :: r_p(:)
   577)   type(grid_type), pointer :: grid
   578)   type(option_type), pointer :: option
   579)   type(field_type), pointer :: field
   580)   type(patch_type), pointer :: patch
   581)   type(TH_auxvar_type), pointer :: TH_auxvars(:)
   582)   type(global_auxvar_type), pointer :: global_auxvars(:)  
   583)   type(TH_parameter_type), pointer :: TH_parameter
   584)   type(sec_heat_type), pointer :: TH_sec_heat_vars(:)
   585)   class(material_auxvar_type), pointer :: material_auxvars(:)
   586) 
   587)   PetscInt :: local_id, ghosted_id
   588)   PetscReal :: Res(2)
   589)   PetscReal :: inf_norm, global_inf_norm
   590)   PetscReal :: vol_frac_prim
   591)   PetscInt :: istart, iend
   592)   
   593)   patch => this%realization%patch
   594)   grid => patch%grid
   595)   option => this%realization%option
   596)   field => this%realization%field
   597)   TH_auxvars => patch%aux%TH%auxvars
   598)   TH_parameter => patch%aux%TH%TH_parameter
   599)   global_auxvars => patch%aux%Global%auxvars
   600)   TH_sec_heat_vars => patch%aux%SC_heat%sec_heat_vars
   601)   material_auxvars => patch%aux%Material%auxvars
   602)   
   603)   dX_changed = PETSC_FALSE
   604)   X1_changed = PETSC_FALSE
   605)   
   606)   if (this%check_post_convergence) then
   607)     call VecGetArrayF90(dX,dX_p,ierr);CHKERRQ(ierr)
   608)     call VecGetArrayF90(X1,X1_p,ierr);CHKERRQ(ierr)
   609)     call VecGetArrayF90(field%ithrm_loc,ithrm_loc_p,ierr);CHKERRQ(ierr)
   610)     call VecGetArrayF90(field%flow_r,r_p,ierr);CHKERRQ(ierr)
   611)     
   612)     inf_norm = 0.d0
   613)     vol_frac_prim = 1.d0
   614)     
   615)     do local_id = 1, grid%nlmax
   616)       ghosted_id = grid%nL2G(local_id)
   617)       if (patch%imat(ghosted_id) <= 0) cycle
   618)     
   619)       iend = local_id*option%nflowdof
   620)       istart = iend-option%nflowdof+1
   621)       
   622)       if (option%use_mc) then
   623)         vol_frac_prim = TH_sec_heat_vars(local_id)%epsilon
   624)       endif
   625) 
   626)       call THAccumulation(TH_auxvars(ghosted_id), &
   627)                            global_auxvars(ghosted_id), &
   628)                            material_auxvars(ghosted_id), &
   629)                            TH_parameter%dencpr(int(ithrm_loc_p(ghosted_id))), &
   630)                            option,vol_frac_prim,Res)
   631)                                                         
   632)       inf_norm = max(inf_norm,min(dabs(dX_p(istart)/X1_p(istart)), &
   633)                                   dabs(r_p(istart)/Res(1)), &
   634)                                   dabs(dX_p(iend)/X1_p(iend)), &
   635)                                   dabs(r_p(iend)/Res(2))))
   636)     enddo
   637)     call MPI_Allreduce(inf_norm,global_inf_norm,ONE_INTEGER_MPI, &
   638)                        MPI_DOUBLE_PRECISION, &
   639)                        MPI_MAX,option%mycomm,ierr)
   640)     option%converged = PETSC_TRUE
   641)     if (global_inf_norm > th_itol_scaled_res) &
   642)       option%converged = PETSC_FALSE
   643)     call VecRestoreArrayF90(dX,dX_p,ierr);CHKERRQ(ierr)
   644)     call VecRestoreArrayF90(X1,X1_p,ierr);CHKERRQ(ierr)
   645)     call VecGetArrayF90(field%flow_r,r_p,ierr);CHKERRQ(ierr)
   646)   endif
   647)   
   648) end subroutine PMTHCheckUpdatePost
   649) 
   650) ! ************************************************************************** !
   651) 
   652) subroutine PMTHTimeCut(this)
   653)   ! 
   654)   ! This routine
   655)   ! 
   656)   ! Author: Gautam Bisht, LBNL
   657)   ! Date: 03/90/13
   658)   ! 
   659) 
   660)   use TH_module, only : THTimeCut
   661) 
   662)   implicit none
   663)   
   664)   class(pm_th_type) :: this
   665)   
   666)   call PMSubsurfaceFlowTimeCut(this)
   667)   call THTimeCut(this%realization)
   668) 
   669) end subroutine PMTHTimeCut
   670) 
   671) ! ************************************************************************** !
   672) 
   673) subroutine PMTHUpdateSolution(this)
   674)   ! 
   675)   ! This routine
   676)   ! 
   677)   ! Author: Gautam Bisht, LBNL
   678)   ! Date: 03/90/13
   679)   ! 
   680) 
   681)   use TH_module, only : THUpdateSolution, THUpdateSurfaceBC
   682) 
   683)   implicit none
   684)   
   685)   class(pm_th_type) :: this
   686)   
   687)   call PMSubsurfaceFlowUpdateSolution(this)
   688)   call THUpdateSolution(this%realization)
   689)   if (this%option%surf_flow_on) &
   690)     call THUpdateSurfaceBC(this%realization)
   691) 
   692) end subroutine PMTHUpdateSolution     
   693) 
   694) ! ************************************************************************** !
   695) 
   696) subroutine PMTHUpdateAuxVars(this)
   697)   ! 
   698)   ! Author: Glenn Hammond
   699)   ! Date: 04/21/14
   700) 
   701)   use TH_module, only : THUpdateAuxVars
   702)   
   703)   implicit none
   704)   
   705)   class(pm_th_type) :: this
   706) 
   707)   call THUpdateAuxVars(this%realization)
   708) 
   709) end subroutine PMTHUpdateAuxVars   
   710) 
   711) ! ************************************************************************** !
   712) 
   713) subroutine PMTHMaxChange(this)
   714)   ! 
   715)   ! This routine
   716)   ! 
   717)   ! Author: Gautam Bisht, LBNL
   718)   ! Date: 03/90/13
   719)   ! 
   720) 
   721)   use TH_module, only : THMaxChange
   722) 
   723)   implicit none
   724)   
   725)   class(pm_th_type) :: this
   726)   
   727)   call THMaxChange(this%realization,this%max_pressure_change, &
   728)                    this%max_temperature_change)
   729)     if (this%option%print_screen_flag) then
   730)     write(*,'("  --> max chng: dpmx= ",1pe12.4," dtmpmx= ",1pe12.4)') &
   731)       this%max_pressure_change,this%max_temperature_change
   732)   endif
   733)   if (this%option%print_file_flag) then
   734)     write(this%option%fid_out,'("  --> max chng: dpmx= ",1pe12.4, &
   735)       & " dtmpmx= ",1pe12.4)') &
   736)       this%max_pressure_change,this%max_temperature_change
   737)   endif 
   738) 
   739) end subroutine PMTHMaxChange
   740) 
   741) ! ************************************************************************** !
   742) 
   743) subroutine PMTHComputeMassBalance(this,mass_balance_array)
   744)   ! 
   745)   ! This routine
   746)   ! 
   747)   ! Author: Gautam Bisht, LBNL
   748)   ! Date: 03/90/13
   749)   ! 
   750) 
   751)   use TH_module, only : THComputeMassBalance
   752) 
   753)   implicit none
   754)   
   755)   class(pm_th_type) :: this
   756)   PetscReal :: mass_balance_array(:)
   757)   
   758)   call THComputeMassBalance(this%realization,mass_balance_array)
   759) 
   760) end subroutine PMTHComputeMassBalance
   761) 
   762) ! ************************************************************************** !
   763) 
   764) subroutine PMTHInputRecord(this)
   765)   ! 
   766)   ! Writes ingested information to the input record file.
   767)   ! 
   768)   ! Author: Jenn Frederick, SNL
   769)   ! Date: 03/21/2016
   770)   ! 
   771)   
   772)   implicit none
   773)   
   774)   class(pm_th_type) :: this
   775) 
   776)   character(len=MAXWORDLENGTH) :: word
   777)   PetscInt :: id
   778) 
   779)   id = INPUT_RECORD_UNIT
   780) 
   781)   write(id,'(a29)',advance='no') 'pm: '
   782)   write(id,'(a)') this%name
   783)   write(id,'(a29)',advance='no') 'mode: '
   784)   write(id,'(a)') 'thermo-hydro'
   785) 
   786) end subroutine PMTHInputRecord
   787) 
   788) ! ************************************************************************** !
   789) 
   790) subroutine PMTHDestroy(this)
   791)   ! 
   792)   ! This routine
   793)   ! 
   794)   ! Author: Gautam Bisht, LBNL
   795)   ! Date: 03/90/13
   796)   ! 
   797) 
   798)   use TH_module, only : THDestroy
   799) 
   800)   implicit none
   801)   
   802)   class(pm_th_type) :: this
   803)   
   804)   if (associated(this%next)) then
   805)     call this%next%Destroy()
   806)   endif
   807) 
   808)   call this%commN%Destroy()
   809) 
   810)   ! preserve this ordering
   811)   call THDestroy(this%realization%patch)
   812)   call PMSubsurfaceFlowDestroy(this)
   813) 
   814) end subroutine PMTHDestroy
   815) 
   816) end module PM_TH_class

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