v0.15.0
Loading...
Searching...
No Matches
Public Types | Public Member Functions | Public Attributes | List of all members
MoFEM::TSMethod Struct Reference

Data structure for TS (time stepping) context. More...

#include "src/interfaces/LoopMethods.hpp"

Inheritance diagram for MoFEM::TSMethod:
[legend]
Collaboration diagram for MoFEM::TSMethod:
[legend]

Public Types

enum  TSContext {
  CTX_TSSETRHSFUNCTION , CTX_TSSETRHSJACOBIAN , CTX_TSSETIFUNCTION , CTX_TSSETIJACOBIAN ,
  CTX_TSTSMONITORSET , CTX_TSNONE
}
 Context enumeration for TS solver phases. More...
 

Public Member Functions

MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 Query interface for type casting.
 
 TSMethod ()
 Default constructor.
 
virtual ~TSMethod ()=default
 Virtual destructor.
 
MoFEMErrorCode copyTs (const TSMethod &ts)
 Copy TS solver data from another instance.
 

Public Attributes

TS ts
 PETSc time stepping solver object.
 
TSContext ts_ctx
 Current TS computation context.
 
PetscInt ts_step
 Current time step number.
 
PetscReal ts_a
 Shift parameter for U_t (see PETSc Time Solver documentation)
 
PetscReal ts_aa
 Shift parameter for U_tt (second time derivative)
 
PetscReal ts_t
 Current time value.
 
PetscReal ts_dt
 Current time step size.
 
Vec & ts_u
 Reference to current state vector U(t)
 
Vec & ts_u_t
 Reference to first time derivative of state vector dU/dt.
 
Vec & ts_u_tt
 Reference to second time derivative of state vector d²U/dt²
 
Vec & ts_F
 Reference to residual vector F(t,U,U_t,U_tt)
 
Mat & ts_A
 Reference to Jacobian matrix: dF/dU + a*dF/dU_t + aa*dF/dU_tt.
 
Mat & ts_B
 Reference to preconditioner matrix for ts_A.
 

Additional Inherited Members

- Protected Types inherited from MoFEM::PetscData
enum  DataContext {
  CTX_SET_NONE = 0 , CTX_SET_F = 1 << 0 , CTX_SET_A = 1 << 1 , CTX_SET_B = 1 << 2 ,
  CTX_SET_X = 1 << 3 , CTX_SET_DX = 1 << 4 , CTX_SET_X_T = 1 << 5 , CTX_SET_X_TT = 1 << 6 ,
  CTX_SET_TIME = 1 << 7
}
 Enumeration for data context flags. More...
 
using Switches = std::bitset< 8 >
 Bitset type for context switches.
 
- Protected Member Functions inherited from MoFEM::PetscData
 PetscData ()
 Default constructor.
 
virtual ~PetscData ()=default
 Virtual destructor.
 
MoFEMErrorCode copyPetscData (const PetscData &petsc_data)
 Copy PETSc data from another instance.
 
- Protected Member Functions inherited from MoFEM::UnknownInterface
template<class IFACE >
MoFEMErrorCode registerInterface (bool error_if_registration_failed=true)
 Register interface.
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE *&iface) const
 Get interface reference to pointer of interface.
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE **const iface) const
 Get interface pointer to pointer of interface.
 
template<class IFACE , typename boost::enable_if< boost::is_pointer< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get interface pointer to pointer of interface.
 
template<class IFACE , typename boost::enable_if< boost::is_reference< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get reference to interface.
 
template<class IFACE >
IFACE * getInterface () const
 Function returning pointer to interface.
 
virtual ~UnknownInterface ()=default
 
- Static Protected Member Functions inherited from MoFEM::UnknownInterface
static MoFEMErrorCode getLibVersion (Version &version)
 Get library version.
 
static MoFEMErrorCode getFileVersion (moab::Interface &moab, Version &version)
 Get database major version.
 
static MoFEMErrorCode setFileVersion (moab::Interface &moab, Version version=Version(MoFEM_VERSION_MAJOR, MoFEM_VERSION_MINOR, MoFEM_VERSION_BUILD))
 Get database major version.
 
static MoFEMErrorCode getInterfaceVersion (Version &version)
 Get database major version.
 
- Protected Attributes inherited from MoFEM::PetscData
Switches data_ctx
 Current data context switches.
 
Vec f
 PETSc residual vector.
 
Mat A
 PETSc Jacobian matrix.
 
Mat B
 PETSc preconditioner matrix.
 
Vec x
 PETSc solution vector.
 
Vec dx
 PETSc solution increment vector.
 
Vec x_t
 PETSc first time derivative vector.
 
Vec x_tt
 PETSc second time derivative vector.
 
- Static Protected Attributes inherited from MoFEM::PetscData
static constexpr Switches CtxSetNone = PetscData::Switches(CTX_SET_NONE)
 No data switch.
 
static constexpr Switches CtxSetF = PetscData::Switches(CTX_SET_F)
 Residual vector switch.
 
static constexpr Switches CtxSetA = PetscData::Switches(CTX_SET_A)
 Jacobian matrix switch.
 
static constexpr Switches CtxSetB = PetscData::Switches(CTX_SET_B)
 Preconditioner matrix switch.
 
static constexpr Switches CtxSetX = PetscData::Switches(CTX_SET_X)
 Solution vector switch.
 
static constexpr Switches CtxSetDX = PetscData::Switches(CTX_SET_DX)
 Solution increment switch.
 
static constexpr Switches CtxSetX_T = PetscData::Switches(CTX_SET_X_T)
 First time derivative switch.
 
static constexpr Switches CtxSetX_TT = PetscData::Switches(CTX_SET_X_TT)
 Second time derivative switch.
 
static constexpr Switches CtxSetTime = PetscData::Switches(CTX_SET_TIME)
 Time value switch.
 

Detailed Description

Data structure for TS (time stepping) context.

This structure extends PetscData to provide context and data management specifically for PETSc TS (Time Stepping) solvers. It manages time-dependent computations including time derivatives, time step information, and temporal integration contexts for transient finite element analyses.

Definition at line 261 of file LoopMethods.hpp.

Member Enumeration Documentation

◆ TSContext

Context enumeration for TS solver phases.

Indicates the current phase of time stepping computation, determining which temporal operations are being performed and what data is active.

Enumerator
CTX_TSSETRHSFUNCTION 

Setting up right-hand side function.

CTX_TSSETRHSJACOBIAN 

Setting up right-hand side Jacobian.

CTX_TSSETIFUNCTION 

Setting up implicit function.

CTX_TSSETIJACOBIAN 

Setting up implicit Jacobian.

CTX_TSTSMONITORSET 

Setting up time step monitoring.

CTX_TSNONE 

No specific TS context.

Definition at line 278 of file LoopMethods.hpp.

278 {
279 CTX_TSSETRHSFUNCTION, ///< Setting up right-hand side function
280 CTX_TSSETRHSJACOBIAN, ///< Setting up right-hand side Jacobian
281 CTX_TSSETIFUNCTION, ///< Setting up implicit function
282 CTX_TSSETIJACOBIAN, ///< Setting up implicit Jacobian
283 CTX_TSTSMONITORSET, ///< Setting up time step monitoring
284 CTX_TSNONE ///< No specific TS context
285 };
@ CTX_TSSETIFUNCTION
Setting up implicit function.
@ CTX_TSSETRHSFUNCTION
Setting up right-hand side function.
@ CTX_TSSETRHSJACOBIAN
Setting up right-hand side Jacobian.
@ CTX_TSTSMONITORSET
Setting up time step monitoring.
@ CTX_TSSETIJACOBIAN
Setting up implicit Jacobian.
@ CTX_TSNONE
No specific TS context.

Constructor & Destructor Documentation

◆ TSMethod()

MoFEM::TSMethod::TSMethod ( )

Default constructor.

Definition at line 112 of file LoopMethods.cpp.

Vec f
PETSc residual vector.
Vec x_t
PETSc first time derivative vector.
Vec x_tt
PETSc second time derivative vector.
Vec x
PETSc solution vector.
Mat B
PETSc preconditioner matrix.
Mat A
PETSc Jacobian matrix.
PetscReal ts_t
Current time value.
Mat & ts_A
Reference to Jacobian matrix: dF/dU + a*dF/dU_t + aa*dF/dU_tt.
Vec & ts_F
Reference to residual vector F(t,U,U_t,U_tt)
Vec & ts_u_tt
Reference to second time derivative of state vector d²U/dt²
TSContext ts_ctx
Current TS computation context.
PetscReal ts_a
Shift parameter for U_t (see PETSc Time Solver documentation)
Vec & ts_u_t
Reference to first time derivative of state vector dU/dt.
PetscInt ts_step
Current time step number.
Vec & ts_u
Reference to current state vector U(t)
Mat & ts_B
Reference to preconditioner matrix for ts_A.

◆ ~TSMethod()

virtual MoFEM::TSMethod::~TSMethod ( )
virtualdefault

Virtual destructor.

Member Function Documentation

◆ copyTs()

MoFEMErrorCode MoFEM::TSMethod::copyTs ( const TSMethod ts)

Copy TS solver data from another instance.

Parameters
tsSource TS method
Returns
Error code

Definition at line 117 of file LoopMethods.cpp.

117 {
120 this->ts_ctx = ts.ts_ctx;
121 this->ts = ts.ts;
122 this->ts_step = ts.ts_step;
123 this->ts_a = ts.ts_a;
124 this->ts_t = ts.ts_t;
125 this->ts_dt = ts.ts_dt;
127}
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
MoFEMErrorCode copyPetscData(const PetscData &petsc_data)
Copy PETSc data from another instance.
TS ts
PETSc time stepping solver object.
PetscReal ts_dt
Current time step size.

◆ query_interface()

MoFEMErrorCode MoFEM::TSMethod::query_interface ( boost::typeindex::type_index  type_index,
UnknownInterface **  iface 
) const
virtual

Query interface for type casting.

Parameters
type_indexType information for interface querying
ifacePointer to interface object
Returns
Error code

Reimplemented from MoFEM::PetscData.

Definition at line 106 of file LoopMethods.cpp.

107 {
108 *iface = const_cast<TSMethod *>(this);
109 return 0;
110}
TSMethod()
Default constructor.

Member Data Documentation

◆ ts

TS MoFEM::TSMethod::ts

◆ ts_a

PetscReal MoFEM::TSMethod::ts_a

Shift parameter for U_t (see PETSc Time Solver documentation)

Examples
mofem/tutorials/cor-0to1/src/UnsaturatedFlow.hpp, and mofem/tutorials/cor-9/reaction_diffusion.cpp.

Definition at line 309 of file LoopMethods.hpp.

◆ ts_A

Mat& MoFEM::TSMethod::ts_A

Reference to Jacobian matrix: dF/dU + a*dF/dU_t + aa*dF/dU_tt.

Definition at line 319 of file LoopMethods.hpp.

◆ ts_aa

PetscReal MoFEM::TSMethod::ts_aa

Shift parameter for U_tt (second time derivative)

Definition at line 310 of file LoopMethods.hpp.

◆ ts_B

Mat& MoFEM::TSMethod::ts_B

Reference to preconditioner matrix for ts_A.

Examples
mofem/tutorials/cor-0to1/src/UnsaturatedFlow.hpp.

Definition at line 320 of file LoopMethods.hpp.

◆ ts_ctx

TSContext MoFEM::TSMethod::ts_ctx

Current TS computation context.

Examples
mofem/users_modules/basic_finite_elements/nonlinear_elasticity/nonlinear_dynamics.cpp.

Definition at line 306 of file LoopMethods.hpp.

◆ ts_dt

PetscReal MoFEM::TSMethod::ts_dt

Current time step size.

Definition at line 312 of file LoopMethods.hpp.

◆ ts_F

Vec& MoFEM::TSMethod::ts_F

Reference to residual vector F(t,U,U_t,U_tt)

Examples
mofem/tutorials/cor-0to1/src/UnsaturatedFlow.hpp.

Definition at line 317 of file LoopMethods.hpp.

◆ ts_step

PetscInt MoFEM::TSMethod::ts_step

◆ ts_t

PetscReal MoFEM::TSMethod::ts_t

◆ ts_u

Vec& MoFEM::TSMethod::ts_u

Reference to current state vector U(t)

Examples
mofem/tutorials/vec-4/shallow_wave.cpp.

Definition at line 314 of file LoopMethods.hpp.

◆ ts_u_t

Vec& MoFEM::TSMethod::ts_u_t

Reference to first time derivative of state vector dU/dt.

Definition at line 315 of file LoopMethods.hpp.

◆ ts_u_tt

Vec& MoFEM::TSMethod::ts_u_tt

Reference to second time derivative of state vector d²U/dt²

Definition at line 316 of file LoopMethods.hpp.


The documentation for this struct was generated from the following files: