pm_general.F90       coverage:  75.00 %func     77.06 %block


     1) module PM_General_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_general_type
    21)     PetscInt, pointer :: max_change_ivar(:)
    22)     PetscInt, pointer :: max_change_isubvar(:)
    23)   contains
    24)     procedure, public :: Read => PMGeneralRead
    25)     procedure, public :: InitializeRun => PMGeneralInitializeRun
    26)     procedure, public :: InitializeTimestep => PMGeneralInitializeTimestep
    27)     procedure, public :: Residual => PMGeneralResidual
    28)     procedure, public :: Jacobian => PMGeneralJacobian
    29)     procedure, public :: UpdateTimestep => PMGeneralUpdateTimestep
    30)     procedure, public :: PreSolve => PMGeneralPreSolve
    31)     procedure, public :: PostSolve => PMGeneralPostSolve
    32)     procedure, public :: CheckUpdatePre => PMGeneralCheckUpdatePre
    33)     procedure, public :: CheckUpdatePost => PMGeneralCheckUpdatePost
    34)     procedure, public :: TimeCut => PMGeneralTimeCut
    35)     procedure, public :: UpdateSolution => PMGeneralUpdateSolution
    36)     procedure, public :: UpdateAuxVars => PMGeneralUpdateAuxVars
    37)     procedure, public :: MaxChange => PMGeneralMaxChange
    38)     procedure, public :: ComputeMassBalance => PMGeneralComputeMassBalance
    39)     procedure, public :: InputRecord => PMGeneralInputRecord
    40)     procedure, public :: CheckpointBinary => PMGeneralCheckpointBinary
    41)     procedure, public :: RestartBinary => PMGeneralRestartBinary
    42)     procedure, public :: Destroy => PMGeneralDestroy
    43)   end type pm_general_type
    44)   
    45)   public :: PMGeneralCreate
    46)   
    47) contains
    48) 
    49) ! ************************************************************************** !
    50) 
    51) function PMGeneralCreate()
    52)   ! 
    53)   ! Creates General process models shell
    54)   ! 
    55)   ! Author: Glenn Hammond
    56)   ! Date: 03/14/13
    57)   ! 
    58)   use Variables_module, only : LIQUID_PRESSURE, GAS_PRESSURE, AIR_PRESSURE, &
    59)                                LIQUID_MOLE_FRACTION, TEMPERATURE, &
    60)                                GAS_SATURATION
    61)   implicit none
    62)   
    63)   class(pm_general_type), pointer :: PMGeneralCreate
    64) 
    65)   class(pm_general_type), pointer :: general_pm
    66)   
    67) #ifdef PM_GENERAL_DEBUG  
    68)   print *, 'PMGeneralCreate()'
    69) #endif  
    70) 
    71)   allocate(general_pm)
    72)   allocate(general_pm%max_change_ivar(6))
    73)   general_pm%max_change_ivar = [LIQUID_PRESSURE, GAS_PRESSURE, AIR_PRESSURE, &
    74)                                 LIQUID_MOLE_FRACTION, TEMPERATURE, &
    75)                                 GAS_SATURATION]
    76)   allocate(general_pm%max_change_isubvar(6))
    77)                                    ! UNINITIALIZED_INTEGER avoids zeroing of 
    78)                                    ! pressures not represented in phase
    79)                                        ! 2 = air in xmol(air,liquid)
    80)   general_pm%max_change_isubvar = [0,0,0,2,0,0]
    81)   
    82)   call PMSubsurfaceFlowCreate(general_pm)
    83)   general_pm%name = 'PMGeneral'
    84) 
    85)   PMGeneralCreate => general_pm
    86)   
    87) end function PMGeneralCreate
    88) 
    89) ! ************************************************************************** !
    90) 
    91) subroutine PMGeneralRead(this,input)
    92)   ! 
    93)   ! Sets up SNES solvers.
    94)   ! 
    95)   ! Author: Glenn Hammond
    96)   ! Date: 03/04/15
    97)   !
    98)   use General_module
    99)   use General_Aux_module
   100)   use Input_Aux_module
   101)   use String_module
   102)   use Option_module
   103) 
   104)   implicit none
   105)   
   106)   type(input_type), pointer :: input
   107)   
   108)   character(len=MAXWORDLENGTH) :: keyword, word
   109)   class(pm_general_type) :: this
   110)   type(option_type), pointer :: option
   111)   PetscReal :: tempreal
   112)   character(len=MAXSTRINGLENGTH) :: error_string
   113)   PetscBool :: found
   114) 
   115)   option => this%option
   116) 
   117)   error_string = 'General Options'
   118)   
   119)   input%ierr = 0
   120)   do
   121)   
   122)     call InputReadPflotranString(input,option)
   123) 
   124)     if (InputCheckExit(input,option)) exit  
   125) 
   126)     call InputReadWord(input,option,keyword,PETSC_TRUE)
   127)     call InputErrorMsg(input,option,'keyword',error_string)
   128)     call StringToUpper(keyword)
   129)     
   130)     found = PETSC_FALSE
   131)     call PMSubsurfaceFlowReadSelectCase(this,input,keyword,found,option)    
   132)     if (found) cycle
   133)     
   134)     select case(trim(keyword))
   135)       case('ITOL_SCALED_RESIDUAL')
   136)         call InputReadDouble(input,option,general_itol_scaled_res)
   137)         call InputDefaultMsg(input,option,'general_itol_scaled_res')
   138)         this%check_post_convergence = PETSC_TRUE
   139)       case('ITOL_RELATIVE_UPDATE')
   140)         call InputReadDouble(input,option,general_itol_rel_update)
   141)         call InputDefaultMsg(input,option,'general_itol_rel_update')
   142)         this%check_post_convergence = PETSC_TRUE        
   143)       case('TOUGH2_ITOL_SCALED_RESIDUAL')
   144)         ! since general_tough2_itol_scaled_res_e1 is an array, we must read
   145)         ! the tolerance into a double and copy it to the array.
   146)         tempreal = UNINITIALIZED_DOUBLE
   147)         call InputReadDouble(input,option,tempreal)
   148)         ! tempreal will remain uninitialized if the read fails.
   149)         call InputDefaultMsg(input,option,'tough_itol_scaled_residual_e1')
   150)         if (Initialized(tempreal)) then
   151)           general_tough2_itol_scaled_res_e1 = tempreal
   152)         endif
   153)         call InputReadDouble(input,option,general_tough2_itol_scaled_res_e2)
   154)         call InputDefaultMsg(input,option,'tough_itol_scaled_residual_e2')
   155)         general_tough2_conv_criteria = PETSC_TRUE
   156)         this%check_post_convergence = PETSC_TRUE
   157)       case('T2_ITOL_SCALED_RESIDUAL_TEMP')
   158)         call InputReadDouble(input,option,tempreal)
   159)         call InputErrorMsg(input,option, &
   160)                            'tough_itol_scaled_residual_e1 for temperature', &
   161)                            error_string)
   162)         general_tough2_itol_scaled_res_e1(3,:) = tempreal
   163)       case('WINDOW_EPSILON') 
   164)         call InputReadDouble(input,option,window_epsilon)
   165)         call InputErrorMsg(input,option,'window epsilon',error_string)
   166)       case('GAS_COMPONENT_FORMULA_WEIGHT')
   167)         !geh: assuming gas component is index 2
   168)         call InputReadDouble(input,option,fmw_comp(2))
   169)         call InputErrorMsg(input,option,'gas component formula wt.', &
   170)                            error_string)
   171)       case('TWO_PHASE_ENERGY_DOF')
   172)         call InputReadWord(input,option,word,PETSC_TRUE)
   173)         call InputErrorMsg(input,option,'two_phase_energy_dof',error_string)
   174)         call GeneralAuxSetEnergyDOF(word,option)
   175)       case('ISOTHERMAL')
   176)         general_isothermal = PETSC_TRUE
   177)       case('NO_AIR')
   178)         general_no_air = PETSC_TRUE
   179)       case('MAXIMUM_PRESSURE_CHANGE')
   180)         call InputReadDouble(input,option,general_max_pressure_change)
   181)         call InputErrorMsg(input,option,'maximum pressure change', &
   182)                            error_string)
   183)       case('MAX_ITERATION_BEFORE_DAMPING')
   184)         call InputReadInt(input,option,general_max_it_before_damping)
   185)         call InputErrorMsg(input,option,'maximum iteration before damping', &
   186)                            error_string)
   187)       case('DAMPING_FACTOR')
   188)         call InputReadDouble(input,option,general_damping_factor)
   189)         call InputErrorMsg(input,option,'damping factor',error_string)
   190) #if 0        
   191)       case('GOVERN_MAXIMUM_PRESSURE_CHANGE')
   192)         call InputReadDouble(input,option,this%dPmax_allowable)
   193)         call InputErrorMsg(input,option,'maximum allowable pressure change', &
   194)                            error_string)
   195)       case('GOVERN_MAXIMUM_TEMPERATURE_CHANGE')
   196)         call InputReadDouble(input,option,this%dTmax_allowable)
   197)         call InputErrorMsg(input,option, &
   198)                            'maximum allowable temperature change', &
   199)                            error_string)
   200)       case('GOVERN_MAXIMUM_SATURATION_CHANGE')
   201)         call InputReadDouble(input,option,this%dSmax_allowable)
   202)         call InputErrorMsg(input,option,'maximum allowable saturation change', &
   203)                            error_string)
   204)       case('GOVERN_MAXIMUM_MOLE_FRACTION_CHANGE')
   205)         call InputReadDouble(input,option,this%dXmax_allowable)
   206)         call InputErrorMsg(input,option, &
   207)                            'maximum allowable mole fraction change', &
   208)                            error_string)
   209) #endif
   210)       case('DEBUG_CELL')
   211)         call InputReadInt(input,option,general_debug_cell_id)
   212)         call InputErrorMsg(input,option,'debug cell id',error_string)
   213)       case('NO_TEMP_DEPENDENT_DIFFUSION')
   214)         general_temp_dep_gas_air_diff = PETSC_FALSE
   215)       case('DIFFUSE_XMASS')
   216)         general_diffuse_xmol = PETSC_FALSE
   217)       case('HARMONIC_GAS_DIFFUSIVE_DENSITY')
   218)         general_harmonic_diff_density = PETSC_TRUE
   219)       case('ARITHMETIC_GAS_DIFFUSIVE_DENSITY')
   220)         general_harmonic_diff_density = PETSC_FALSE
   221)       case('ANALYTICAL_DERIVATIVES')
   222)         general_analytical_derivatives = PETSC_TRUE
   223)       case default
   224)         call InputKeywordUnrecognized(keyword,'GENERAL Mode',option)
   225)     end select
   226)     
   227)   enddo  
   228)   
   229)   if (general_isothermal .and. &
   230)       general_2ph_energy_dof == GENERAL_AIR_PRESSURE_INDEX) then
   231)     option%io_buffer = 'Isothermal GENERAL mode may only be run with ' // &
   232)                        'temperature as the two phase energy dof.'
   233)     call printErrMsg(option)
   234)   endif
   235) 
   236) end subroutine PMGeneralRead
   237) 
   238) ! ************************************************************************** !
   239) 
   240) recursive subroutine PMGeneralInitializeRun(this)
   241)   ! 
   242)   ! Initializes the time stepping
   243)   ! 
   244)   ! Author: Glenn Hammond
   245)   ! Date: 04/21/14 
   246) 
   247)   use Realization_Base_class
   248)   
   249)   implicit none
   250)   
   251)   class(pm_general_type) :: this
   252)   
   253)   PetscInt :: i
   254)   PetscErrorCode :: ierr
   255) 
   256)   ! need to allocate vectors for max change
   257)   call VecDuplicateVecsF90(this%realization%field%work,SIX_INTEGER, &
   258)                            this%realization%field%max_change_vecs, &
   259)                            ierr);CHKERRQ(ierr)
   260)   ! set initial values
   261)   do i = 1, 6
   262)     call RealizationGetVariable(this%realization, &
   263)                                 this%realization%field%max_change_vecs(i), &
   264)                                 this%max_change_ivar(i), &
   265)                                 this%max_change_isubvar(i))
   266)   enddo
   267) 
   268)   ! call parent implementation
   269)   call PMSubsurfaceFlowInitializeRun(this)
   270) 
   271) end subroutine PMGeneralInitializeRun
   272) 
   273) ! ************************************************************************** !
   274) 
   275) subroutine PMGeneralInitializeTimestep(this)
   276)   ! 
   277)   ! Should not need this as it is called in PreSolve.
   278)   ! 
   279)   ! Author: Glenn Hammond
   280)   ! Date: 03/14/13
   281)   ! 
   282) 
   283)   use General_module, only : GeneralInitializeTimestep
   284)   use Global_module
   285)   use Variables_module, only : TORTUOSITY
   286)   use Material_module, only : MaterialAuxVarCommunicate
   287)   
   288)   implicit none
   289)   
   290)   class(pm_general_type) :: this
   291) 
   292)   call PMSubsurfaceFlowInitializeTimestepA(this)                                 
   293) !geh:remove   everywhere                                
   294)   call MaterialAuxVarCommunicate(this%comm1, &
   295)                                  this%realization%patch%aux%Material, &
   296)                                  this%realization%field%work_loc,TORTUOSITY, &
   297)                                  ZERO_INTEGER)
   298)                                  
   299)   if (this%option%print_screen_flag) then
   300)     write(*,'(/,2("=")," GENERAL FLOW ",64("="))')
   301)   endif
   302)   
   303)   call GeneralInitializeTimestep(this%realization)
   304)   call PMSubsurfaceFlowInitializeTimestepB(this)                                 
   305)   
   306) end subroutine PMGeneralInitializeTimestep
   307) 
   308) ! ************************************************************************** !
   309) 
   310) subroutine PMGeneralPreSolve(this)
   311)   ! 
   312)   ! Author: Glenn Hammond
   313)   ! Date: 03/14/13
   314) 
   315)   implicit none
   316) 
   317)   class(pm_general_type) :: this
   318) 
   319) end subroutine PMGeneralPreSolve
   320) 
   321) ! ************************************************************************** !
   322) 
   323) subroutine PMGeneralPostSolve(this)
   324)   ! 
   325)   ! Author: Glenn Hammond
   326)   ! Date: 03/14/13
   327) 
   328)   implicit none
   329) 
   330)   class(pm_general_type) :: this
   331) 
   332) end subroutine PMGeneralPostSolve
   333) 
   334) ! ************************************************************************** !
   335) 
   336) subroutine PMGeneralUpdateTimestep(this,dt,dt_min,dt_max,iacceleration, &
   337)                                     num_newton_iterations,tfac)
   338)   ! 
   339)   ! Author: Glenn Hammond
   340)   ! Date: 03/14/13
   341)   ! 
   342) 
   343)   use Realization_Base_class, only : RealizationGetVariable
   344)   use Field_module
   345)   use Global_module, only : GlobalSetAuxVarVecLoc
   346)   use Variables_module, only : LIQUID_SATURATION, GAS_SATURATION
   347) 
   348)   implicit none
   349)   
   350)   class(pm_general_type) :: this
   351)   PetscReal :: dt
   352)   PetscReal :: dt_min,dt_max
   353)   PetscInt :: iacceleration
   354)   PetscInt :: num_newton_iterations
   355)   PetscReal :: tfac(:)
   356)   
   357)   PetscReal :: fac
   358)   PetscInt :: ifac
   359)   PetscReal :: up, ut, ux, us, umin
   360)   PetscReal :: dtt
   361)   type(field_type), pointer :: field
   362)   
   363) #ifdef PM_GENERAL_DEBUG  
   364)   call printMsg(this%option,'PMGeneral%UpdateTimestep()')
   365) #endif
   366)   
   367)   fac = 0.5d0
   368)   if (num_newton_iterations >= iacceleration) then
   369)     fac = 0.33d0
   370)     umin = 0.d0
   371)   else
   372)     up = this%pressure_change_governor/(this%max_pressure_change+0.1)
   373)     ut = this%temperature_change_governor/(this%max_temperature_change+1.d-5)
   374)     ux = this%xmol_change_governor/(this%max_xmol_change+1.d-5)
   375)     us = this%saturation_change_governor/(this%max_saturation_change+1.d-5)
   376)     umin = min(up,ut,ux,us)
   377)   endif
   378)   ifac = max(min(num_newton_iterations,size(tfac)),1)
   379)   dtt = fac * dt * (1.d0 + umin)
   380)   dt = min(dtt,tfac(ifac)*dt,dt_max)
   381)   dt = max(dt,dt_min)
   382) 
   383)   if (Initialized(this%cfl_governor)) then
   384)     ! Since saturations are not stored in global_auxvar for general mode, we
   385)     ! must copy them over for the CFL check
   386)     ! liquid saturation
   387)     field => this%realization%field
   388)     call RealizationGetVariable(this%realization,field%work, &
   389)                                 LIQUID_SATURATION,ZERO_INTEGER)
   390)     call this%realization%comm1%GlobalToLocal(field%work,field%work_loc)
   391)     call GlobalSetAuxVarVecLoc(this%realization,field%work_loc, &
   392)                                LIQUID_SATURATION,TIME_NULL)
   393)     call RealizationGetVariable(this%realization,field%work, &
   394)                                 GAS_SATURATION,ZERO_INTEGER)
   395)     call this%realization%comm1%GlobalToLocal(field%work,field%work_loc)
   396)     call GlobalSetAuxVarVecLoc(this%realization,field%work_loc, &
   397)                                GAS_SATURATION,TIME_NULL)
   398)     call PMSubsurfaceFlowLimitDTByCFL(this,dt)
   399)   endif
   400) 
   401) end subroutine PMGeneralUpdateTimestep
   402) 
   403) ! ************************************************************************** !
   404) 
   405) subroutine PMGeneralResidual(this,snes,xx,r,ierr)
   406)   ! 
   407)   ! Author: Glenn Hammond
   408)   ! Date: 03/14/13
   409)   ! 
   410) 
   411)   use General_module, only : GeneralResidual
   412) 
   413)   implicit none
   414)   
   415)   class(pm_general_type) :: this
   416)   SNES :: snes
   417)   Vec :: xx
   418)   Vec :: r
   419)   PetscErrorCode :: ierr
   420)   
   421)   call PMSubsurfaceFlowUpdatePropertiesNI(this)
   422)   call GeneralResidual(snes,xx,r,this%realization,ierr)
   423) 
   424) end subroutine PMGeneralResidual
   425) 
   426) ! ************************************************************************** !
   427) 
   428) subroutine PMGeneralJacobian(this,snes,xx,A,B,ierr)
   429)   ! 
   430)   ! Author: Glenn Hammond
   431)   ! Date: 03/14/13
   432)   ! 
   433) 
   434)   use General_module, only : GeneralJacobian
   435) 
   436)   implicit none
   437)   
   438)   class(pm_general_type) :: this
   439)   SNES :: snes
   440)   Vec :: xx
   441)   Mat :: A, B
   442)   PetscErrorCode :: ierr
   443)   
   444)   call GeneralJacobian(snes,xx,A,B,this%realization,ierr)
   445) 
   446) end subroutine PMGeneralJacobian
   447) 
   448) ! ************************************************************************** !
   449) 
   450) subroutine PMGeneralCheckUpdatePre(this,line_search,X,dX,changed,ierr)
   451)   ! 
   452)   ! Author: Glenn Hammond
   453)   ! Date: 03/14/13
   454)   ! 
   455)   use Realization_Subsurface_class
   456)   use Grid_module
   457)   use Field_module
   458)   use Option_module
   459)   use Saturation_Function_module
   460)   use Patch_module
   461)   use General_Aux_module
   462)   use Global_Aux_module
   463)   
   464)   implicit none
   465)   
   466)   class(pm_general_type) :: this
   467)   SNESLineSearch :: line_search
   468)   Vec :: X
   469)   Vec :: dX
   470)   PetscBool :: changed
   471)   PetscErrorCode :: ierr
   472)   
   473)   PetscReal, pointer :: X_p(:)
   474)   PetscReal, pointer :: dX_p(:)
   475)   PetscReal, pointer :: r_p(:)
   476)   type(grid_type), pointer :: grid
   477)   type(option_type), pointer :: option
   478)   type(patch_type), pointer :: patch
   479)   type(field_type), pointer :: field
   480)   type(general_auxvar_type), pointer :: gen_auxvars(:,:)
   481)   type(global_auxvar_type), pointer :: global_auxvars(:)  
   482) #ifdef DEBUG_GENERAL_INFO
   483)   character(len=MAXSTRINGLENGTH) :: string, string2, string3
   484)   character(len=MAXWORDLENGTH) :: cell_id_word
   485)   PetscInt, parameter :: max_cell_id = 10
   486)   PetscInt :: cell_id, cell_locator(0:max_cell_id)
   487)   PetscInt :: i, ii
   488) #endif
   489)   PetscInt :: local_id, ghosted_id
   490)   PetscInt :: offset
   491)   PetscInt :: liquid_pressure_index, gas_pressure_index, air_pressure_index
   492)   PetscInt :: saturation_index, temperature_index, xmol_index
   493)   PetscInt :: lid, gid, apid, cpid, vpid, spid
   494)   PetscInt :: pgas_index
   495)   PetscReal :: liquid_pressure0, liquid_pressure1, del_liquid_pressure
   496)   PetscReal :: gas_pressure0, gas_pressure1, del_gas_pressure
   497)   PetscReal :: air_pressure0, air_pressure1, del_air_pressure
   498)   PetscReal :: temperature0, temperature1, del_temperature
   499)   PetscReal :: saturation0, saturation1, del_saturation
   500)   PetscReal :: xmol0, xmol1, del_xmol
   501)   PetscReal :: max_saturation_change = 0.125d0
   502)   PetscReal :: max_temperature_change = 10.d0
   503)   PetscReal :: min_pressure
   504)   PetscReal :: scale, temp_scale, temp_real
   505)   PetscReal, parameter :: tolerance = 0.99d0
   506)   PetscReal, parameter :: initial_scale = 1.d0
   507)   SNES :: snes
   508)   PetscInt :: newton_iteration
   509)   
   510)   grid => this%realization%patch%grid
   511)   option => this%realization%option
   512)   field => this%realization%field
   513)   gen_auxvars => this%realization%patch%aux%General%auxvars
   514)   global_auxvars => this%realization%patch%aux%Global%auxvars
   515) 
   516)   patch => this%realization%patch
   517) 
   518)   spid = option%saturation_pressure_id
   519)   apid = option%air_pressure_id
   520) 
   521)   call SNESLineSearchGetSNES(line_search,snes,ierr)
   522)   call SNESGetIterationNumber(snes,newton_iteration,ierr)
   523) 
   524)   call VecGetArrayF90(dX,dX_p,ierr);CHKERRQ(ierr)
   525)   call VecGetArrayReadF90(X,X_p,ierr);CHKERRQ(ierr)
   526) 
   527)   changed = PETSC_TRUE
   528) 
   529)   ! truncation
   530)   ! air mole fraction in liquid phase must be truncated.  we do not use scaling
   531)   ! here because of the very small values.  just truncation.
   532)   do local_id = 1, grid%nlmax
   533)     ghosted_id = grid%nL2G(local_id)
   534)     if (patch%imat(ghosted_id) <= 0) cycle
   535)     offset = (local_id-1)*option%nflowdof
   536)     select case(global_auxvars(ghosted_id)%istate)
   537)       case(LIQUID_STATE)
   538)         xmol_index = offset + GENERAL_LIQUID_STATE_X_MOLE_DOF
   539)         if (X_p(xmol_index) - dX_p(xmol_index) < 0.d0) then
   540)           ! we use 1.d-10 since cancelation can occur with smaller values
   541)           dX_p(xmol_index) = X_p(xmol_index)
   542)         endif
   543)       case(TWO_PHASE_STATE)
   544)         pgas_index = offset + GENERAL_GAS_PRESSURE_DOF
   545)         if (X_p(pgas_index) - dX_p(pgas_index) < &
   546)             gen_auxvars(ZERO_INTEGER,ghosted_id)% &
   547)               pres(option%saturation_pressure_id)) then
   548)           dX_p(pgas_index) = X_p(pgas_index) - &
   549)             gen_auxvars(ZERO_INTEGER,ghosted_id)% &
   550)               pres(option%saturation_pressure_id)
   551)         endif
   552)     end select
   553)   enddo
   554) 
   555)   scale = initial_scale
   556)   if (general_max_it_before_damping > 0 .and. &
   557)       newton_iteration > general_max_it_before_damping) then
   558)     scale = general_damping_factor
   559)   endif
   560) 
   561) #define LIMIT_MAX_PRESSURE_CHANGE
   562) #define LIMIT_MAX_SATURATION_CHANGE
   563) !!#define LIMIT_MAX_TEMPERATURE_CHANGE
   564) !#define TRUNCATE_LIQUID_PRESSURE
   565) !! TRUNCATE_GAS/AIR_PRESSURE is needed for times when the solve wants
   566) !! to pull them negative.
   567) !#define TRUNCATE_GAS_PRESSURE
   568) !#define TRUNCATE_AIR_PRESSURE
   569) 
   570) #ifdef DEBUG_GENERAL_INFO
   571)   cell_locator = 0
   572) #endif
   573) 
   574)   ! scaling
   575)   do local_id = 1, grid%nlmax
   576)     ghosted_id = grid%nL2G(local_id)
   577)     offset = (local_id-1)*option%nflowdof
   578)     temp_scale = 1.d0
   579) #ifdef DEBUG_GENERAL_INFO
   580)     cell_id = grid%nG2A(ghosted_id)
   581)     write(cell_id_word,*) cell_id
   582)     cell_id_word = '(Cell ' // trim(adjustl(cell_id_word)) // '): '
   583) #endif
   584)     select case(global_auxvars(ghosted_id)%istate)
   585)       case(LIQUID_STATE)
   586)         liquid_pressure_index  = offset + GENERAL_LIQUID_PRESSURE_DOF
   587)         temperature_index  = offset + GENERAL_ENERGY_DOF
   588)         dX_p(liquid_pressure_index) = dX_p(liquid_pressure_index) * &
   589)                                       general_pressure_scale
   590)         temp_scale = 1.d0
   591)         del_liquid_pressure = dX_p(liquid_pressure_index)
   592)         liquid_pressure0 = X_p(liquid_pressure_index)
   593)         liquid_pressure1 = liquid_pressure0 - del_liquid_pressure
   594)         del_temperature = dX_p(temperature_index)
   595)         temperature0 = X_p(temperature_index)
   596)         temperature1 = temperature0 - del_temperature
   597) #ifdef LIMIT_MAX_PRESSURE_CHANGE
   598)         if (dabs(del_liquid_pressure) > general_max_pressure_change) then
   599)           temp_real = dabs(general_max_pressure_change/del_liquid_pressure)
   600) #ifdef DEBUG_GENERAL_INFO
   601)           if (cell_locator(0) < max_cell_id) then
   602)             cell_locator(0) = cell_locator(0) + 1
   603)             cell_locator(cell_locator(0)) = ghosted_id
   604)           endif
   605)           string = trim(cell_id_word) // &
   606)             'Liquid pressure change scaled to truncate at max_pressure_change: '
   607)           call printMsg(option,string)
   608)           write(string2,*) liquid_pressure0
   609)           string = '  Liquid Pressure 0    : ' // adjustl(string2)
   610)           call printMsg(option,string)
   611)           write(string2,*) liquid_pressure1
   612)           string = '  Liquid Pressure 1    : ' // adjustl(string2)
   613)           call printMsg(option,string)
   614)           write(string2,*) -1.d0*del_liquid_pressure
   615)           string = 'Liquid Pressure change : ' // adjustl(string2)
   616)           call printMsg(option,string)
   617)           write(string2,*) temp_real
   618)           string = '          scaling  : ' // adjustl(string2)
   619)           call printMsg(option,string)
   620) #endif
   621)           temp_scale = min(temp_scale,temp_real)
   622)         endif
   623) #endif 
   624) !LIMIT_MAX_PRESSURE_CHANGE
   625) #ifdef TRUNCATE_LIQUID_PRESSURE
   626)         ! truncate liquid pressure change to prevent liquid pressure from 
   627)         ! dropping below the air pressure while in the liquid state
   628)         min_pressure = gen_auxvars(ZERO_INTEGER,ghosted_id)%pres(apid) + &
   629)                        gen_auxvars(ZERO_INTEGER,ghosted_id)%pres(spid)
   630)         if (liquid_pressure1 < min_pressure) then
   631)           temp_real = tolerance * (liquid_pressure0 - min_pressure)
   632)           if (dabs(del_liquid_pressure) > 1.d-40) then
   633)             temp_real = dabs(temp_real / del_liquid_pressure)
   634)           else
   635)             option%io_buffer = 'Something is seriously wrong in ' // &
   636)               'GeneralCheckUpdatePre(liquid<min).  Contact Glenn through ' // &
   637)               'pflotran-dev@googlegroups.com.'
   638)             call printErrMsg(option)
   639)           endif
   640) #ifdef DEBUG_GENERAL_INFO
   641)           if (cell_locator(0) < max_cell_id) then
   642)             cell_locator(0) = cell_locator(0) + 1
   643)             cell_locator(cell_locator(0)) = ghosted_id
   644)           endif
   645)           string = trim(cell_id_word) // &
   646)             'Liquid pressure change scaled to prevent liquid ' // &
   647)             'pressure from dropping below air pressure: '
   648)           call printMsg(option,string)
   649)           write(string2,*) liquid_pressure0
   650)           string = '  Liquid pressure 0: ' // adjustl(string2)
   651)           call printMsg(option,string)
   652)           write(string2,*) liquid_pressure1
   653)           string = '  Liquid pressure 1: ' // adjustl(string2)
   654)           call printMsg(option,string)
   655)           write(string2,*) -1.d0*del_liquid_pressure
   656)           string = '  pressure change  : ' // adjustl(string2)
   657)           call printMsg(option,string)
   658)           write(string2,*) temp_real
   659)           string = '          scaling  : ' // adjustl(string2)
   660)           call printMsg(option,string)
   661) #endif
   662)           temp_scale = min(temp_scale,temp_real)
   663)         endif
   664) #endif 
   665) !TRUNCATE_LIQUID_PRESSURE  
   666) #ifdef LIMIT_MAX_TEMPERATURE_CHANGE
   667)         if (dabs(del_temperature) > max_temperature_change) then
   668)           temp_real = dabs(max_temperature_change/del_temperature)
   669) #ifdef DEBUG_GENERAL_INFO
   670)           if (cell_locator(0) < max_cell_id) then
   671)             cell_locator(0) = cell_locator(0) + 1
   672)             cell_locator(cell_locator(0)) = ghosted_id
   673)           endif
   674)           string = trim(cell_id_word) // &
   675)             'Temperature change scaled to truncate at max_temperature_change: '
   676)           call printMsg(option,string)
   677)           write(string2,*) temperature0
   678)           string = '  Temperature 0    : ' // adjustl(string2)
   679)           call printMsg(option,string)
   680)           write(string2,*) temperature1
   681)           string = '  Temperature 1    : ' // adjustl(string2)
   682)           call printMsg(option,string)
   683)           write(string2,*) -1.d0*del_temperature
   684)           string = 'Temperature change : ' // adjustl(string2)
   685)           call printMsg(option,string)
   686)           write(string2,*) temp_real
   687)           string = '          scaling  : ' // adjustl(string2)
   688)           call printMsg(option,string)
   689) #endif
   690)           temp_scale = min(temp_scale,temp_real)
   691)         endif
   692) #endif 
   693) !LIMIT_MAX_TEMPERATURE_CHANGE        
   694)       case(TWO_PHASE_STATE)
   695)         gas_pressure_index = offset + GENERAL_GAS_PRESSURE_DOF
   696) !        air_pressure_index = offset + 2
   697)         saturation_index = offset + GENERAL_GAS_SATURATION_DOF
   698)         temperature_index  = offset + GENERAL_ENERGY_DOF
   699)         dX_p(gas_pressure_index) = dX_p(gas_pressure_index) * &
   700)                                    general_pressure_scale
   701)         if (general_2ph_energy_dof == GENERAL_AIR_PRESSURE_INDEX) then                                   
   702)           air_pressure_index = offset + GENERAL_ENERGY_DOF
   703)           dX_p(air_pressure_index) = dX_p(air_pressure_index) * &
   704)                                      general_pressure_scale
   705)           del_air_pressure = dX_p(air_pressure_index)
   706)           air_pressure0 = X_p(air_pressure_index)
   707)           air_pressure1 = air_pressure0 - del_air_pressure
   708)         endif
   709)         temp_scale = 1.d0
   710)         del_gas_pressure = dX_p(gas_pressure_index)
   711)         gas_pressure0 = X_p(gas_pressure_index)
   712)         gas_pressure1 = gas_pressure0 - del_gas_pressure
   713)         del_saturation = dX_p(saturation_index)
   714)         saturation0 = X_p(saturation_index)
   715)         saturation1 = saturation0 - del_saturation
   716) #ifdef LIMIT_MAX_PRESSURE_CHANGE
   717)         if (dabs(del_gas_pressure) > general_max_pressure_change) then
   718)           temp_real = dabs(general_max_pressure_change/del_gas_pressure)
   719) #ifdef DEBUG_GENERAL_INFO
   720)           if (cell_locator(0) < max_cell_id) then
   721)             cell_locator(0) = cell_locator(0) + 1
   722)             cell_locator(cell_locator(0)) = ghosted_id
   723)           endif
   724)           string = trim(cell_id_word) // &
   725)             'Gas pressure change scaled to truncate at max_pressure_change: '
   726)           call printMsg(option,string)
   727)           write(string2,*) gas_pressure0
   728)           string = '  Gas Pressure 0    : ' // adjustl(string2)
   729)           call printMsg(option,string)
   730)           write(string2,*) gas_pressure1
   731)           string = '  Gas Pressure 1    : ' // adjustl(string2)
   732)           call printMsg(option,string)
   733)           write(string2,*) -1.d0*del_gas_pressure
   734)           string = 'Gas Pressure change : ' // adjustl(string2)
   735)           call printMsg(option,string)
   736)           write(string2,*) temp_real
   737)           string = '          scaling  : ' // adjustl(string2)
   738)           call printMsg(option,string)
   739) #endif
   740)           temp_scale = min(temp_scale,temp_real)
   741)         endif
   742) #endif
   743) #ifdef TRUNCATE_GAS_PRESSURE
   744)         if (gas_pressure1 <= 0.d0) then
   745)           if (dabs(del_gas_pressure) > 1.d-40) then
   746)             temp_real = tolerance * dabs(gas_pressure0 / del_gas_pressure)
   747) #ifdef DEBUG_GENERAL_INFO
   748)             if (cell_locator(0) < max_cell_id) then
   749)               cell_locator(0) = cell_locator(0) + 1
   750)               cell_locator(cell_locator(0)) = ghosted_id
   751)             endif
   752)             string = trim(cell_id_word) // &
   753)               'Gas pressure change scaled to prevent gas ' // &
   754)               'pressure from dropping below zero: '
   755)             call printMsg(option,string)
   756)             write(string2,*) gas_pressure0
   757)             string = '  Gas pressure 0   : ' // adjustl(string2)
   758)             call printMsg(option,string)
   759)             write(string2,*) gas_pressure1
   760)             string = '  Gas pressure 1   : ' // adjustl(string2)
   761)             call printMsg(option,string)
   762)             write(string2,*) -1.d0*del_gas_pressure
   763)             string = '  pressure change  : ' // adjustl(string2)
   764)             call printMsg(option,string)
   765)             write(string2,*) temp_real
   766)             string = '          scaling  : ' // adjustl(string2)
   767)             call printMsg(option,string)
   768) #endif
   769)             temp_scale = min(temp_scale,temp_real)
   770)           endif
   771)         endif
   772) #endif 
   773) !TRUNCATE_GAS_PRESSURE
   774) #ifdef TRUNCATE_AIR_PRESSURE
   775)         if (air_pressure1 <= 0.d0) then
   776)           if (dabs(del_air_pressure) > 1.d-40) then
   777)             temp_real = tolerance * dabs(air_pressure0 / del_air_pressure)
   778) #ifdef DEBUG_GENERAL_INFO
   779)             if (cell_locator(0) < max_cell_id) then
   780)               cell_locator(0) = cell_locator(0) + 1
   781)               cell_locator(cell_locator(0)) = ghosted_id
   782)             endif
   783)             string = trim(cell_id_word) // &
   784)               'Air pressure change scaled to prevent air ' // &
   785)               'pressure from dropping below zero: '
   786)             call printMsg(option,string)
   787)             write(string2,*) air_pressure0
   788)             string = '  Air pressure 0   : ' // adjustl(string2)
   789)             call printMsg(option,string)
   790)             write(string2,*) air_pressure1
   791)             string = '  Air pressure 1   : ' // adjustl(string2)
   792)             call printMsg(option,string)
   793)             write(string2,*) -1.d0*del_air_pressure
   794)             string = '  pressure change  : ' // adjustl(string2)
   795)             call printMsg(option,string)
   796)             write(string2,*) temp_real
   797)             string = '          scaling  : ' // adjustl(string2)
   798)             call printMsg(option,string)
   799) #endif
   800)             temp_scale = min(temp_scale,temp_real)
   801)           endif
   802)         endif
   803) #endif 
   804) !TRUNCATE_AIR_PRESSURE
   805) #if defined(TRUNCATE_GAS_PRESSURE) && defined(TRUNCATE_AIR_PRESSURE)
   806)         ! have to factor in scaled update from previous conditionals
   807)         gas_pressure1 = gas_pressure0 - temp_scale * del_gas_pressure
   808)         air_pressure1 = air_pressure0 - temp_scale * del_air_pressure
   809)         if (gas_pressure1 <= air_pressure1) then
   810) !          temp_real = (air_pressure0 - gas_pressure0) / &
   811) !                      (temp_scale * (del_air_pressure - del_gas_pressure))
   812) !          temp_real = temp_real * tolerance * temp_scale
   813)           ! refactored to prevent divide by 0
   814)           temp_real = temp_scale * (del_air_pressure - del_gas_pressure)
   815)           if (temp_real > 0.d0) then
   816)             temp_real = (air_pressure0 - gas_pressure0) / temp_real
   817)             temp_real = temp_real * tolerance * temp_scale
   818)           endif
   819)           if (temp_real <= 0.d0) then
   820)             ! add info on pressures here
   821)             option%io_buffer = 'Something is seriously wrong in ' // &
   822)               'GeneralCheckUpdatePre(gas<=air).  Contact Glenn through ' // &
   823)               'pflotran-dev@googlegroups.com.'
   824)             call printErrMsg(option)
   825)           endif
   826) #ifdef DEBUG_GENERAL_INFO
   827)           if (cell_locator(0) < max_cell_id) then
   828)             cell_locator(0) = cell_locator(0) + 1
   829)             cell_locator(cell_locator(0)) = ghosted_id
   830)           endif
   831)           string = trim(cell_id_word) // &
   832)             'Gas/Air pressure change scaled again to prevent gas ' // &
   833)             'pressure from dropping below air pressure: '
   834)           call printMsg(option,string)
   835)           write(string2,*) gas_pressure0
   836)           string = '  Gas pressure 0       : ' // adjustl(string2)
   837)           call printMsg(option,string)
   838)           write(string2,*) gas_pressure1
   839)           string = '  Gas pressure 1       : ' // adjustl(string2)
   840)           call printMsg(option,string)
   841)           write(string2,*) -1.d0*temp_real*del_gas_pressure
   842)           string = '  Gas pressure change  : ' // adjustl(string2)
   843)           call printMsg(option,string)
   844)           write(string2,*) air_pressure0
   845)           string = '  Air pressure 0       : ' // adjustl(string2)
   846)           call printMsg(option,string)
   847)           write(string2,*) air_pressure1
   848)           string = '  Air pressure 1       : ' // adjustl(string2)
   849)           call printMsg(option,string)
   850)           write(string2,*) -1.d0*temp_real*del_air_pressure
   851)           string = '  Air pressure change  : ' // adjustl(string2)
   852)           write(string2,*) temp_real
   853)           string = '          scaling  : ' // adjustl(string2)
   854)           call printMsg(option,string)
   855) #endif
   856)           temp_scale = min(temp_scale,temp_real)
   857)         endif
   858) #endif 
   859) !TRUNCATE_GAS_PRESSURE && TRUNCATE_AIR_PRESSURE
   860) #ifdef LIMIT_MAX_SATURATION_CHANGE
   861)         if (dabs(del_saturation) > max_saturation_change) then
   862)           temp_real = dabs(max_saturation_change/del_saturation)
   863) #ifdef DEBUG_GENERAL_INFO
   864)           if (cell_locator(0) < max_cell_id) then
   865)             cell_locator(0) = cell_locator(0) + 1
   866)             cell_locator(cell_locator(0)) = ghosted_id
   867)           endif
   868)           string = trim(cell_id_word) // &
   869)             'Gas saturation change scaled to truncate at ' // &
   870)             'max_saturation_change: '
   871)           call printMsg(option,string)
   872)           write(string2,*) saturation0
   873)           string = '  Saturation 0    : ' // adjustl(string2)
   874)           call printMsg(option,string)
   875)           write(string2,*) saturation1
   876)           string = '  Saturation 1    : ' // adjustl(string2)
   877)           call printMsg(option,string)
   878)           write(string2,*) -1.d0*del_saturation
   879)           string = 'Saturation change : ' // adjustl(string2)
   880)           call printMsg(option,string)
   881)           write(string2,*) temp_real
   882)           string = '          scaling  : ' // adjustl(string2)
   883)           call printMsg(option,string)
   884) #endif
   885)           temp_scale = min(temp_scale,temp_real)
   886)         endif
   887) #endif 
   888) !LIMIT_MAX_SATURATION_CHANGE        
   889)       case(GAS_STATE) 
   890)         gas_pressure_index = offset + GENERAL_GAS_PRESSURE_DOF
   891)         air_pressure_index = offset + GENERAL_GAS_STATE_AIR_PRESSURE_DOF
   892)         dX_p(gas_pressure_index) = dX_p(gas_pressure_index) * &
   893)                                    general_pressure_scale
   894)         dX_p(air_pressure_index) = dX_p(air_pressure_index) * &
   895)                                    general_pressure_scale
   896)     end select
   897)     scale = min(scale,temp_scale) 
   898)   enddo
   899) 
   900)   temp_scale = scale
   901)   call MPI_Allreduce(temp_scale,scale,ONE_INTEGER_MPI, &
   902)                      MPI_DOUBLE_PRECISION, &
   903)                      MPI_MIN,option%mycomm,ierr)
   904) 
   905)   if (scale < 0.9999d0) then
   906) #ifdef DEBUG_GENERAL_INFO
   907)     string  = '++++++++++++++++++++++++++++++++++++++++++++++++++++++'
   908)     call printMsg(option,string)
   909)     write(string2,*) scale, (grid%nG2A(cell_locator(i)),i=1,cell_locator(0))
   910)     string = 'Final scaling: : ' // adjustl(string2)
   911)     call printMsg(option,string)
   912)     do i = 1, cell_locator(0)
   913)       ghosted_id = cell_locator(i)
   914)       offset = (ghosted_id-1)*option%nflowdof
   915)       write(string2,*) grid%nG2A(ghosted_id)
   916)       string = 'Cell ' // trim(adjustl(string2))
   917)       write(string2,*) global_auxvars(ghosted_id)%istate
   918)       string = trim(string) // ' (State = ' // trim(adjustl(string2)) // ') '
   919)       call printMsg(option,string)
   920)       ! for some reason cannot perform array operation on dX_p(:)
   921)       write(string2,*) (X_p(offset+ii),ii=1,3)
   922)       string = '   Orig. Solution: ' // trim(adjustl(string2))
   923)       call printMsg(option,string)
   924)       write(string2,*) (X_p(offset+ii)-dX_p(offset+ii),ii=1,3)
   925)       string = '  Solution before: ' // trim(adjustl(string2))
   926)       call printMsg(option,string)
   927)       write(string2,*) (X_p(offset+ii)-scale*dX_p(offset+ii),ii=1,3)
   928)       string = '   Solution after: ' // trim(adjustl(string2))
   929)       call printMsg(option,string)
   930)     enddo
   931)     string  = '++++++++++++++++++++++++++++++++++++++++++++++++++++++'
   932)     call printMsg(option,string)
   933) #endif
   934)     dX_p = scale*dX_p
   935)   endif
   936) 
   937)   call VecRestoreArrayF90(dX,dX_p,ierr);CHKERRQ(ierr)
   938)   call VecRestoreArrayReadF90(X,X_p,ierr);CHKERRQ(ierr)
   939) 
   940) end subroutine PMGeneralCheckUpdatePre
   941) 
   942) ! ************************************************************************** !
   943) 
   944) subroutine PMGeneralCheckUpdatePost(this,line_search,X0,dX,X1,dX_changed, &
   945)                                     X1_changed,ierr)
   946)   ! 
   947)   ! Author: Glenn Hammond
   948)   ! Date: 03/14/13
   949)   ! 
   950)   use General_Aux_module
   951)   use Global_Aux_module
   952)   use Grid_module
   953)   use Option_module
   954)   use Realization_Subsurface_class
   955)   use Grid_module
   956)   use Field_module
   957)   use Patch_module
   958)   use Option_module
   959)   use Material_Aux_class  
   960)   use Output_EKG_module
   961)   
   962)   implicit none
   963)   
   964)   class(pm_general_type) :: this
   965)   SNESLineSearch :: line_search
   966)   Vec :: X0
   967)   Vec :: dX
   968)   Vec :: X1
   969)   PetscBool :: dX_changed
   970)   PetscBool :: X1_changed
   971)   PetscErrorCode :: ierr
   972) 
   973)   PetscReal, pointer :: X0_p(:)
   974)   PetscReal, pointer :: X1_p(:)
   975)   PetscReal, pointer :: dX_p(:)
   976)   PetscReal, pointer :: r_p(:)
   977)   PetscReal, pointer :: accum_p(:), accum_p2(:)
   978)   type(grid_type), pointer :: grid
   979)   type(option_type), pointer :: option
   980)   type(field_type), pointer :: field
   981)   type(patch_type), pointer :: patch
   982)   type(general_auxvar_type), pointer :: general_auxvars(:,:)
   983)   type(global_auxvar_type), pointer :: global_auxvars(:)  
   984)   class(material_auxvar_type), pointer :: material_auxvars(:)  
   985)   type(material_parameter_type), pointer :: material_parameter
   986)   PetscInt :: local_id, ghosted_id
   987)   PetscInt :: offset , ival, idof
   988) #ifdef DEBUG_GENERAL_INFO
   989)   PetscInt :: icell_max_rel_update(3), icell_max_scaled_residual(3)
   990)   PetscInt :: istate_max_rel_update(3), istate_max_scaled_residual(3)
   991)   character(len=2) :: state_char
   992) #endif
   993)   PetscReal :: dX_X0, R_A, R
   994)   PetscReal :: inf_norm_rel_update(3,3), global_inf_norm_rel_update(3,3)
   995)   PetscReal :: inf_norm_scaled_residual(3,3), global_inf_norm_scaled_residual(3,3)
   996)   PetscReal :: inf_norm_update(3,3), global_inf_norm_update(3,3)
   997)   PetscReal :: inf_norm_residual(3,3), global_inf_norm_residual(3,3)
   998)   PetscReal :: two_norm_residual(3,3), global_two_norm_residual(3,3)
   999)   PetscReal, parameter :: inf_pres_tol = 1.d-1
  1000)   PetscReal, parameter :: inf_temp_tol = 1.d-5
  1001)   PetscReal, parameter :: inf_sat_tol = 1.d-6
  1002)   PetscReal, parameter :: inf_xmol_tol = 1.d-6
  1003)   !geh: note the scaling by 0.d0 several lines down which prevent false 
  1004)   !     convergence
  1005)   PetscReal, parameter :: inf_norm_update_tol(3,3) = &
  1006)     reshape([inf_pres_tol,inf_xmol_tol,inf_temp_tol, &
  1007)              inf_pres_tol,inf_pres_tol,inf_temp_tol, &
  1008)              inf_pres_tol,inf_pres_tol,inf_sat_tol], &
  1009)             shape(inf_norm_update_tol)) * &
  1010)             0.d0
  1011)   PetscReal :: temp(3,12), global_temp(3,12)
  1012)   PetscMPIInt :: mpi_int
  1013)   PetscBool :: converged_abs_update
  1014)   PetscBool :: converged_rel_update
  1015)   PetscBool :: converged_scaled_residual
  1016)   PetscInt :: istate
  1017)   PetscReal :: t_over_v
  1018)   PetscReal :: two_norm
  1019)   
  1020)   grid => this%realization%patch%grid
  1021)   option => this%realization%option
  1022)   field => this%realization%field
  1023)   patch => this%realization%patch
  1024)   general_auxvars => patch%aux%General%auxvars
  1025)   global_auxvars => patch%aux%Global%auxvars
  1026)   material_auxvars => patch%aux%Material%auxvars
  1027)   material_parameter => patch%aux%Material%material_parameter
  1028)   
  1029)   dX_changed = PETSC_FALSE
  1030)   X1_changed = PETSC_FALSE
  1031)   
  1032)   option%converged = PETSC_FALSE
  1033)   if (this%check_post_convergence) then
  1034)     call VecGetArrayReadF90(dX,dX_p,ierr);CHKERRQ(ierr)
  1035)     call VecGetArrayReadF90(X0,X0_p,ierr);CHKERRQ(ierr)
  1036)     call VecGetArrayReadF90(field%flow_r,r_p,ierr);CHKERRQ(ierr)
  1037)     call VecGetArrayReadF90(field%flow_accum,accum_p,ierr);CHKERRQ(ierr)
  1038)     call VecGetArrayReadF90(field%flow_accum2,accum_p2,ierr);CHKERRQ(ierr)
  1039) #ifdef DEBUG_GENERAL_INFO
  1040)     icell_max_rel_update = 0
  1041)     istate_max_rel_update = 0
  1042)     icell_max_scaled_residual = 0
  1043)     istate_max_scaled_residual = 0
  1044) #endif
  1045)     inf_norm_update(:,:) = -1.d20
  1046)     inf_norm_rel_update(:,:) = -1.d20
  1047)     inf_norm_scaled_residual(:,:) = -1.d20
  1048)     inf_norm_residual(:,:) = -1.d20
  1049)     two_norm_residual(:,:) = 0.d0
  1050)     do local_id = 1, grid%nlmax
  1051)       offset = (local_id-1)*option%nflowdof
  1052)       ghosted_id = grid%nL2G(local_id)
  1053)       if (patch%imat(ghosted_id) <= 0) cycle
  1054)       istate = global_auxvars(ghosted_id)%istate
  1055)       do idof = 1, option%nflowdof
  1056)         ival = offset+idof
  1057)         R = r_p(ival)
  1058) #ifdef DEBUG_GENERAL_INFO
  1059)         two_norm_residual(idof,istate) = two_norm_residual(idof,istate) + R*R
  1060) #endif
  1061)         inf_norm_residual(idof,istate) = max(inf_norm_residual(idof,istate), &
  1062)                                              dabs(R))
  1063)         if (general_tough2_conv_criteria) then
  1064)           !geh: scale by t_over_v to match TOUGH2 residual units. see equation
  1065)           !     B.5 of TOUGH2 user manual (LBNL-43134)
  1066)           t_over_v = option%flow_dt/material_auxvars(ghosted_id)%volume
  1067)           if (accum_p2(ival)*t_over_v < general_tough2_itol_scaled_res_e2) then
  1068)             R_A = dabs(R*t_over_v)
  1069)           else
  1070)             R_A = dabs(R/accum_p2(ival))
  1071)           endif
  1072)         else
  1073)           R_A = dabs(R/accum_p(ival))
  1074)         endif
  1075)         dX_X0 = dabs(dX_p(ival)/X0_p(ival))
  1076)         inf_norm_update(idof,istate) = max(inf_norm_update(idof,istate), &
  1077)                                            dabs(dX_p(ival)))
  1078)         if (inf_norm_rel_update(idof,istate) < dX_X0) then
  1079) #ifdef DEBUG_GENERAL_INFO
  1080)           if (maxval(inf_norm_rel_update(idof,:)) < dX_X0) then
  1081)             icell_max_rel_update(idof) = grid%nG2A(ghosted_id)
  1082)             istate_max_rel_update(idof) = global_auxvars(ghosted_id)%istate
  1083)           endif
  1084) #endif
  1085)           inf_norm_rel_update(idof,istate) = dX_X0
  1086)         endif
  1087)         if (inf_norm_scaled_residual(idof,istate) < R_A) then
  1088) #ifdef DEBUG_GENERAL_INFO
  1089)           if (maxval(inf_norm_scaled_residual(idof,:)) < R_A) then
  1090)             icell_max_scaled_residual(idof) = grid%nG2A(ghosted_id)
  1091)             istate_max_scaled_residual(idof) = global_auxvars(ghosted_id)%istate
  1092)           endif
  1093) #endif
  1094)           inf_norm_scaled_residual(idof,istate) = R_A
  1095)         endif
  1096)       enddo
  1097)     enddo
  1098)     temp(1:3,1:3) = inf_norm_update(:,:)
  1099)     temp(1:3,4:6) = inf_norm_rel_update(:,:)
  1100)     temp(1:3,7:9) = inf_norm_scaled_residual(:,:)
  1101)     temp(1:3,10:12) = inf_norm_residual(:,:)
  1102)     mpi_int = 36
  1103)     call MPI_Allreduce(temp,global_temp,mpi_int, &
  1104)                        MPI_DOUBLE_PRECISION,MPI_MAX,option%mycomm,ierr)
  1105)     global_inf_norm_update(:,:) = global_temp(1:3,1:3)
  1106)     global_inf_norm_rel_update(:,:) = global_temp(1:3,4:6)
  1107)     global_inf_norm_scaled_residual(:,:) = global_temp(1:3,7:9)
  1108)     global_inf_norm_residual(:,:) = global_temp(1:3,10:12)
  1109) 
  1110) #ifdef DEBUG_GENERAL_INFO
  1111)     mpi_int = 9
  1112)     call MPI_Allreduce(two_norm_residual,global_two_norm_residual,mpi_int, &
  1113)                        MPI_DOUBLE_PRECISION,MPI_SUM,option%mycomm,ierr)
  1114)     global_two_norm_residual = sqrt(global_two_norm_residual)
  1115) #endif
  1116) 
  1117)     converged_abs_update = PETSC_TRUE
  1118)     converged_scaled_residual = PETSC_TRUE
  1119)     do istate = 1, 3
  1120)       do idof = 1, option%nflowdof
  1121)         if (global_inf_norm_update(idof,istate) > &
  1122)           inf_norm_update_tol(idof,istate)) then
  1123)           converged_abs_update = PETSC_FALSE
  1124)         endif
  1125)         if (general_tough2_conv_criteria) then
  1126)           if (global_inf_norm_scaled_residual(idof,istate) > &
  1127)             general_tough2_itol_scaled_res_e1(idof,istate)) then
  1128)             converged_scaled_residual = PETSC_FALSE
  1129)           endif
  1130)         endif
  1131)       enddo  
  1132)     enddo  
  1133)     converged_rel_update = maxval(global_inf_norm_rel_update) < &
  1134)                                   general_itol_rel_update
  1135)     if (.not.general_tough2_conv_criteria) then
  1136)       converged_scaled_residual = maxval(global_inf_norm_scaled_residual) < &
  1137)                                   general_itol_scaled_res
  1138)     endif
  1139) #if 0
  1140)     do idof = 1, option%nflowdof
  1141)       if (global_inf_norm(idof) > option%flow%post_convergence_tol) then
  1142)         converged_rel_update = PETSC_FALSE
  1143) #ifdef DEBUG_GENERAL_INFO
  1144)         select case(istate_max(idof))
  1145)           case(1)
  1146)             state_char = 'L'
  1147)           case(2)
  1148)             state_char = 'G'
  1149)           case(3)
  1150)             state_char = '2P'
  1151)         end select
  1152)         write(*,'(''-+ '',a3,i2,''('',i5,''):'',es12.4, &
  1153)                  &'' dX_X/dX/X:'',3es12.4, &
  1154)                  &'' R_A/R/A:'',3es12.4)') state_char,idof, &
  1155)            icell_max(idof),global_inf_norm(idof), &
  1156)            dX_X1_max(idof), dX_max(idof),  X1_max(idof), &
  1157)            R_A_max(idof), R_max(idof), A_max(idof)
  1158) #endif
  1159)       endif
  1160)     enddo
  1161) #endif
  1162) #ifdef DEBUG_GENERAL_INFO
  1163)     write(*,'(4x,''-+  dpl:'',es12.4,''  dxa:'',es12.4,''  dt:'',es12.4)') &
  1164)       (max(global_inf_norm_update(idof,1),0.d0),idof=1,3)
  1165)     write(*,'(4x,''-+  dpg:'',es12.4,''  dpa:'',es12.4,''  dt:'',es12.4)') &
  1166)       (max(global_inf_norm_update(idof,2),0.d0),idof=1,3)
  1167)     if (general_2ph_energy_dof == GENERAL_TEMPERATURE_INDEX) then
  1168)       write(*,'(4x,''-+  dpg:'',es12.4,''  dsg:'',es12.4,''  dt:'',es12.4)') &
  1169)         (max(global_inf_norm_update(idof,3),0.d0),idof=1,3)
  1170)     else
  1171)       write(*,'(4x,''-+  dpg:'',es12.4,''  dsg:'',es12.4,'' dpa:'',es12.4)') &
  1172)         (max(global_inf_norm_update(idof,3),0.d0),idof=1,3)
  1173)     endif
  1174)     write(*,'(4x,''-+ rupl:'',es12.4,'' ruxa:'',es12.4,'' rut:'',es12.4)') &
  1175)       (max(global_inf_norm_rel_update(idof,1),0.d0),idof=1,3)
  1176)     write(*,'(4x,''-+ rupg:'',es12.4,'' rupa:'',es12.4,'' rut:'',es12.4)') &
  1177)       (max(global_inf_norm_rel_update(idof,2),0.d0),idof=1,3)
  1178)     write(*,'(4x,''-+ rupg:'',es12.4,'' rusg:'',es12.4,'' rut:'',es12.4)') &
  1179)       (max(global_inf_norm_rel_update(idof,3),0.d0),idof=1,3)
  1180)     write(*,'(4x,''-+  srl:'',es12.4,''  srg:'',es12.4,'' sre:'',es12.4)') &
  1181)       (max(global_inf_norm_scaled_residual(idof,1),0.d0),idof=1,3)
  1182)     write(*,'(4x,''-+  srl:'',es12.4,''  srg:'',es12.4,'' sre:'',es12.4)') &
  1183)       (max(global_inf_norm_scaled_residual(idof,2),0.d0),idof=1,3)
  1184)     write(*,'(4x,''-+  srl:'',es12.4,''  srg:'',es12.4,'' sre:'',es12.4)') &
  1185)       (max(global_inf_norm_scaled_residual(idof,3),0.d0),idof=1,3)
  1186)     write(*,'(4x,''-+ ru1 icell:'',i7,''  st:'',i3,''  X:'',es11.3, &
  1187)               &''  dX:'',es11.3,''  R:'',es11.3)') &
  1188)       icell_max_rel_update(1), istate_max_rel_update(1), &
  1189)       X0_p((icell_max_rel_update(1)-1)*3+1), &
  1190)       -1.d0*dX_p((icell_max_rel_update(1)-1)*3+1), &
  1191)       r_p((icell_max_rel_update(1)-1)*3+1)
  1192)     write(*,'(4x,''-+ ru2 icell:'',i7,''  st:'',i3,''  X:'',es11.3, &
  1193)               &''  dX:'',es11.3,''  R:'',es11.3)') &
  1194)       icell_max_rel_update(2), istate_max_rel_update(2), &
  1195)       X0_p((icell_max_rel_update(2)-1)*3+2), &
  1196)       -1.d0*dX_p((icell_max_rel_update(2)-1)*3+2), &
  1197)       r_p((icell_max_rel_update(2)-1)*3+2)
  1198)     write(*,'(4x,''-+ ru3 icell:'',i7,''  st:'',i3,''  X:'',es11.3, &
  1199)               &''  dX:'',es11.3,''  R:'',es11.3)') &
  1200)       icell_max_rel_update(3), istate_max_rel_update(3), &
  1201)       X0_p((icell_max_rel_update(3)-1)*3+3), &
  1202)       -1.d0*dX_p((icell_max_rel_update(3)-1)*3+3), &
  1203)       r_p((icell_max_rel_update(3)-1)*3+3)
  1204)     write(*,'(4x,''-+ sr1 icell:'',i7,''  st:'',i3,''  X:'',es11.3, &
  1205)               &''  dX:'',es11.3,''  R:'',es11.3)') &
  1206)       icell_max_scaled_residual(1), istate_max_scaled_residual(1), &
  1207)       X0_p((icell_max_scaled_residual(1)-1)*3+1), &
  1208)       -1.d0*dX_p((icell_max_scaled_residual(1)-1)*3+1), &
  1209)       r_p((icell_max_scaled_residual(1)-1)*3+1)
  1210)     write(*,'(4x,''-+ sr2 icell:'',i7,''  st:'',i3,''  X:'',es11.3, &
  1211)               &''  dX:'',es11.3,''  R:'',es11.3)') &
  1212)       icell_max_scaled_residual(2), istate_max_scaled_residual(2), &
  1213)       X0_p((icell_max_scaled_residual(2)-1)*3+2), &
  1214)       -1.d0*dX_p((icell_max_scaled_residual(2)-1)*3+2), &
  1215)       r_p((icell_max_scaled_residual(2)-1)*3+2)
  1216)     write(*,'(4x,''-+ sr3 icell:'',i7,''  st:'',i3,''  X:'',es11.3, &
  1217)               &''  dX:'',es11.3,''  R:'',es11.3)') &
  1218)       icell_max_scaled_residual(3), istate_max_scaled_residual(3), &
  1219)       X0_p((icell_max_scaled_residual(3)-1)*3+3), &
  1220)       -1.d0*dX_p((icell_max_scaled_residual(3)-1)*3+3), &
  1221)       r_p((icell_max_scaled_residual(3)-1)*3+3)
  1222) #endif
  1223)     option%converged = PETSC_FALSE
  1224)     if (converged_abs_update .or. converged_rel_update .or. &
  1225)         converged_scaled_residual) then
  1226)       option%converged = PETSC_TRUE
  1227)     endif
  1228)     call VecRestoreArrayReadF90(dX,dX_p,ierr);CHKERRQ(ierr)
  1229)     call VecRestoreArrayReadF90(X0,X0_p,ierr);CHKERRQ(ierr)
  1230)     call VecRestoreArrayReadF90(field%flow_r,r_p,ierr);CHKERRQ(ierr)
  1231)     call VecRestoreArrayReadF90(field%flow_accum,accum_p,ierr);CHKERRQ(ierr)
  1232)     call VecRestoreArrayReadF90(field%flow_accum2,accum_p2,ierr);CHKERRQ(ierr)
  1233)                                
  1234)     if (this%print_ekg) then
  1235)       call VecNorm(field%flow_r,NORM_2,two_norm,ierr);CHKERRQ(ierr)
  1236)       if (OptionPrintToFile(option)) then
  1237)   100 format("GENERAL NEWTON_ITERATION ",100es11.3)
  1238)         write(IUNIT_EKG,100) &
  1239)           global_inf_norm_update(:,:), &
  1240)           global_inf_norm_rel_update(:,:), &
  1241)           global_inf_norm_scaled_residual(:,:), &
  1242)           global_inf_norm_residual(:,:), &
  1243)           two_norm
  1244)       endif    
  1245)     endif                               
  1246)   endif                               
  1247) 
  1248) end subroutine PMGeneralCheckUpdatePost
  1249) 
  1250) ! ************************************************************************** !
  1251) 
  1252) subroutine PMGeneralTimeCut(this)
  1253)   ! 
  1254)   ! Author: Glenn Hammond
  1255)   ! Date: 03/14/13
  1256)   ! 
  1257) 
  1258)   use General_module, only : GeneralTimeCut
  1259) 
  1260)   implicit none
  1261)   
  1262)   class(pm_general_type) :: this
  1263)   
  1264)   call PMSubsurfaceFlowTimeCut(this)
  1265)   call GeneralTimeCut(this%realization)
  1266) 
  1267) end subroutine PMGeneralTimeCut
  1268) 
  1269) ! ************************************************************************** !
  1270) 
  1271) subroutine PMGeneralUpdateSolution(this)
  1272)   ! 
  1273)   ! Author: Glenn Hammond
  1274)   ! Date: 03/14/13
  1275)   ! 
  1276) 
  1277)   use General_module, only : GeneralUpdateSolution, &
  1278)                              GeneralMapBCAuxVarsToGlobal
  1279) 
  1280)   implicit none
  1281)   
  1282)   class(pm_general_type) :: this
  1283)   
  1284)   call PMSubsurfaceFlowUpdateSolution(this)
  1285)   call GeneralUpdateSolution(this%realization)
  1286)   call GeneralMapBCAuxVarsToGlobal(this%realization)
  1287) 
  1288) end subroutine PMGeneralUpdateSolution     
  1289) 
  1290) ! ************************************************************************** !
  1291) 
  1292) subroutine PMGeneralUpdateAuxVars(this)
  1293)   ! 
  1294)   ! Author: Glenn Hammond
  1295)   ! Date: 04/21/14
  1296)   use General_module, only : GeneralUpdateAuxVars
  1297) 
  1298)   implicit none
  1299)   
  1300)   class(pm_general_type) :: this
  1301) 
  1302)   call GeneralUpdateAuxVars(this%realization,PETSC_FALSE)
  1303) 
  1304) end subroutine PMGeneralUpdateAuxVars   
  1305) 
  1306) ! ************************************************************************** !
  1307) 
  1308) subroutine PMGeneralMaxChange(this)
  1309)   ! 
  1310)   ! Not needed given GeneralMaxChange is called in PostSolve
  1311)   ! 
  1312)   ! Author: Glenn Hammond
  1313)   ! Date: 03/14/13
  1314)   ! 
  1315) 
  1316)   use Realization_Base_class
  1317)   use Realization_Subsurface_class
  1318)   use Option_module
  1319)   use Field_module
  1320)   use Grid_module
  1321)   use Global_Aux_module
  1322)   use General_Aux_module
  1323)   use Variables_module, only : LIQUID_PRESSURE, LIQUID_MOLE_FRACTION, &
  1324)                                TEMPERATURE, GAS_PRESSURE, AIR_PRESSURE, &
  1325)                                GAS_SATURATION
  1326)   implicit none
  1327)   
  1328)   class(pm_general_type) :: this
  1329)   
  1330)   class(realization_subsurface_type), pointer :: realization
  1331)   type(option_type), pointer :: option
  1332)   type(field_type), pointer :: field
  1333)   type(grid_type), pointer :: grid
  1334)   PetscReal, pointer :: vec_ptr(:), vec_ptr2(:)
  1335)   PetscReal :: max_change_local(6)
  1336)   PetscReal :: max_change_global(6)
  1337)   PetscReal :: max_change
  1338)   PetscInt :: i, j
  1339)   PetscInt :: local_id, ghosted_id
  1340) 
  1341)   
  1342)   PetscErrorCode :: ierr
  1343)   
  1344)   realization => this%realization
  1345)   option => realization%option
  1346)   field => realization%field
  1347)   grid => realization%patch%grid
  1348) 
  1349)   max_change_global = 0.d0
  1350)   max_change_local = 0.d0
  1351)   
  1352)   ! max change variables: [LIQUID_PRESSURE, GAS_PRESSURE, AIR_PRESSURE, &
  1353)   !                        LIQUID_MOLE_FRACTION, TEMPERATURE, GAS_SATURATION]
  1354)   do i = 1, 6
  1355)     call RealizationGetVariable(realization,field%work, &
  1356)                                 this%max_change_ivar(i), &
  1357)                                 this%max_change_isubvar(i))
  1358)     ! yes, we could use VecWAXPY and a norm here, but we need the ability
  1359)     ! to customize
  1360)     call VecGetArrayF90(field%work,vec_ptr,ierr);CHKERRQ(ierr)
  1361)     call VecGetArrayF90(field%max_change_vecs(i),vec_ptr2,ierr);CHKERRQ(ierr)
  1362)     max_change = 0.d0
  1363)     do j = 1, grid%nlmax
  1364)       ! have to weed out cells that changed state
  1365)       if (dabs(vec_ptr(j)) > 1.d-40 .and. dabs(vec_ptr2(j)) > 1.d-40) then
  1366)         max_change = max(max_change,dabs(vec_ptr(j)-vec_ptr2(j)))
  1367)       endif
  1368)     enddo
  1369)     max_change_local(i) = max_change
  1370)     call VecRestoreArrayF90(field%work,vec_ptr,ierr);CHKERRQ(ierr)
  1371)     call VecRestoreArrayF90(field%max_change_vecs(i),vec_ptr2, &
  1372)                             ierr);CHKERRQ(ierr)
  1373)     call VecCopy(field%work,field%max_change_vecs(i),ierr);CHKERRQ(ierr)
  1374)   enddo
  1375)   call MPI_Allreduce(max_change_local,max_change_global,SIX_INTEGER, &
  1376)                       MPI_DOUBLE_PRECISION,MPI_MAX,option%mycomm,ierr)
  1377)   ! print them out
  1378)   if (OptionPrintToScreen(option)) then
  1379)     write(*,'("  --> max chng: dpl= ",1pe12.4, " dpg= ",1pe12.4,&
  1380)       & " dpa= ",1pe12.4,/,15x," dxa= ",1pe12.4,"  dt= ",1pe12.4,&
  1381)       & " dsg= ",1pe12.4)') &
  1382)       max_change_global(1:6)
  1383)   endif
  1384)   if (OptionPrintToFile(option)) then
  1385)     write(option%fid_out,'("  --> max chng: dpl= ",1pe12.4, " dpg= ",1pe12.4,&
  1386)       & " dpa= ",1pe12.4,/,15x," dxa= ",1pe12.4,"  dt= ",1pe12.4, &
  1387)       & " dsg= ",1pe12.4)') &
  1388)       max_change_global(1:6)
  1389)   endif
  1390) 
  1391)   ! max change variables: [LIQUID_PRESSURE, GAS_PRESSURE, AIR_PRESSURE, &
  1392)   !                        LIQUID_MOLE_FRACTION, TEMPERATURE, GAS_SATURATION]
  1393)   ! ignore air pressure as it jumps during phase change
  1394)   this%max_pressure_change = maxval(max_change_global(1:2))
  1395)   this%max_xmol_change = max_change_global(4)
  1396)   this%max_temperature_change = max_change_global(5)
  1397)   this%max_saturation_change = max_change_global(6)
  1398)   
  1399) end subroutine PMGeneralMaxChange
  1400) 
  1401) ! ************************************************************************** !
  1402) 
  1403) subroutine PMGeneralComputeMassBalance(this,mass_balance_array)
  1404)   ! 
  1405)   ! Author: Glenn Hammond
  1406)   ! Date: 03/14/13
  1407)   ! 
  1408) 
  1409)   use General_module, only : GeneralComputeMassBalance
  1410) 
  1411)   implicit none
  1412)   
  1413)   class(pm_general_type) :: this
  1414)   PetscReal :: mass_balance_array(:)
  1415)   
  1416)   call GeneralComputeMassBalance(this%realization,mass_balance_array)
  1417) 
  1418) end subroutine PMGeneralComputeMassBalance
  1419) 
  1420) ! ************************************************************************** !
  1421) 
  1422) subroutine PMGeneralInputRecord(this)
  1423)   ! 
  1424)   ! Writes ingested information to the input record file.
  1425)   ! 
  1426)   ! Author: Jenn Frederick, SNL
  1427)   ! Date: 03/21/2016
  1428)   ! 
  1429)   
  1430)   implicit none
  1431)   
  1432)   class(pm_general_type) :: this
  1433) 
  1434)   character(len=MAXWORDLENGTH) :: word
  1435)   PetscInt :: id
  1436) 
  1437)   id = INPUT_RECORD_UNIT
  1438) 
  1439)   write(id,'(a29)',advance='no') 'pm: '
  1440)   write(id,'(a)') this%name
  1441)   write(id,'(a29)',advance='no') 'mode: '
  1442)   write(id,'(a)') 'general'
  1443)   if (this%check_post_convergence) then
  1444)     write(id,'(a29)',advance='no') 'ITOL_SCALED_RESIDUAL: '
  1445)     write(id,'(a)') 'ON'
  1446)     write(id,'(a29)',advance='no') 'ITOL_RELATIVE_RESIDUAL: '
  1447)     write(id,'(a)') 'ON'
  1448)   endif
  1449) 
  1450) end subroutine PMGeneralInputRecord
  1451) 
  1452) ! ************************************************************************** !
  1453) 
  1454) subroutine PMGeneralCheckpointBinary(this,viewer)
  1455)   ! 
  1456)   ! Checkpoints data associated with General PM
  1457)   ! 
  1458)   ! Author: Glenn Hammond
  1459)   ! Date: 02/18/15
  1460) 
  1461)   use Checkpoint_module
  1462)   use Global_module
  1463)   use Variables_module, only : STATE
  1464) 
  1465)   implicit none
  1466) #include "petsc/finclude/petscviewer.h"      
  1467) 
  1468)   class(pm_general_type) :: this
  1469)   PetscViewer :: viewer
  1470)   
  1471)   call GlobalGetAuxVarVecLoc(this%realization, &
  1472)                              this%realization%field%iphas_loc, &
  1473)                              STATE,ZERO_INTEGER)
  1474)   call PMSubsurfaceFlowCheckpointBinary(this,viewer)
  1475)   
  1476) end subroutine PMGeneralCheckpointBinary
  1477) 
  1478) ! ************************************************************************** !
  1479) 
  1480) subroutine PMGeneralRestartBinary(this,viewer)
  1481)   ! 
  1482)   ! Restarts data associated with General PM
  1483)   ! 
  1484)   ! Author: Glenn Hammond
  1485)   ! Date: 02/18/15
  1486) 
  1487)   use Checkpoint_module
  1488)   use Global_module
  1489)   use Variables_module, only : STATE
  1490) 
  1491)   implicit none
  1492) #include "petsc/finclude/petscviewer.h"      
  1493) 
  1494)   class(pm_general_type) :: this
  1495)   PetscViewer :: viewer
  1496)   
  1497)   call PMSubsurfaceFlowRestartBinary(this,viewer)
  1498)   call GlobalSetAuxVarVecLoc(this%realization, &
  1499)                              this%realization%field%iphas_loc, &
  1500)                              STATE,ZERO_INTEGER)
  1501)   
  1502) end subroutine PMGeneralRestartBinary
  1503) ! ************************************************************************** !
  1504) 
  1505) subroutine PMGeneralDestroy(this)
  1506)   ! 
  1507)   ! Destroys General process model
  1508)   ! 
  1509)   ! Author: Glenn Hammond
  1510)   ! Date: 03/14/13
  1511)   ! 
  1512) 
  1513)   use General_module, only : GeneralDestroy
  1514) 
  1515)   implicit none
  1516)   
  1517)   class(pm_general_type) :: this
  1518)   
  1519)   if (associated(this%next)) then
  1520)     call this%next%Destroy()
  1521)   endif
  1522) 
  1523)   deallocate(this%max_change_ivar)
  1524)   nullify(this%max_change_ivar)
  1525)   deallocate(this%max_change_isubvar)
  1526)   nullify(this%max_change_isubvar)
  1527) 
  1528)   ! preserve this ordering
  1529)   call GeneralDestroy(this%realization)
  1530)   call PMSubsurfaceFlowDestroy(this)
  1531)   
  1532) end subroutine PMGeneralDestroy
  1533)   
  1534) end module PM_General_class

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