20 #ifndef GEOS_PHYSICSSOLVERS_PHYSICSSOLVERBASE_HPP_
21 #define GEOS_PHYSICSSOLVERS_PHYSICSSOLVERBASE_HPP_
23 #include "codingUtilities/traits.hpp"
30 #include "physicsSolvers/NonlinearSolverParameters.hpp"
32 #include "physicsSolvers/SolverStatistics.hpp"
40 class DomainPartition;
60 Group *
const parent );
121 real64 const eventProgress,
131 real64 const eventProgress,
407 bool const setSparsity =
true );
466 integer const nonlinearIteration,
480 integer const nonlinearIteration,
537 real64 const scalingFactor );
575 real64 const scalingFactor,
825 template<
typename LAMBDA >
828 for(
auto const & target: m_meshTargets )
830 string const meshBodyName = target.first.first;
831 string const meshLevelName = target.first.second;
836 if( meshLevelPtr==
nullptr )
840 lambda( meshBodyName, *meshLevelPtr, regionNames );
851 template<
typename LAMBDA >
854 for(
auto const & target: m_meshTargets )
856 string const meshBodyName = target.first.first;
857 string const meshLevelName = target.first.second;
862 if( meshLevelPtr==
nullptr )
866 lambda( meshBodyName, *meshLevelPtr, regionNames );
884 virtual bool registerCallback(
void * func,
const std::type_info & funcType )
final override;
898 #if defined(GEOS_USE_PYGEOSX)
903 virtual PyTypeObject * getPythonType()
const override;
912 return m_meshTargets;
942 real64 const oldNewtonNorm,
952 template<
typename CONSTITUTIVE_BASE_TYPE >
962 template<
typename CONSTITUTIVE_BASE_TYPE >
982 template<
typename BASETYPE = constitutive::ConstitutiveBase,
typename LOOKUP_TYPE >
986 return constitutiveModels.
getGroup< BASETYPE >( key );
997 template<
typename BASETYPE = constitutive::ConstitutiveBase,
typename LOOKUP_TYPE >
1001 return constitutiveModels.
getGroup< BASETYPE >( key );
1037 std::unique_ptr< PreconditionerBase< LAInterface > >
m_precond;
1061 std::map< std::string, std::chrono::system_clock::duration >
m_timers;
1086 bool solveNonlinearSystem(
real64 const & time_n,
1097 void logEndOfCycleInformation(
integer const cycleNumber,
1099 std::vector< real64 >
const & subStepDt )
const;
1103 template<
typename CONSTITUTIVE_BASE_TYPE >
1111 GEOS_ERROR_IF( !validName.empty(),
"A valid constitutive model was already found." );
1112 validName = model.getName();
1117 template<
typename CONSTITUTIVE_BASE_TYPE >
1125 GEOS_ERROR_IF( !validName.empty(),
"A valid constitutive model was already found." );
1126 validName = model.getName();
#define GEOS_UNUSED_VAR(...)
Mark an unused variable and silence compiler warnings.
#define GEOS_ERROR_IF(EXP, msg)
Conditionally raise a hard error and terminate the program.
The DoFManager is responsible for allocating global dofs, constructing sparsity patterns,...
Partition of the decomposed physical domain. It also manages the connexion information to its neighbo...
dataRepository::Group const & getConstitutiveModels() const
Get the group in which the constitutive models of this subregion are registered.
The class is used to manage mesh body.
Group & getMeshLevels()
Get the meshLevels group.
Class facilitating the representation of a multi-level discretization of a MeshBody.
dataRepository::Group const & getConstitutiveModels() const
Get the group in which the constitutive models of this subregion are registered.
Base class for all physics solvers.
SolverStatistics const & getSolverStatistics() const
const accessor for the solver statistics.
void setSystemSetupTimestamp(Timestamp timestamp)
set the timestamp of the system setup
virtual void registerDataOnMesh(Group &MeshBodies) override
Register wrappers that contain data on the mesh objects.
PhysicsSolverBase(PhysicsSolverBase const &)=delete
Deleted copy constructor.
static BASETYPE const & getConstitutiveModel(dataRepository::Group const &dataGroup, LOOKUP_TYPE const &key)
Get the Constitutive Model object.
virtual string getCatalogName() const =0
DofManager & getDofManager()
Getter for degree-of-freedom manager.
SolverStatistics m_solverStatistics
Solver statistics.
real64 m_cflFactor
Courant–Friedrichs–Lewy factor for the timestep.
string m_discretizationName
name of the FV discretization object in the data repository
map< std::pair< string, string >, array1d< string > > const & getMeshTargets() const
accessor for m_meshTargets
string getDiscretizationName() const
return the name of the discretization object
ParallelVector const & getSystemRhs() const
Getter for system rhs vector.
ParallelMatrix m_matrix
System matrix.
std::map< std::string, std::chrono::system_clock::duration > m_timers
Timers for the aggregate profiling of the solver.
LinearSolverParameters & getLinearSolverParameters()
accessor for the linear solver parameters.
static BASETYPE & getConstitutiveModel(dataRepository::Group &dataGroup, LOOKUP_TYPE const &key)
Get the Constitutive Model object.
void forDiscretizationOnMeshTargets(Group const &meshBodies, LAMBDA &&lambda) const
Loop over the target discretization on all mesh targets and apply callback.
virtual bool checkSequentialSolutionIncrements(DomainPartition &domain) const
Check if the solution increments are ok to use.
PhysicsSolverBase(string const &name, Group *const parent)
Constructor for PhysicsSolverBase.
ParallelMatrix & getSystemMatrix()
Getter for system matrix.
static string getConstitutiveName(ElementSubRegionBase const &subRegion)
Get the Constitutive Name object.
integer m_numTimestepsSinceLastDtCut
Number of cycles since last timestep cut.
virtual void setConstitutiveNamesCallSuper(ElementSubRegionBase &subRegion) const
This function sets constitutive name fields on an ElementSubRegionBase, and calls the base function i...
ParallelVector & getSystemRhs()
Getter for system rhs vector.
virtual void cleanup(real64 const time_n, integer const cycleNumber, integer const eventCounter, real64 const eventProgress, DomainPartition &domain) override
Called as the code exits the main run loop.
std::unique_ptr< PreconditionerBase< LAInterface > > m_precond
Custom preconditioner for the "native" iterative solver.
Timestamp getMeshModificationTimestamp(DomainPartition &domain) const
getter for the timestamp of the mesh modification on the mesh levels
void forDiscretizationOnMeshTargets(Group &meshBodies, LAMBDA &&lambda) const
Loop over the target discretization on all mesh targets and apply callback.
localIndex targetRegionIndex(string const ®ionName) const
Get position of a given region within solver's target region list.
ParallelMatrix const & getSystemMatrix() const
Getter for system rhs vector.
integer m_writeLinearSystem
flag for debug output of matrix, rhs, and solution
PhysicsSolverBase()=delete
Deleted constructor.
LinearSolverParameters const & getLinearSolverParameters() const
const accessor for the linear solver parameters.
LinearSolverResult m_linearSolverResult
Result of the last linear solve.
virtual ~PhysicsSolverBase() override
Destructor for PhysicsSolverBase.
CRSMatrix< real64, globalIndex > m_localMatrix
Local system matrix and rhs.
ParallelVector const & getSystemSolution() const
Getter for system solution vector.
virtual void synchronizeNonlinearSolverParameters()
syncronize the nonlinear solver parameters.
virtual bool registerCallback(void *func, const std::type_info &funcType) final override
function to set the value of m_assemblyCallback
SolverStatistics & getSolverStatistics()
accessor for the solver statistics.
virtual bool execute(real64 const time_n, real64 const dt, integer const cycleNumber, integer const eventCounter, real64 const eventProgress, DomainPartition &domain) override
Main extension point of executable targets.
virtual Group * createChild(string const &childKey, string const &childName) override
creates a child group of of this PhysicsSolverBase instantiation
virtual void saveSequentialIterationState(DomainPartition &domain)
Save the state of the solver for sequential iteration.
PhysicsSolverBase(PhysicsSolverBase &&)=default
Move constructor for PhysicsSolverBase.
real64 getTimestepRequest()
getter for the next timestep size
PhysicsSolverBase & operator=(PhysicsSolverBase const &)=delete
Deleted copy assignment operator.
static CatalogInterface::CatalogType & getCatalog()
Get the singleton catalog for PhysicsSolverBase.
PhysicsSolverBase & operator=(PhysicsSolverBase &&)=delete
Deleted move assignment operator.
NonlinearSolverParameters const & getNonlinearSolverParameters() const
const accessor for the nonlinear solver parameters.
DofManager m_dofManager
Data structure to handle degrees of freedom.
R1Tensor const gravityVector() const
return the value of the gravity vector specified in PhysicsSolverManager
virtual void initialize_postMeshGeneration() override
Initialization tasks after mesh generation is completed.
real64 m_maxStableDt
maximum stable time step
real64 eisenstatWalker(real64 const newNewtonNorm, real64 const oldNewtonNorm, LinearSolverParameters::Krylov const &krylovParams)
Eisenstat-Walker adaptive tolerance.
ParallelVector m_solution
System solution vector.
ParallelVector m_rhs
System right-hand side vector.
std::function< void(CRSMatrix< real64, globalIndex >, array1d< real64 >) > m_assemblyCallback
Callback function for assembly step.
CRSMatrixView< real64 const, globalIndex const > getLocalMatrix() const
Getter for local matrix.
Timestamp getSystemSetupTimestamp() const
getter for the timestamp of the system setup
ParallelVector & getSystemSolution()
Getter for system solution vector.
NonlinearSolverParameters & getNonlinearSolverParameters()
accessor for the nonlinear solver parameters.
Timestamp m_systemSetupTimestamp
Timestamp of the last call to setup system.
LinearSolverParametersInput m_linearSolverParameters
Linear solver parameters.
CRSMatrix< real64, globalIndex > & getLocalMatrix()
Getter for local matrix.
NonlinearSolverParameters m_nonlinearSolverParameters
Nonlinear solver parameters.
real64 m_nextDt
timestep of the next cycle
void generateMeshTargetsFromTargetRegions(Group const &meshBodies)
Generate mesh targets from target regions.
DofManager const & getDofManager() const
Getter for degree-of-freedom manager.
This class is used to log the solver statistics.
This class provides the base class/interface for the catalog value objects.
std::unordered_map< std::string, std::unique_ptr< CatalogInterface< BASETYPE, ARGS... > > > CatalogType
This is the type that will be used for the catalog. The catalog is actually instantiated in the BASET...
T * getGroupPointer(KEY const &key)
Return a pointer to a sub-group of the current Group.
T & getGroup(KEY const &key)
Return a reference to a sub-group of the current Group.
void forSubGroups(LAMBDA &&lambda)
Apply the given functor to subgroups that can be casted to one of specified types.
Base template for ordered and unordered maps.
virtual void updateState(DomainPartition &domain)
Recompute all dependent quantities from primary variables (including constitutive models)
virtual void applySystemSolution(DofManager const &dofManager, arrayView1d< real64 const > const &localSolution, real64 const scalingFactor, real64 const dt, DomainPartition &domain)
Function to apply the solution vector to the state.
virtual bool lineSearch(real64 const &time_n, real64 const &dt, integer const cycleNumber, DomainPartition &domain, DofManager const &dofManager, CRSMatrixView< real64, globalIndex const > const &localMatrix, ParallelVector &rhs, ParallelVector &solution, real64 const scaleFactor, real64 &lastResidual)
Function to perform line search.
virtual real64 setNextDtBasedOnStateChange(real64 const ¤tDt, DomainPartition &domain)
function to set the next dt based on state change
virtual void resetStateToBeginningOfStep(DomainPartition &domain)
reset state of physics back to the beginning of the step.
virtual void assembleSystem(real64 const time, real64 const dt, DomainPartition &domain, DofManager const &dofManager, CRSMatrixView< real64, globalIndex const > const &localMatrix, arrayView1d< real64 > const &localRhs)
function to assemble the linear system matrix and rhs
virtual bool resetConfigurationToDefault(DomainPartition &domain) const
resets the configuration to the default value.
virtual real64 linearImplicitStep(real64 const &time_n, real64 const &dt, integer const cycleNumber, DomainPartition &domain)
Function for a linear implicit integration step.
virtual bool updateConfiguration(DomainPartition &domain)
updates the configuration (if needed) based on the state after a converged Newton loop.
virtual void applyBoundaryConditions(real64 const time, real64 const dt, DomainPartition &domain, DofManager const &dofManager, CRSMatrixView< real64, globalIndex const > const &localMatrix, arrayView1d< real64 > const &localRhs)
apply boundary condition to system
virtual real64 nonlinearImplicitStep(real64 const &time_n, real64 const &dt, integer const cycleNumber, DomainPartition &domain)
Function for a nonlinear implicit integration step.
virtual real64 scalingForSystemSolution(DomainPartition &domain, DofManager const &dofManager, arrayView1d< real64 const > const &localSolution)
Function to determine if the solution vector should be scaled back in order to maintain a known const...
virtual void implicitStepSetup(real64 const &time_n, real64 const &dt, DomainPartition &domain)
function to perform setup for implicit timestep
void debugOutputSolution(real64 const &time, integer const cycleNumber, integer const nonlinearIteration, ParallelVector const &solution) const
Output the linear system solution for debug purposes.
virtual void setupSystem(DomainPartition &domain, DofManager &dofManager, CRSMatrix< real64, globalIndex > &localMatrix, ParallelVector &rhs, ParallelVector &solution, bool const setSparsity=true)
Set up the linear system (DOF indices and sparsity patterns)
virtual bool checkSystemSolution(DomainPartition &domain, DofManager const &dofManager, arrayView1d< real64 const > const &localSolution, real64 const scalingFactor)
Function to check system solution for physical consistency and constraint violation.
virtual real64 calculateResidualNorm(real64 const &time, real64 const &dt, DomainPartition const &domain, DofManager const &dofManager, arrayView1d< real64 const > const &localRhs)
calculate the norm of the global system residual
virtual real64 setNextDtBasedOnCFL(real64 const ¤tDt, DomainPartition &domain)
function to set the next dt based on state change
virtual real64 solverStep(real64 const &time_n, real64 const &dt, integer const cycleNumber, DomainPartition &domain)
entry function to perform a solver step
virtual void setupDofs(DomainPartition const &domain, DofManager &dofManager) const
Populate degree-of-freedom manager with fields relevant to this solver.
virtual void implicitStepComplete(real64 const &time, real64 const &dt, DomainPartition &domain)
perform cleanup for implicit timestep
virtual real64 setNextDt(real64 const ¤tDt, DomainPartition &domain)
function to set the next time step size
virtual void outputConfigurationStatistics(DomainPartition const &domain) const
virtual real64 setNextDtBasedOnNewtonIter(real64 const ¤tDt)
function to set the next time step size based on Newton convergence
void debugOutputSystem(real64 const &time, integer const cycleNumber, integer const nonlinearIteration, ParallelMatrix const &matrix, ParallelVector const &rhs) const
Output the assembled linear system for debug purposes.
virtual void solveLinearSystem(DofManager const &dofManager, ParallelMatrix &matrix, ParallelVector &rhs, ParallelVector &solution)
function to apply a linear system solver to the assembled system.
virtual real64 getTimestepRequest(real64 const) override
getter for the next timestep size
virtual bool lineSearchWithParabolicInterpolation(real64 const &time_n, real64 const &dt, integer const cycleNumber, DomainPartition &domain, DofManager const &dofManager, CRSMatrixView< real64, globalIndex const > const &localMatrix, ParallelVector &rhs, ParallelVector &solution, real64 const scaleFactor, real64 &lastResidual, real64 &residualNormT)
Function to perform line search using a parabolic interpolation to find the scaling factor.
virtual void resetConfigurationToBeginningOfStep(DomainPartition &domain)
resets the configuration to the beginning of the time-step.
virtual real64 explicitStep(real64 const &time_n, real64 const &dt, integer const cycleNumber, DomainPartition &domain)
Entry function for an explicit time integration step.
ArrayView< T, 1 > arrayView1d
Alias for 1D array view.
unsigned long long int Timestamp
Timestamp type (used to perform actions such a sparsity pattern computation after mesh modifications)
LvArray::CRSMatrixView< T, COL_INDEX, localIndex const, LvArray::ChaiBuffer > CRSMatrixView
Alias for CRS Matrix View.
double real64
64-bit floating point type.
GEOS_LOCALINDEX_TYPE localIndex
Local index type (for indexing objects within an MPI partition).
std::int32_t integer
Signed integer type.
LvArray::CRSMatrix< T, COL_INDEX, localIndex, LvArray::ChaiBuffer > CRSMatrix
Alias for CRS Matrix class.
LAInterface::ParallelMatrix ParallelMatrix
Alias for ParallelMatrix.
Array< T, 1 > array1d
Alias for 1D array.
LAInterface::ParallelVector ParallelVector
Alias for ParallelVector.
static constexpr auto constitutiveModelsString()
Krylov-method parameters.
Set of parameters for a linear solver or preconditioner.
Results/stats of a linear solve.
static constexpr char const * baseDiscretizationString()
Structure to hold scoped key names.
static constexpr char const * nonlinearSolverParametersString()
static constexpr char const * solverStatisticsString()
static constexpr char const * linearSolverParametersString()
Structure to hold scoped key names.
static constexpr char const * cflFactorString()
static constexpr char const * discretizationString()
static constexpr char const * writeLinearSystemString()
static constexpr char const * targetRegionsString()
static constexpr char const * maxStableDtString()
static constexpr char const * minDtIncreaseIntervalString()
static constexpr char const * initialDtString()
static constexpr char const * meshTargetsString()