implicit_PicardBDE_integrator_procedures Submodule

Contains module procedures associated with the implicit BDE integrator class



Functions

function checkConvergence(oldVars, newVars, indicesToCheck, nonlinTol, absTol, use2Norm, convergenceCounter) result(conv)

Checks whether all variables determined by indicesToCheck have converged based on a given nonlinear tolerance

Arguments

Type IntentOptional Attributes Name
type(RealArray), intent(in), dimension(:) :: oldVars

Previous variable values

type(RealArray), intent(in), dimension(:) :: newVars

New variable values

integer(kind=ik), intent(in), dimension(:) :: indicesToCheck

Indices of oldVars/newVars to check for convergence

real(kind=rk), intent(in) :: nonlinTol

Relative convergence tolerances for each variable

real(kind=rk), intent(in) :: absTol

Absolute tolerance in epsilon units for each variable

logical, intent(in) :: use2Norm

True if this should use 2-norm instead of local inf-norm

integer(kind=ik), intent(inout), dimension(:) :: convergenceCounter

Incremented if the variable with a given index has not yet converged

Return Value logical


Subroutines

subroutine tryIntegrate(this, manipulatedModeller, outputVars, inputVars, numSteps, dt, solveSuccess)

Arguments

Type IntentOptional Attributes Name
type(PicardBDEIntegrator), intent(inout) :: this
class(ModellerSurrogate), intent(inout) :: manipulatedModeller

Modeller to be used in callbacks during integration

class(VariableContainer), intent(inout) :: outputVars

VariableContainer object to store the integration output

class(VariableContainer), intent(in), value :: inputVars

VariableContainer object housing input data for the integration routine

integer(kind=ik), intent(in) :: numSteps
real(kind=rk), intent(in), dimension(:) :: dt
logical, intent(inout) :: solveSuccess

Module Functions

pure module function getNonlinTol(this) result(nonlinTol)

Getter for nonlinTol value

Arguments

Type IntentOptional Attributes Name
class(PicardBDEIntegrator), intent(in) :: this

Return Value real(kind=rk)

pure module function getMaxIterations(this) result(maxIter)

Geter for maxIterations value

Arguments

Type IntentOptional Attributes Name
class(PicardBDEIntegrator), intent(in) :: this

Return Value integer(kind=ik)

pure module function getConvergenceIndices(this) result(convVars)

Getter for convergenceTestVars array

Arguments

Type IntentOptional Attributes Name
class(PicardBDEIntegrator), intent(in) :: this

Return Value integer(kind=ik), dimension(:), allocatable


Module Subroutines

module subroutine initBDEIntegrator(this, indexingObj, procRank, nonlinTol, absTol, maxIters, convergenceIndices, modelList, termGroups, evolvesTimeVar, dtController, initialTimestep, use2Norm, petscGroup, intContOptions, integratorName)

BDE integrator constructor

Arguments

Type IntentOptional Attributes Name
class(PicardBDEIntegrator), intent(inout) :: this
type(Indexing), intent(in) :: indexingObj

Indexing object to be used in initializing the implicit vectors

integer(kind=ik), intent(in) :: procRank

Rank of current processors

real(kind=rk), intent(in), optional :: nonlinTol

Picard iteration relative tolerance

real(kind=rk), intent(in), optional :: absTol

Picard iteration absolute tolerance in epsilon units

integer(kind=ik), intent(in), optional :: maxIters

Maximum number of Picard iterations

integer(kind=ik), intent(in), optional, dimension(:) :: convergenceIndices

Variable indices in VariableContainer used to determine convergence of Picard iterations

integer(kind=ik), intent(in), optional, dimension(:) :: modelList

List of models this integrator is responsible for

type(IntArray), intent(in), optional, dimension(:) :: termGroups

Term groups this integrator is responsible for - should conform with modelList

logical, intent(in), optional :: evolvesTimeVar

Set to true if this integrator is allowed to change the "time" varible (if present in passed variable container)

class(TimestepController), intent(in), optional :: dtController

User-supplied timestep controller object

real(kind=rk), intent(in), optional :: initialTimestep

Default timestep

logical, intent(in), optional :: use2Norm

Set to true if this integrator uses the 2-norm instead of the local inf-norm for checking individual variable convergence

integer(kind=ik), intent(in), optional :: petscGroup

PETSc obj group this solver should interact with (defaults to 1)

type(InternalControllerOptions), intent(in), optional :: intContOptions

InternalControlOptions (if present turns on internal control)

character(len=*), intent(in), optional :: integratorName

Name of integrator used in printing

module subroutine integrateBDE(this, manipulatedModeller, outputVars, inputVars)

Integration routine for BDE integrator - implementation of abstract manipulate routine

Arguments

Type IntentOptional Attributes Name
class(PicardBDEIntegrator), intent(inout) :: this
class(ModellerSurrogate), intent(inout) :: manipulatedModeller

Modeller to be used in callbacks during integration

class(VariableContainer), intent(inout) :: outputVars

VariableContainer object to store the integration output

class(VariableContainer), intent(in) :: inputVars

VariableContainer object housing input data for the integration routine

pure module subroutine setNonlinTol(this, nonlinTol)

Setter for nonlinTol value

Arguments

Type IntentOptional Attributes Name
class(PicardBDEIntegrator), intent(inout) :: this
real(kind=rk), intent(in) :: nonlinTol

Picard iteration tolerance

pure module subroutine setMaxIterations(this, maxIter)

Setter for maxIterations value

Arguments

Type IntentOptional Attributes Name
class(PicardBDEIntegrator), intent(inout) :: this
integer(kind=ik), intent(in) :: maxIter

pure module subroutine setConvergenceIndices(this, convVars)

Setter for convergenceTestVars array

Arguments

Type IntentOptional Attributes Name
class(PicardBDEIntegrator), intent(inout) :: this
integer(kind=ik), intent(in), dimension(:) :: convVars