v0.9.0
Public Types | Public Member Functions | Public Attributes | Protected Member Functions | Protected Attributes | Private Member Functions | Private Attributes | Friends | List of all members
MoFEM::ForcesAndSourcesCore::UserDataOperator Struct Reference

Data operator to do calculations at integration points.Is inherited and implemented by user to do calculations. It can be used in many different ways but typically is used to integrate matrices (f.e. stiffness matrix) and the right hand vector (f.e. force vector). More...

#include <src/finite_elements/ForcesAndSourcesCore.hpp>

Inheritance diagram for MoFEM::ForcesAndSourcesCore::UserDataOperator:
[legend]
Collaboration diagram for MoFEM::ForcesAndSourcesCore::UserDataOperator:
[legend]

Public Types

enum  OpType { OPROW = 1 << 0, OPCOL = 1 << 1, OPROWCOL = 1 << 2, OPLAST = 1 << 3 }
 Controls loop over entities on element. More...
 

Public Member Functions

 UserDataOperator (const FieldSpace space, const char type=OPLAST, const bool symm=true)
 
 UserDataOperator (const std::string &field_name, const char type, const bool symm=true)
 
 UserDataOperator (const std::string &row_field_name, const std::string &col_field_name, const char type, const bool symm=true)
 
boost::shared_ptr< const NumeredEntFiniteElementgetNumeredEntFiniteElementPtr () const
 Return raw pointer to NumeredEntFiniteElement. More...
 
EntityHandle getFEEntityHandle () const
 Return finite element entity handle. More...
 
boost::weak_ptr< SideNumber > getSideNumberPtr (const int side_number, const EntityType type)
 Get the side number pointer. More...
 
EntityHandle getSideEntity (const int side_number, const EntityType type)
 Get the side entity. More...
 
int getNumberOfNodesOnElement ()
 Get the number of nodes on finite element. More...
 
MoFEMErrorCode getProblemRowIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get row indices. More...
 
MoFEMErrorCode getProblemColIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get col indices. More...
 
const FEMethodgetFEMethod () const
 Return raw pointer to Finite Element Method object. More...
 
int getOpType () const
 Get operator types. More...
 
void setOpType (const OpType type)
 Set operator type. More...
 
void addOpType (const OpType type)
 Add operator type. More...
 
int getNinTheLoop () const
 get number of finite element in the loop More...
 
int getLoopSize () const
 get size of elements in the loop More...
 
const std::string & getFEName () const
 Get name of the element. More...
 
Accessing SNES
Vec getSnesF () const
 
Vec getSnesX () const
 
Mat getSnesA () const
 
Mat getSnesB () const
 
Accessing TS
Vec getTSu () const
 
Vec getTSu_t () const
 
Vec getTSf () const
 
Mat getTSA () const
 
Mat getTSB () const
 
int getTSstep () const
 
double getTStime () const
 
double getTSa () const
 
Base funtions and integration points
MatrixDoublegetGaussPts ()
 matrix of integration (Gauss) points for Volume Element More...
 
auto getFTensor0IntegrationWeight ()
 Get integration weights. More...
 
Deprecated (do not use)
DEPRECATED MoFEMErrorCode getPorblemRowIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 
- Public Member Functions inherited from MoFEM::DataOperator
 DataOperator (const bool symm=true, const bool do_vertices=true, const bool do_edges=true, const bool do_quads=true, const bool do_tris=true, const bool do_tets=true, const bool do_prisms=true)
 
virtual ~DataOperator ()
 
virtual MoFEMErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, DataForcesAndSourcesCore::EntData &row_data, DataForcesAndSourcesCore::EntData &col_data)
 Operator for bi-linear form, usually to calculate values on left hand side. More...
 
virtual MoFEMErrorCode opLhs (DataForcesAndSourcesCore &row_data, DataForcesAndSourcesCore &col_data, bool symm=true)
 
virtual MoFEMErrorCode opLhs (DataForcesAndSourcesCore &row_data, DataForcesAndSourcesCore &col_data)
 
virtual MoFEMErrorCode doWork (int side, EntityType type, DataForcesAndSourcesCore::EntData &data)
 Operator for linear form, usually to calculate values on right hand side. More...
 
virtual MoFEMErrorCode opRhs (DataForcesAndSourcesCore &data, const bool do_vertices, const bool do_edges, const bool do_quads, const bool do_tris, const bool do_tets, const bool do_prisms, const bool error_if_no_base=true)
 
virtual MoFEMErrorCode opRhs (DataForcesAndSourcesCore &data, const bool error_if_no_base=true)
 
bool getSymm () const
 Get if operator uses symmetry of DOFs or not. More...
 
void setSymm ()
 set if operator is executed taking in account symmetry More...
 
void unSetSymm ()
 unset if operator is executed for non symmetric problem More...
 

Public Attributes

char opType
 
std::string rowFieldName
 
std::string colFieldName
 
FieldSpace sPace
 
- Public Attributes inherited from MoFEM::DataOperator
bool sYmm
 If true assume that matrix is symmetric structure. More...
 
bool doVertices
 If false skip vertices. More...
 
bool doEdges
 If false skip edges. More...
 
bool doQuads
 
bool doTris
 
bool doTets
 
bool doPrisms
 

Protected Member Functions

MoFEMErrorCode setPtrFE (ForcesAndSourcesCore *ptr)
 
ForcesAndSourcesCoregetPtrFE () const
 
ForcesAndSourcesCoregetSidePtrFE () const
 

Protected Attributes

ForcesAndSourcesCoreptrFE
 

Private Member Functions

MoFEMErrorCode loopSide (const string &fe_name, ForcesAndSourcesCore *side_fe, const size_t dim)
 User call this function to loop over elements on the side of face. This function calls finite element with is operator to do calculations. More...
 

Private Attributes

friend ForcesAndSourcesCore
 

Friends

class EdgeElementForcesAndSourcesCoreBase
 
class FaceElementForcesAndSourcesCoreBase
 

Detailed Description

Data operator to do calculations at integration points.

Is inherited and implemented by user to do calculations. It can be used in many different ways but typically is used to integrate matrices (f.e. stiffness matrix) and the right hand vector (f.e. force vector).

Note: It is assumed that operator is executed for symmetric problem. That means that is executed for not repeating entities on finite element. For example on triangle we have nodes, 3 edges and one face. Because of symmetry calculations are for: nodes-nodes, nodes-edge0, nodes-edge_1, nodes-edge_2, nodes-face, edges_1-edges_1, edges_1-edge_1, edge_1-edge_2, edge_1-edge_1, edge_1-edge_2, edge_2-edge_2, edge_1-face, edge_1-face, edges_3-face, face - face

In case of non symmetric problem in addition calculations of lower off diagonal terms. F.e. edge_1-edge_0, esges_3-edge_0, edge_3-edge_1,

In that case class variable UserDataOperator::sYmm = false;

NoteL: By default sYmm is set for symmetric problems

Examples
hello_world.cpp, and mesh_smoothing.cpp.

Definition at line 84 of file ForcesAndSourcesCore.hpp.

Member Enumeration Documentation

◆ OpType

Controls loop over entities on element.

OPROW is used if row vector is assembled OPCOL is usually used if column vector is assembled OPROWCOL is usually used for assemble matrices.

For typical problem like Bubnov-Galerkin OPROW and OPCOL are the same. In more general case for example for non-square matrices columns and rows could have different numeration and/or different set of fields.

Enumerator
OPROW 
OPCOL 
OPROWCOL 
OPLAST 

Definition at line 98 of file ForcesAndSourcesCore.hpp.

Constructor & Destructor Documentation

◆ UserDataOperator() [1/3]

MoFEM::ForcesAndSourcesCore::UserDataOperator::UserDataOperator ( const FieldSpace  space,
const char  type = OPLAST,
const bool  symm = true 
)

This Constructor is used typically when some modification base shape functions on some approx. space is applied. Operator is run for all data on space.

User has no access to field data from this operator.

Examples
ElasticityMixedFormulation.hpp, and SmallStrainPlasticity.hpp.

Definition at line 1386 of file ForcesAndSourcesCore.cpp.

1389  : DataOperator(symm), opType(type), sPace(space), ptrFE(nullptr) {}
DataOperator(const bool symm=true, const bool do_vertices=true, const bool do_edges=true, const bool do_quads=true, const bool do_tris=true, const bool do_tets=true, const bool do_prisms=true)

◆ UserDataOperator() [2/3]

MoFEM::ForcesAndSourcesCore::UserDataOperator::UserDataOperator ( const std::string &  field_name,
const char  type,
const bool  symm = true 
)

Definition at line 1391 of file ForcesAndSourcesCore.cpp.

1393  : DataOperator(symm), opType(type), rowFieldName(field_name),
1394  colFieldName(field_name), sPace(LASTSPACE), ptrFE(nullptr) {}
FieldSpace in [ 0, LASTSPACE )
Definition: definitions.h:175
DataOperator(const bool symm=true, const bool do_vertices=true, const bool do_edges=true, const bool do_quads=true, const bool do_tris=true, const bool do_tets=true, const bool do_prisms=true)

◆ UserDataOperator() [3/3]

MoFEM::ForcesAndSourcesCore::UserDataOperator::UserDataOperator ( const std::string &  row_field_name,
const std::string &  col_field_name,
const char  type,
const bool  symm = true 
)

Definition at line 1396 of file ForcesAndSourcesCore.cpp.

1399  : DataOperator(symm), opType(type), rowFieldName(row_field_name),
1400  colFieldName(col_field_name), sPace(LASTSPACE), ptrFE(nullptr) {}
FieldSpace in [ 0, LASTSPACE )
Definition: definitions.h:175
DataOperator(const bool symm=true, const bool do_vertices=true, const bool do_edges=true, const bool do_quads=true, const bool do_tris=true, const bool do_tets=true, const bool do_prisms=true)

Member Function Documentation

◆ addOpType()

void MoFEM::ForcesAndSourcesCore::UserDataOperator::addOpType ( const OpType  type)

Add operator type.

Definition at line 884 of file ForcesAndSourcesCore.hpp.

◆ getFEEntityHandle()

EntityHandle MoFEM::ForcesAndSourcesCore::UserDataOperator::getFEEntityHandle ( ) const

Return finite element entity handle.

Returns
Finite element entity handle
Examples
UnsaturatedFlow.hpp.

Definition at line 842 of file ForcesAndSourcesCore.hpp.

842  {
843  return getNumeredEntFiniteElementPtr()->getEnt();
844 }
boost::shared_ptr< const NumeredEntFiniteElement > getNumeredEntFiniteElementPtr() const
Return raw pointer to NumeredEntFiniteElement.

◆ getFEMethod()

const FEMethod * MoFEM::ForcesAndSourcesCore::UserDataOperator::getFEMethod ( ) const

Return raw pointer to Finite Element Method object.

Examples
UnsaturatedFlow.hpp.

Definition at line 874 of file ForcesAndSourcesCore.hpp.

874  {
875  return ptrFE;
876 }

◆ getFEName()

const std::string & MoFEM::ForcesAndSourcesCore::UserDataOperator::getFEName ( ) const

Get name of the element.

Definition at line 896 of file ForcesAndSourcesCore.hpp.

896  {
897  return getFEMethod()->feName;
898 }
std::string feName
Name of finite element.
const FEMethod * getFEMethod() const
Return raw pointer to Finite Element Method object.

◆ getFTensor0IntegrationWeight()

auto MoFEM::ForcesAndSourcesCore::UserDataOperator::getFTensor0IntegrationWeight ( )

Get integration weights.

for(int gg = 0; gg!=getGaussPts.size2(); ++gg) {
// integrate something
++t_w;
}
Returns
FTensor::Tensor0<FTensor::PackPtr<double *, 1>>
Examples
HookeElement.cpp, prism_polynomial_approximation.cpp, quad_polynomial_approximation.cpp, and UnsaturatedFlow.hpp.

Definition at line 952 of file ForcesAndSourcesCore.hpp.

952  {
954  &(getGaussPts()(getGaussPts().size1() - 1, 0)));
955 }
MatrixDouble & getGaussPts()
matrix of integration (Gauss) points for Volume Element

◆ getGaussPts()

MatrixDouble & MoFEM::ForcesAndSourcesCore::UserDataOperator::getGaussPts ( )

matrix of integration (Gauss) points for Volume Element

For triangle: columns 0,1 are x,y coordinates respectively and column 2 is a weight value for example getGaussPts()(1,13) returns y coordinate of 13th Gauss point on particular volume element

For tetrahedron: columns 0,1,2 are x,y,z coordinates respectively and column 3 is a weight value for example getGaussPts()(1,13) returns y coordinate of 13th Gauss point on particular volume element

Examples
HookeElement.cpp, MagneticElement.hpp, and UnsaturatedFlow.hpp.

Definition at line 948 of file ForcesAndSourcesCore.hpp.

948  {
949  return static_cast<ForcesAndSourcesCore *>(ptrFE)->gaussPts;
950 }

◆ getLoopSize()

int MoFEM::ForcesAndSourcesCore::UserDataOperator::getLoopSize ( ) const

get size of elements in the loop

Returns
loop size

Definition at line 892 of file ForcesAndSourcesCore.hpp.

892  {
893  return getFEMethod()->getLoopSize();
894 }
int getLoopSize() const
get loop size
const FEMethod * getFEMethod() const
Return raw pointer to Finite Element Method object.

◆ getNinTheLoop()

int MoFEM::ForcesAndSourcesCore::UserDataOperator::getNinTheLoop ( ) const

get number of finite element in the loop

Returns
number of finite element

Definition at line 888 of file ForcesAndSourcesCore.hpp.

888  {
889  return getFEMethod()->getNinTheLoop();
890 }
int getNinTheLoop() const
get number of evaluated element in the loop
const FEMethod * getFEMethod() const
Return raw pointer to Finite Element Method object.

◆ getNumberOfNodesOnElement()

int MoFEM::ForcesAndSourcesCore::UserDataOperator::getNumberOfNodesOnElement ( )

Get the number of nodes on finite element.

Returns
int

Definition at line 868 of file ForcesAndSourcesCore.hpp.

868  {
869  int num_nodes;
870  CHKERR ptrFE->getNumberOfNodes(num_nodes);
871  return num_nodes;
872 }
#define CHKERR
Inline error check.
Definition: definitions.h:596
MoFEMErrorCode getNumberOfNodes(int &num_nodes) const
Get number of DOFs on element.

◆ getNumeredEntFiniteElementPtr()

boost::shared_ptr< const NumeredEntFiniteElement > MoFEM::ForcesAndSourcesCore::UserDataOperator::getNumeredEntFiniteElementPtr ( ) const

Return raw pointer to NumeredEntFiniteElement.

Examples
ElasticityMixedFormulation.hpp, HookeElement.cpp, and UnsaturatedFlow.hpp.

Definition at line 838 of file ForcesAndSourcesCore.hpp.

838  {
840 };
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr

◆ getOpType()

int MoFEM::ForcesAndSourcesCore::UserDataOperator::getOpType ( ) const

Get operator types.

Returns
Return operator type

Definition at line 878 of file ForcesAndSourcesCore.hpp.

◆ getPorblemRowIndices()

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::UserDataOperator::getPorblemRowIndices ( const std::string  filed_name,
const EntityType  type,
const int  side,
VectorInt indices 
) const

Definition at line 957 of file ForcesAndSourcesCore.hpp.

959  {
960  return getProblemRowIndices(filed_name, type, side, indices);
961 }
MoFEMErrorCode getProblemRowIndices(const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
Get row indices.

◆ getProblemColIndices()

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::UserDataOperator::getProblemColIndices ( const std::string  filed_name,
const EntityType  type,
const int  side,
VectorInt indices 
) const

Get col indices.

Field could be or not declared for this element but is declared for problem

Parameters
field_name
typeentity type
sideside number, any number if type is MBVERTEX
Returns
indices

NOTE: Using those indices to assemble matrix will result in error if new non-zero values need to be created.

Definition at line 1291 of file ForcesAndSourcesCore.cpp.

1293  {
1295  if (ptrFE == NULL)
1296  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
1297 
1298  switch (type) {
1299  case MBVERTEX:
1300  CHKERR ptrFE->getProblemNodesColIndices(field_name, indices);
1301  break;
1302  default:
1303  CHKERR ptrFE->getProblemTypeColIndices(field_name, type, side, indices);
1304  }
1306 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
MoFEMErrorCode getProblemNodesColIndices(const std::string &field_name, VectorInt &nodes_indices) const
MoFEMErrorCode getProblemTypeColIndices(const std::string &field_name, EntityType type, int side_number, VectorInt &indices) const
#define CHKERR
Inline error check.
Definition: definitions.h:596
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407

◆ getProblemRowIndices()

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::UserDataOperator::getProblemRowIndices ( const std::string  filed_name,
const EntityType  type,
const int  side,
VectorInt indices 
) const

Get row indices.

Field could be or not declared for this element but is declared for problem

Parameters
field_name
typeentity type
sideside number, any number if type is MBVERTEX
Returns
indices

NOTE: Using those indices to assemble matrix will result in error if new non-zero values need to be created.

Definition at line 1274 of file ForcesAndSourcesCore.cpp.

1276  {
1278  if (ptrFE == NULL)
1279  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
1280 
1281  switch (type) {
1282  case MBVERTEX:
1283  CHKERR ptrFE->getProblemNodesRowIndices(field_name, indices);
1284  break;
1285  default:
1286  CHKERR ptrFE->getProblemTypeRowIndices(field_name, type, side, indices);
1287  }
1289 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
MoFEMErrorCode getProblemNodesRowIndices(const std::string &field_name, VectorInt &nodes_indices) const
MoFEMErrorCode getProblemTypeRowIndices(const std::string &field_name, EntityType type, int side_number, VectorInt &indices) const
#define CHKERR
Inline error check.
Definition: definitions.h:596
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407

◆ getPtrFE()

ForcesAndSourcesCore * MoFEM::ForcesAndSourcesCore::UserDataOperator::getPtrFE ( ) const
protected

Definition at line 970 of file ForcesAndSourcesCore.hpp.

970  {
971  return ptrFE;
972 }

◆ getSideEntity()

EntityHandle MoFEM::ForcesAndSourcesCore::UserDataOperator::getSideEntity ( const int  side_number,
const EntityType  type 
)

Get the side entity.

Note
For vertex is expection. Side basses in argument of function doWork is zero. For other entity types side can be used as argument of this function.
MoFEMErrorCode doWork(int side, EntityType type,
if (type == MBVERTEX) {
for (int n = 0; n != number_of_nodes; ++n)
EntityHandle ent = getSideEntity(n, type);
// Do somthing
} else {
EntityHandle ent = getSideEntity(side, type);
// Do somthing
}
}
Parameters
side_number
type

Definition at line 860 of file ForcesAndSourcesCore.hpp.

861  {
862  if (auto side_ptr = getSideNumberPtr(side_number, type).lock())
863  return side_ptr->ent;
864  else
865  return 0;
866 }
boost::weak_ptr< SideNumber > getSideNumberPtr(const int side_number, const EntityType type)
Get the side number pointer.

◆ getSideNumberPtr()

boost::weak_ptr< SideNumber > MoFEM::ForcesAndSourcesCore::UserDataOperator::getSideNumberPtr ( const int  side_number,
const EntityType  type 
)

Get the side number pointer.

Note
For vertex is expection. Side basses in argument of function doWork is zero. For other entity types side can be used as argument of this function.
Parameters
side_number
type
Returns
boost::weak_ptr<SideNumber>

Definition at line 847 of file ForcesAndSourcesCore.hpp.

848  {
849  auto &side_table_by_side_and_type =
850  ptrFE->numeredEntFiniteElementPtr->getSideNumberTable().get<1>();
851  auto side_it =
852  side_table_by_side_and_type.find(boost::make_tuple(type, side_number));
853  if (side_it != side_table_by_side_and_type.end())
854  return *side_it;
855  else
856  return boost::weak_ptr<SideNumber>();
857 }
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr

◆ getSidePtrFE()

ForcesAndSourcesCore * MoFEM::ForcesAndSourcesCore::UserDataOperator::getSidePtrFE ( ) const
protected

Definition at line 975 of file ForcesAndSourcesCore.hpp.

975  {
976  return ptrFE->sidePtrFE;
977 }
ForcesAndSourcesCore * sidePtrFE
Element to integrate on the sides.

◆ getSnesA()

Mat MoFEM::ForcesAndSourcesCore::UserDataOperator::getSnesA ( ) const

Definition at line 908 of file ForcesAndSourcesCore.hpp.

908  {
909  return getFEMethod()->snes_A;
910 }
Mat snes_A
jacobian matrix
const FEMethod * getFEMethod() const
Return raw pointer to Finite Element Method object.

◆ getSnesB()

Mat MoFEM::ForcesAndSourcesCore::UserDataOperator::getSnesB ( ) const

Definition at line 912 of file ForcesAndSourcesCore.hpp.

912  {
913  return getFEMethod()->snes_B;
914 }
const FEMethod * getFEMethod() const
Return raw pointer to Finite Element Method object.
Mat snes_B
preconditioner of jacobian matrix

◆ getSnesF()

Vec MoFEM::ForcesAndSourcesCore::UserDataOperator::getSnesF ( ) const

Definition at line 900 of file ForcesAndSourcesCore.hpp.

900  {
901  return getFEMethod()->snes_f;
902 }
Vec snes_f
residual
const FEMethod * getFEMethod() const
Return raw pointer to Finite Element Method object.

◆ getSnesX()

Vec MoFEM::ForcesAndSourcesCore::UserDataOperator::getSnesX ( ) const

Definition at line 904 of file ForcesAndSourcesCore.hpp.

904  {
905  return getFEMethod()->snes_x;
906 }
const FEMethod * getFEMethod() const
Return raw pointer to Finite Element Method object.
Vec snes_x
state vector

◆ getTSA()

Mat MoFEM::ForcesAndSourcesCore::UserDataOperator::getTSA ( ) const

Definition at line 928 of file ForcesAndSourcesCore.hpp.

928  {
929  return getFEMethod()->ts_A;
930 }
const FEMethod * getFEMethod() const
Return raw pointer to Finite Element Method object.

◆ getTSa()

double MoFEM::ForcesAndSourcesCore::UserDataOperator::getTSa ( ) const

Definition at line 944 of file ForcesAndSourcesCore.hpp.

944  {
945  return getFEMethod()->ts_a;
946 }
const FEMethod * getFEMethod() const
Return raw pointer to Finite Element Method object.
PetscReal ts_a
shift for U_tt (see PETSc Time Solver)

◆ getTSB()

Mat MoFEM::ForcesAndSourcesCore::UserDataOperator::getTSB ( ) const

Definition at line 932 of file ForcesAndSourcesCore.hpp.

932  {
933  return getFEMethod()->ts_B;
934 }
const FEMethod * getFEMethod() const
Return raw pointer to Finite Element Method object.
Mat ts_B
Preconditioner for ts_A.

◆ getTSf()

Vec MoFEM::ForcesAndSourcesCore::UserDataOperator::getTSf ( ) const

Definition at line 924 of file ForcesAndSourcesCore.hpp.

924  {
925  return getFEMethod()->ts_F;
926 }
const FEMethod * getFEMethod() const
Return raw pointer to Finite Element Method object.
Vec ts_F
residual vector

◆ getTSstep()

int MoFEM::ForcesAndSourcesCore::UserDataOperator::getTSstep ( ) const

Definition at line 936 of file ForcesAndSourcesCore.hpp.

936  {
937  return getFEMethod()->ts_step;
938 }
const FEMethod * getFEMethod() const
Return raw pointer to Finite Element Method object.
PetscInt ts_step
time step

◆ getTStime()

double MoFEM::ForcesAndSourcesCore::UserDataOperator::getTStime ( ) const

Definition at line 940 of file ForcesAndSourcesCore.hpp.

940  {
941  return getFEMethod()->ts_t;
942 }
PetscReal ts_t
time
const FEMethod * getFEMethod() const
Return raw pointer to Finite Element Method object.

◆ getTSu()

Vec MoFEM::ForcesAndSourcesCore::UserDataOperator::getTSu ( ) const

Definition at line 916 of file ForcesAndSourcesCore.hpp.

916  {
917  return getFEMethod()->ts_u;
918 }
const FEMethod * getFEMethod() const
Return raw pointer to Finite Element Method object.
Vec ts_u
state vector

◆ getTSu_t()

Vec MoFEM::ForcesAndSourcesCore::UserDataOperator::getTSu_t ( ) const

Definition at line 920 of file ForcesAndSourcesCore.hpp.

920  {
921  return getFEMethod()->ts_u_t;
922 }
Vec ts_u_t
time derivative of state vector
const FEMethod * getFEMethod() const
Return raw pointer to Finite Element Method object.

◆ loopSide()

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::UserDataOperator::loopSide ( const string &  fe_name,
ForcesAndSourcesCore side_fe,
const size_t  dim 
)
private

User call this function to loop over elements on the side of face. This function calls finite element with is operator to do calculations.

Parameters
fe_name
side_fe
dim
Returns
MoFEMErrorCode

Definition at line 1316 of file ForcesAndSourcesCore.cpp.

1318  {
1320  const EntityHandle ent = getNumeredEntFiniteElementPtr()->getEnt();
1321  const Problem *problem_ptr = getFEMethod()->problemPtr;
1322  Range adjacent_ents;
1323  CHKERR ptrFE->mField.getInterface<BitRefManager>()->getAdjacenciesAny(
1324  ent, side_dim, adjacent_ents);
1325  typedef NumeredEntFiniteElement_multiIndex::index<
1326  Composite_Name_And_Ent_mi_tag>::type FEByComposite;
1327  FEByComposite &numered_fe =
1328  problem_ptr->numeredFiniteElements->get<Composite_Name_And_Ent_mi_tag>();
1329 
1330  side_fe->feName = fe_name;
1331 
1332  CHKERR side_fe->setSideFEPtr(ptrFE);
1333  CHKERR side_fe->copyBasicMethod(*getFEMethod());
1334  CHKERR side_fe->copyKsp(*getFEMethod());
1335  CHKERR side_fe->copySnes(*getFEMethod());
1336  CHKERR side_fe->copyTs(*getFEMethod());
1337 
1338  CHKERR side_fe->preProcess();
1339 
1340  int nn = 0;
1341  side_fe->loopSize = adjacent_ents.size();
1342  for (Range::iterator vit = adjacent_ents.begin(); vit != adjacent_ents.end();
1343  vit++) {
1344  FEByComposite::iterator miit =
1345  numered_fe.find(boost::make_tuple(fe_name, *vit));
1346  if (miit != numered_fe.end()) {
1347  side_fe->nInTheLoop = nn++;
1348  side_fe->numeredEntFiniteElementPtr = *miit;
1349  side_fe->dataFieldEntsPtr = (*miit)->sPtr->data_field_ents_view;
1350  side_fe->rowFieldEntsPtr = (*miit)->sPtr->row_field_ents_view;
1351  side_fe->colFieldEntsPtr = (*miit)->sPtr->col_field_ents_view;
1352  side_fe->dataPtr = (*miit)->sPtr->data_dofs;
1353  side_fe->rowPtr = (*miit)->rows_dofs;
1354  side_fe->colPtr = (*miit)->cols_dofs;
1355  CHKERR (*side_fe)();
1356  }
1357  }
1358 
1359  CHKERR side_fe->postProcess();
1360 
1362 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
const Problem * problemPtr
raw pointer to problem
const FEMethod * getFEMethod() const
Return raw pointer to Finite Element Method object.
boost::shared_ptr< const NumeredEntFiniteElement > getNumeredEntFiniteElementPtr() const
Return raw pointer to NumeredEntFiniteElement.
MoFEMErrorCode getInterface(const MOFEMuuid &uuid, IFACE *&iface) const
Get interface by uuid and return reference to pointer of interface.
#define CHKERR
Inline error check.
Definition: definitions.h:596
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407

◆ setOpType()

void MoFEM::ForcesAndSourcesCore::UserDataOperator::setOpType ( const OpType  type)

Set operator type.

Parameters
Operatortype

Definition at line 880 of file ForcesAndSourcesCore.hpp.

◆ setPtrFE()

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::UserDataOperator::setPtrFE ( ForcesAndSourcesCore ptr)
protected

Definition at line 964 of file ForcesAndSourcesCore.hpp.

964  {
966  ptrFE = ptr;
968 }
#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

Friends And Related Function Documentation

◆ EdgeElementForcesAndSourcesCoreBase

Definition at line 372 of file ForcesAndSourcesCore.hpp.

◆ FaceElementForcesAndSourcesCoreBase

Definition at line 373 of file ForcesAndSourcesCore.hpp.

Member Data Documentation

◆ colFieldName

std::string MoFEM::ForcesAndSourcesCore::UserDataOperator::colFieldName

Definition at line 107 of file ForcesAndSourcesCore.hpp.

◆ ForcesAndSourcesCore

friend MoFEM::ForcesAndSourcesCore::UserDataOperator::ForcesAndSourcesCore
private

Definition at line 371 of file ForcesAndSourcesCore.hpp.

◆ opType

char MoFEM::ForcesAndSourcesCore::UserDataOperator::opType

Definition at line 105 of file ForcesAndSourcesCore.hpp.

◆ ptrFE

ForcesAndSourcesCore* MoFEM::ForcesAndSourcesCore::UserDataOperator::ptrFE
protected

Definition at line 349 of file ForcesAndSourcesCore.hpp.

◆ rowFieldName

std::string MoFEM::ForcesAndSourcesCore::UserDataOperator::rowFieldName

Definition at line 106 of file ForcesAndSourcesCore.hpp.

◆ sPace

FieldSpace MoFEM::ForcesAndSourcesCore::UserDataOperator::sPace

Definition at line 108 of file ForcesAndSourcesCore.hpp.


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