v0.9.0
Public Types | Public Member Functions | Public Attributes | List of all members
MoFEM::TSMethod Struct Reference

data structure for TS (time stepping) contextStructure stores context data which are set in functions run by PETSc Time Stepping functions. 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
}
 

Public Member Functions

MoFEMErrorCode query_interface (const MOFEMuuid &uuid, UnknownInterface **iface) const
 
 TSMethod ()
 
virtual ~TSMethod ()
 
MoFEMErrorCode setTsCtx (const TSContext &ctx)
 Set Ts context. More...
 
MoFEMErrorCode copyTs (const TSMethod &ts)
 Copy TS solver data. More...
 
MoFEMErrorCode setTs (TS _ts)
 Set TS solver. More...
 
- Public Member Functions inherited from MoFEM::UnknownInterface
template<class IFACE >
MoFEMErrorCode registerInterface (const MOFEMuuid &uuid, bool error_if_registration_failed=true)
 Register interface. More...
 
template<class IFACE , bool VERIFY = false>
MoFEMErrorCode getInterface (const MOFEMuuid &uuid, IFACE *&iface) const
 Get interface by uuid and return reference to pointer of interface. More...
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE *&iface) const
 Get interface refernce to pointer of interface. More...
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE **const iface) const
 Get interface pointer to pointer of interface. More...
 
template<class IFACE , typename boost::enable_if< boost::is_pointer< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get interface pointer to pointer of interface. More...
 
template<class IFACE , typename boost::enable_if< boost::is_reference< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get reference to interface. More...
 
template<class IFACE >
IFACE * getInterface () const
 Function returning pointer to interface. More...
 
virtual ~UnknownInterface ()=default
 
virtual MoFEMErrorCode getLibVersion (Version &version) const
 Get library version. More...
 
virtual const MoFEMErrorCode getFileVersion (moab::Interface &moab, Version &version) const
 Get database major version. More...
 
virtual MoFEMErrorCode getInterfaceVersion (Version &version) const
 Get database major version. More...
 
template<>
MoFEMErrorCode getInterface (const MOFEMuuid &uuid, UnknownInterface *&iface) const
 

Public Attributes

TSContext ts_ctx
 
TS ts
 time solver More...
 
Vec ts_u
 state vector More...
 
Vec ts_u_t
 time derivative of state vector More...
 
Vec ts_u_tt
 second time derivative of state vector More...
 
Vec ts_F
 residual vector More...
 
Mat ts_A
 
Mat ts_B
 Preconditioner for ts_A. More...
 
PetscInt ts_step
 time step More...
 
PetscReal ts_a
 shift for U_tt (see PETSc Time Solver) More...
 
PetscReal ts_v
 shift for U_t shift for U_t More...
 
PetscReal ts_t
 time More...
 

Additional Inherited Members

- Protected Member Functions inherited from MoFEM::UnknownInterface
boost::typeindex::type_index getClassIdx (const MOFEMuuid &uid) const
 Get type name for interface Id. More...
 
MOFEMuuid getUId (const boost::typeindex::type_index &class_idx) const
 Get interface Id for class name. More...
 

Detailed Description

data structure for TS (time stepping) context

Structure stores context data which are set in functions run by PETSc Time Stepping functions.

Definition at line 154 of file LoopMethods.hpp.

Member Enumeration Documentation

◆ TSContext

Enumerator
CTX_TSSETRHSFUNCTION 
CTX_TSSETRHSJACOBIAN 
CTX_TSSETIFUNCTION 
CTX_TSSETIJACOBIAN 
CTX_TSTSMONITORSET 
CTX_TSNONE 

Definition at line 165 of file LoopMethods.hpp.

Constructor & Destructor Documentation

◆ TSMethod()

MoFEM::TSMethod::TSMethod ( )

Definition at line 175 of file LoopMethods.hpp.

176  : ts_ctx(CTX_TSNONE), ts_u(PETSC_NULL), ts_u_t(PETSC_NULL),
177  ts_F(PETSC_NULL), ts_A(PETSC_NULL), ts_B(PETSC_NULL), ts_step(-1),
178  ts_a(0), ts_t(0) {}
Vec ts_u_t
time derivative of state vector
PetscReal ts_t
time
PetscReal ts_a
shift for U_tt (see PETSc Time Solver)
Vec ts_F
residual vector
PetscInt ts_step
time step
Vec ts_u
state vector
Mat ts_B
Preconditioner for ts_A.
TSContext ts_ctx

◆ ~TSMethod()

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

Definition at line 180 of file LoopMethods.hpp.

180 {}

Member Function Documentation

◆ copyTs()

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

Copy TS solver data.

Definition at line 75 of file LoopMethods.cpp.

75  {
77  this->ts_ctx = ts.ts_ctx;
78  this->ts = ts.ts;
79  this->ts_u = ts.ts_u;
80  this->ts_u_t = ts.ts_u_t;
81  this->ts_F = ts.ts_F;
82  this->ts_A = ts.ts_A;
83  this->ts_B = ts.ts_B;
84  this->ts_step = ts.ts_step;
85  this->ts_a = ts.ts_a;
86  this->ts_t = ts.ts_t;
88 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:501
Vec ts_u_t
time derivative of state vector
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:508
PetscReal ts_t
time
PetscReal ts_a
shift for U_tt (see PETSc Time Solver)
Vec ts_F
residual vector
PetscInt ts_step
time step
TS ts
time solver
Vec ts_u
state vector
Mat ts_B
Preconditioner for ts_A.
TSContext ts_ctx

◆ query_interface()

MoFEMErrorCode MoFEM::TSMethod::query_interface ( const MOFEMuuid uuid,
UnknownInterface **  iface 
) const
virtual

Implements MoFEM::UnknownInterface.

Reimplemented in MoFEM::DofMethod, MoFEM::EntityMethod, MoFEM::FEMethod, and MoFEM::BasicMethod.

Definition at line 156 of file LoopMethods.hpp.

157  {
158  if (uuid == IDD_MOFEMTsMethod) {
159  *iface = const_cast<TSMethod *>(this);
161  }
162  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "unknown interface");
163  }
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:508
static const MOFEMuuid IDD_MOFEMTsMethod
Definition: LoopMethods.hpp:28

◆ setTs()

MoFEMErrorCode MoFEM::TSMethod::setTs ( TS  _ts)

Set TS solver.

Definition at line 70 of file LoopMethods.cpp.

70  {
72  ts = _ts;
74 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:501
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:508
TS ts
time solver

◆ setTsCtx()

MoFEMErrorCode MoFEM::TSMethod::setTsCtx ( const TSContext ctx)

Set Ts context.

Definition at line 65 of file LoopMethods.cpp.

65  {
67  ts_ctx = ctx;
69 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:501
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:508
TSContext ts_ctx

Member Data Documentation

◆ ts

TS MoFEM::TSMethod::ts

time solver

Examples
UnsaturatedFlow.hpp.

Definition at line 191 of file LoopMethods.hpp.

◆ ts_A

Mat MoFEM::TSMethod::ts_A

Jacobian of G(U) = F(t,U,W+v*U,W'+a*U), equivalent to dF/dU + v*dF/dU_t + a*dF/dU_tt

Definition at line 197 of file LoopMethods.hpp.

◆ ts_a

PetscReal MoFEM::TSMethod::ts_a

shift for U_tt (see PETSc Time Solver)

Examples
UnsaturatedFlow.hpp.

Definition at line 202 of file LoopMethods.hpp.

◆ ts_B

Mat MoFEM::TSMethod::ts_B

Preconditioner for ts_A.

Examples
UnsaturatedFlow.hpp.

Definition at line 199 of file LoopMethods.hpp.

◆ ts_ctx

TSContext MoFEM::TSMethod::ts_ctx

Definition at line 174 of file LoopMethods.hpp.

◆ ts_F

Vec MoFEM::TSMethod::ts_F

residual vector

Examples
UnsaturatedFlow.hpp.

Definition at line 195 of file LoopMethods.hpp.

◆ ts_step

PetscInt MoFEM::TSMethod::ts_step

time step

Definition at line 201 of file LoopMethods.hpp.

◆ ts_t

PetscReal MoFEM::TSMethod::ts_t

time

Examples
UnsaturatedFlow.hpp.

Definition at line 204 of file LoopMethods.hpp.

◆ ts_u

Vec MoFEM::TSMethod::ts_u

state vector

Definition at line 192 of file LoopMethods.hpp.

◆ ts_u_t

Vec MoFEM::TSMethod::ts_u_t

time derivative of state vector

Definition at line 193 of file LoopMethods.hpp.

◆ ts_u_tt

Vec MoFEM::TSMethod::ts_u_tt

second time derivative of state vector

Definition at line 194 of file LoopMethods.hpp.

◆ ts_v

PetscReal MoFEM::TSMethod::ts_v

shift for U_t shift for U_t

Definition at line 203 of file LoopMethods.hpp.


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