ICFERST  22-06
Reservoir simulator based on DCVFEM, Dynamic Mesh optimisation and Surface-based modelling
solvers_module Module Reference

All the subroutines associated to solving non-linear systems, Schur complement, etc. More...

Functions/Subroutines

subroutine, public boundedsolutioncorrections (state, packed_state, Mdims, CV_funs, small_findrm, small_colm, Field_name, for_sat, min_max_limits)
 The sparcity of the local CV connectivity is in: small_findrm, small_colm. ngl_its=max no of global iterations e.g. 100. error_tol = tolerance on the iterations. More...
 
subroutine, public fpi_backtracking (nphase, Mdims, ndgln, state, packed_state, sat_bak, backtrack_sat, backtrack_par_from_schema, Previous_convergence, satisfactory_convergence, new_backtrack_par, Max_sat_its, its, nonlinear_iteration, useful_sats, res, res_ratio, first_res)
 :In this subroutine we applied some corrections and backtrack_par on the saturations obtained from the saturation equation this subroutine is detailed in Salinas et al. 2017 doi:10.1002/fld.4357 The method ensures convergence "independent" of the time step. More...
 
subroutine, public set_saturation_to_sum_one (mdims, packed_state, state, do_not_update_halos)
 :This subroutines eliminates the oscillations in the saturation that are bigger than a certain tolerance and also sets the saturation to be between bounds More...
 
subroutine, public non_porous_ensure_sum_to_one (Mdims, packed_state, do_not_update_halos)
 : This subroutines eliminates the oscillations in the saturation that are bigger than a certain tolerance and also sets the saturation to be between bounds More...
 
subroutine, public initialise_saturation_sums_one (mdims, packed_state, find_scapegoat_phase)
 :Ensure that the saturations at the beginning sum to one, if they do not all the error is compensated in the scapegoat_phase. Normally the last More...
 
subroutine, public auto_backtracking (Mdims, backtrack_par_factor, courant_number_in, first_time_step, nonlinear_iteration, I_am_temperature)
 : The maximum backtracking factor is calculated based on the Courant number and physical effects ocurring in the domain More...
 
subroutine, public get_anderson_acceleration_new_guess (N, M, NewField, History_field, stored_residuals, max_its, prev_small_matrix, restart_now)
 This subroutine provides a new guess based on previous updates and residuals. Use this subroutine to speed up any system being solved by looping. A new guess is computed and returned Method explained in DOI.10.1137/10078356X. More...
 
subroutine, public scale_petsc_matrix (Mat_petsc)
 In this subroutine the matrix is re-scaled based on the formula D^-0.5 * A * D^-0.5 X'= D^-0.5 b; and next X = D^-0.5 * X'; IMPORTANT: the step X = D^-0.5 * X' needs to be done elsewhere store the diagonal before calling this This should allow to deal with high ranges of viscosity ratio for example A is-written. More...
 
subroutine, public duplicate_petsc_matrix (MAT_A, MAT_B)
 
subroutine, public petsc_stokes_solver (packed_state, Mdims, Mmat, ndgln, Mspars, final_phase, pmat, P_all, deltaP, rhs_p, solver_option_path, Dmat)
 In this subroutine the Schur complement is generated and solved using PETSc to update the pressure field Matrices need to be in petsc format and pmat is the preconditioned matrix, i.e. pmat = A11- A10(AproxA00^-1)A01. More...
 
subroutine convert_c_and_ct_mat_to_petsc_format (packed_state, Mdims, Mmat, ndgln, Mspars, NPHASE)
 : This subroutine converts the C (gradient) and CT (Divergence) matrices into PETSc format Mainly devoted to be used by the PETSc stokes schur solver This can be used when moving from the ICFERST matrix format to PETSc as reference More...
 

Detailed Description

All the subroutines associated to solving non-linear systems, Schur complement, etc.

Function/Subroutine Documentation

◆ auto_backtracking()

subroutine, public solvers_module::auto_backtracking ( type(multi_dimensions), intent(in)  Mdims,
real, intent(inout)  backtrack_par_factor,
real, dimension(:), intent(inout)  courant_number_in,
logical, intent(in)  first_time_step,
integer, intent(in)  nonlinear_iteration,
logical, intent(in), optional  I_am_temperature 
)

: The maximum backtracking factor is calculated based on the Courant number and physical effects ocurring in the domain

Parameters
MdimsData type storing all the dimensions describing the mesh, fields, nodes, etc
backtrack_par_factorINOUT backtracking value to use
courant_number_inCourant number and shock front courant number
first_time_steptrue if the first one
nonlinear_iterationcurrent non-linear iteration
I_am_temperatureIf doing temperature true (DEPRECATED PROBABLY)
Here is the caller graph for this function:

◆ boundedsolutioncorrections()

subroutine, public solvers_module::boundedsolutioncorrections ( type( state_type ), dimension( : ), intent(inout)  state,
type( state_type ), intent(inout)  packed_state,
type(multi_dimensions), intent(in)  Mdims,
type(multi_shape_funs), intent(in)  CV_funs,
integer, dimension( : ), intent(in)  small_findrm,
integer, dimension( : ), intent(in)  small_colm,
character( len=* ), intent(in)  Field_name,
logical, intent(in), optional  for_sat,
real, dimension(2), optional  min_max_limits 
)

The sparcity of the local CV connectivity is in: small_findrm, small_colm. ngl_its=max no of global iterations e.g. 100. error_tol = tolerance on the iterations.

nloc_its: This iteration is very good at avoiding spreading the modifications too far - however it can stagnate. nloc_its2: This iteration is very good at avoiding stagnating but does spread the modifcations far. us a single iteration because of this as default... nits_nod: iterations at a nod - this iteration is very good at avoiding spreading the modifications too far - however it can stagnate.

Parameters
stateLinked list containing all the fields defined in diamond and considered by Fluidity
packed_stateLinked list containing all the fields used by IC-FERST, memory partially shared with state
MdimsData type storing all the dimensions describing the mesh, fields, nodes, etc
CV_funsShape functions for the CV mesh
small_findrmsparcity of the local CV connectivity is in: small_findrm, small_colm.
small_colmsparcity of the local CV connectivity is in: small_findrm, small_colm.
Field_nameName of the field
for_satTrue if doing it for Phase Volume Fraction
min_max_limitsMaximum and minimum values
Here is the call graph for this function:
Here is the caller graph for this function:

◆ convert_c_and_ct_mat_to_petsc_format()

subroutine solvers_module::convert_c_and_ct_mat_to_petsc_format ( type( state_type ), intent(inout)  packed_state,
type(multi_dimensions), intent(in)  Mdims,
type(multi_matrices), intent(inout)  Mmat,
type(multi_ndgln), intent(in)  ndgln,
type(multi_sparsities), intent(in)  Mspars,
integer, intent(in)  NPHASE 
)

: This subroutine converts the C (gradient) and CT (Divergence) matrices into PETSc format Mainly devoted to be used by the PETSc stokes schur solver This can be used when moving from the ICFERST matrix format to PETSc as reference

Parameters
packed_stateLinked list containing all the fields used by IC-FERST, memory partially shared with state
MdimsData type storing all the dimensions describing the mesh, fields, nodes, etc
MmatMatrices for ICFERST
ndglnGlobal to local variables
MsparsSparsity of the matrices
nphasenumber of phases
Here is the caller graph for this function:

◆ duplicate_petsc_matrix()

subroutine, public solvers_module::duplicate_petsc_matrix ( type(petsc_csr_matrix), intent(in)  MAT_A,
type(petsc_csr_matrix), intent(inout)  MAT_B 
)

◆ fpi_backtracking()

subroutine, public solvers_module::fpi_backtracking ( integer, intent(in)  nphase,
type( multi_dimensions ), intent(in)  Mdims,
type(multi_ndgln), intent(in)  ndgln,
type( state_type ), dimension( : ), intent(in)  state,
type( state_type ), intent(inout)  packed_state,
real, dimension(:, :), intent(in)  sat_bak,
real, dimension(:, :), intent(in)  backtrack_sat,
real, intent(in)  backtrack_par_from_schema,
real, intent(inout)  Previous_convergence,
logical, intent(inout)  satisfactory_convergence,
real, intent(inout)  new_backtrack_par,
integer, intent(in)  Max_sat_its,
integer, intent(in)  its,
integer, intent(in)  nonlinear_iteration,
integer, intent(inout)  useful_sats,
real, intent(in)  res,
real, intent(in)  res_ratio,
real, intent(in)  first_res 
)

:In this subroutine we applied some corrections and backtrack_par on the saturations obtained from the saturation equation this subroutine is detailed in Salinas et al. 2017 doi:10.1002/fld.4357 The method ensures convergence "independent" of the time step.

Parameters
nphaseNumber of phases
MdimsData type storing all the dimensions describing the mesh, fields, nodes, etc
ndglnGlobal to local variables
stateLinked list containing all the fields defined in diamond and considered by Fluidity
packed_stateLinked list containing all the fields used by IC-FERST, memory partially shared with state
sat_bakBackup of the saturation field
backtrack_satBacktrack parameter used
backtrack_par_from_schemaobtained from diamond/input file
Previous_convergenceConvergence obtained in the previous attempt
satisfactory_convergenceIf converged true
new_backtrack_parNew obtained backtracking parameter
Max_sat_itsMaximum allowed number of saturation iterations
itsIterations taken?
nonlinear_iterationcurrent non_linear iteration
useful_satsNumber of useful iterations of the saturation iteration method
resresidual obtained using the new sigmas obtained with the new saturation in the transport equation
res_ratioratio of residuals between to iterations
first_resresidual obtained at the first attemp, used as reference
Here is the call graph for this function:
Here is the caller graph for this function:

◆ get_anderson_acceleration_new_guess()

subroutine, public solvers_module::get_anderson_acceleration_new_guess ( integer, intent(in)  N,
integer, intent(in)  M,
real, dimension(n), intent(out)  NewField,
real, dimension(n, m+2), intent(inout)  History_field,
real, dimension(n, m+1), intent(inout)  stored_residuals,
integer, optional  max_its,
real, dimension(:,:), intent(inout), optional, allocatable  prev_small_matrix,
logical, intent(inout)  restart_now 
)

This subroutine provides a new guess based on previous updates and residuals. Use this subroutine to speed up any system being solved by looping. A new guess is computed and returned Method explained in DOI.10.1137/10078356X.

Author
Pablo Salinas
Parameters
NSize if the field of interest. Used also to turn vector/tensor fields into scalar fields internally here
MSize of the field of iterations available - 2
NewFieldNew guess obtained by the Anderson Acceleration method
History_fieldPast results obtained by the outer solver
stored_residualsPast residuals obtained as the (- guessed value + G(guess value) )
max_itsMaximum number of iterations
prev_small_matrixContains the previous A'*A, speeds up the method. On exit is updated (M,M) Prepared to be passed unallocated from M = 1 and allocated internally
Here is the caller graph for this function:

◆ initialise_saturation_sums_one()

subroutine, public solvers_module::initialise_saturation_sums_one ( type( multi_dimensions ), intent(in)  mdims,
type( state_type ), intent(inout)  packed_state,
logical, intent(in), optional  find_scapegoat_phase 
)

:Ensure that the saturations at the beginning sum to one, if they do not all the error is compensated in the scapegoat_phase. Normally the last

Parameters
MdimsData type storing all the dimensions describing the mesh, fields, nodes, etc
find_scapegoat_phaseAgainst which phase we correct the error

◆ non_porous_ensure_sum_to_one()

subroutine, public solvers_module::non_porous_ensure_sum_to_one ( type( multi_dimensions ), intent(in)  Mdims,
type( state_type ), intent(inout)  packed_state,
logical, intent(in), optional  do_not_update_halos 
)

: This subroutines eliminates the oscillations in the saturation that are bigger than a certain tolerance and also sets the saturation to be between bounds

Parameters
MdimsData type storing all the dimensions describing the mesh, fields, nodes, etc
packed_stateLinked list containing all the fields used by IC-FERST, memory partially shared with state
do_not_update_halosIf true do not update halos, to save time in parallel
Here is the caller graph for this function:

◆ petsc_stokes_solver()

subroutine, public solvers_module::petsc_stokes_solver ( type( state_type ), intent(inout)  packed_state,
type(multi_dimensions), intent(in)  Mdims,
type(multi_matrices), intent(inout)  Mmat,
type(multi_ndgln), intent(in)  ndgln,
type(multi_sparsities), intent(in)  Mspars,
integer, intent(in)  final_phase,
type(petsc_csr_matrix), intent(inout)  pmat,
type( tensor_field ), intent(inout)  P_all,
type( vector_field ), intent(inout)  deltaP,
type( vector_field ), intent(in)  rhs_p,
character(len=option_path_len), intent(in)  solver_option_path,
type(petsc_csr_matrix), intent(inout), optional  Dmat 
)

In this subroutine the Schur complement is generated and solved using PETSc to update the pressure field Matrices need to be in petsc format and pmat is the preconditioned matrix, i.e. pmat = A11- A10(AproxA00^-1)A01.

Author
Pablo Salinas
Parameters
packed_stateLinked list containing all the fields used by IC-FERST, memory partially shared with state
MdimsData type storing all the dimensions describing the mesh, fields, nodes, etc
MsparsSparsity of the matrices
ndglnGlobal to local variables
final_phaseThis is the final phase to be assembled, in this way we can assemble from phase 1 to final_phase not necessarily being for all the phases
pmatINOUT preconditioned matrix, i.e. pmat = A11- A10(AproxA00^-1)A01
P_allPressure field
deltaPUpdate of pressure from before solving
rhs_pRHS of the Pressure system
solver_option_pathPath of the solver to be used to solve here the system
Dmat(optional) is the Matrix multipliying pressure in the continuity equation
Here is the call graph for this function:
Here is the caller graph for this function:

◆ scale_petsc_matrix()

subroutine, public solvers_module::scale_petsc_matrix ( type(petsc_csr_matrix), intent(inout)  Mat_petsc)

In this subroutine the matrix is re-scaled based on the formula D^-0.5 * A * D^-0.5 X'= D^-0.5 b; and next X = D^-0.5 * X'; IMPORTANT: the step X = D^-0.5 * X' needs to be done elsewhere store the diagonal before calling this This should allow to deal with high ranges of viscosity ratio for example A is-written.

Author
Pablo Salinas
Here is the caller graph for this function:

◆ set_saturation_to_sum_one()

subroutine, public solvers_module::set_saturation_to_sum_one ( type( multi_dimensions ), intent(in)  mdims,
type( state_type ), intent(inout)  packed_state,
type( state_type ), dimension(:), intent(in)  state,
logical, intent(in), optional  do_not_update_halos 
)

:This subroutines eliminates the oscillations in the saturation that are bigger than a certain tolerance and also sets the saturation to be between bounds

Parameters
MdimsData type storing all the dimensions describing the mesh, fields, nodes, etc
packed_stateLinked list containing all the fields used by IC-FERST, memory partially shared with state
stateLinked list containing all the fields defined in diamond and considered by Fluidity
do_not_update_halosIf true do not update halos, to save time in parallel
Here is the caller graph for this function: