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

#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, OPSPACE = 1 << 3,
  OPLAST = 1 << 3
}
 Controls loop over entities on element. More...
 
- Public Types inherited from MoFEM::DataOperator
using DoWorkLhsHookFunType = boost::function< MoFEMErrorCode(DataOperator *op_ptr, int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)>
 
using DoWorkRhsHookFunType = boost::function< MoFEMErrorCode(DataOperator *op_ptr, int side, EntityType type, EntitiesFieldData::EntData &data)>
 

Public Member Functions

 UserDataOperator (const FieldSpace space, const char type=OPSPACE, 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...
 
int getFEDim () const
 Get dimension of finite element. More...
 
EntityType getFEType () const
 Get dimension of finite element. More...
 
boost::weak_ptr< SideNumbergetSideNumberPtr (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 () const
 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...
 
std::string getFEName () const
 Get name of the element. More...
 
ForcesAndSourcesCoregetPtrFE () const
 
ForcesAndSourcesCoregetSidePtrFE () const
 
ForcesAndSourcesCoregetRefinePtrFE () const
 
Accessing KSP
const PetscData::SwitchesgetDataCtx () const
 
const KspMethod::KSPContext getKSPCtx () const
 
const SnesMethod::SNESContext getSNESCtx () const
 
const TSMethod::TSContext getTSCtx () const
 
Vec getKSPf () const
 
Mat getKSPA () const
 
Mat getKSPB () const
 
Accessing SNES
Vec getSNESf () const
 
Vec getSNESx () const
 
Mat getSNESA () const
 
Mat getSNESB () const
 
Accessing TS
Vec getTSu () const
 
Vec getTSu_t () const
 
Vec getTSu_tt () const
 
Vec getTSf () const
 
Mat getTSA () const
 
Mat getTSB () const
 
int getTSstep () const
 
double getTStime () const
 
double getTStimeStep () const
 
double getTSa () const
 
double getTSaa () const
 
Base funtions and integration points
MatrixDoublegetGaussPts ()
 matrix of integration (Gauss) points for Volume Element More...
 
auto getFTensor0IntegrationWeight ()
 Get integration weights. More...
 
Coordinates and access to internal data
MatrixDoublegetCoordsAtGaussPts ()
 Gauss points and weight, matrix (nb. of points x 3) More...
 
auto getFTensor1CoordsAtGaussPts ()
 Get coordinates at integration points assuming linear geometry. More...
 
Measures (area, volume, length, etc.)
double getMeasure () const
 get measure of element More...
 
doublegetMeasure ()
 get measure of element More...
 
- Public Member Functions inherited from MoFEM::DataOperator
 DataOperator (const bool symm=true)
 
virtual ~DataOperator ()=default
 
virtual MoFEMErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
 Operator for bi-linear form, usually to calculate values on left hand side. More...
 
virtual MoFEMErrorCode opLhs (EntitiesFieldData &row_data, EntitiesFieldData &col_data)
 
virtual MoFEMErrorCode doWork (int side, EntityType type, EntitiesFieldData::EntData &data)
 Operator for linear form, usually to calculate values on right hand side. More...
 
virtual MoFEMErrorCode opRhs (EntitiesFieldData &data, const bool error_if_no_base=false)
 
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
DoWorkLhsHookFunType doWorkLhsHook
 
DoWorkRhsHookFunType doWorkRhsHook
 
bool sYmm
 If true assume that matrix is symmetric structure. More...
 
std::array< bool, MBMAXTYPE > doEntities
 If true operator is executed for entity. More...
 
booldoVertices
 \deprectaed If false skip vertices More...
 
booldoEdges
 \deprectaed If false skip edges More...
 
booldoQuads
 \deprectaed More...
 
booldoTris
 \deprectaed More...
 
booldoTets
 \deprectaed More...
 
booldoPrisms
 \deprectaed More...
 

Static Public Attributes

static const char *const OpTypeNames []
 

Protected Member Functions

virtual MoFEMErrorCode setPtrFE (ForcesAndSourcesCore *ptr)
 

Protected Attributes

ForcesAndSourcesCoreptrFE
 

Friends

class ForcesAndSourcesCore
 
class EdgeElementForcesAndSourcesCore
 
class FaceElementForcesAndSourcesCore
 
class ContactPrismElementForcesAndSourcesCore
 

Loops

}

using AdjCache = std::map< EntityHandle, std::vector< boost::weak_ptr< NumeredEntFiniteElement > >>
 
MoFEMErrorCode loopSide (const string &fe_name, ForcesAndSourcesCore *side_fe, const size_t dim, const EntityHandle ent_for_side=0, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy, AdjCache *adj_cache=nullptr)
 User calls this function to loop over elements on the side of face. This function calls finite element with its operator to do calculations. More...
 
MoFEMErrorCode loopThis (const string &fe_name, ForcesAndSourcesCore *this_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User calls this function to loop over the same element using a different set of integration points. This function calls finite element with its operator to do calculations. More...
 
MoFEMErrorCode loopParent (const string &fe_name, ForcesAndSourcesCore *parent_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User calls this function to loop over parent elements. This function calls finite element with its operator to do calculations. More...
 
MoFEMErrorCode loopChildren (const string &fe_name, ForcesAndSourcesCore *child_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User calls this function to loop over parent elements. This function calls finite element with its operator to do calculations. More...
 

Detailed Description

Examples
analytical_poisson_field_split.cpp, continuity_check_on_contact_prism_side_ele.cpp, continuity_check_on_skeleton_with_simple_2d_for_h1.cpp, continuity_check_on_skeleton_with_simple_2d_for_hcurl.cpp, continuity_check_on_skeleton_with_simple_2d_for_hdiv.cpp, dynamic_first_order_con_law.cpp, ElasticityMixedFormulation.hpp, EshelbianPlasticity.cpp, forces_and_sources_testing_users_base.cpp, hcurl_curl_operator.cpp, hdiv_divergence_operator.cpp, hello_world.cpp, HookeElement.hpp, level_set.cpp, MagneticElement.hpp, mesh_smoothing.cpp, mortar_contact_thermal.cpp, NavierStokesElement.hpp, prism_elements_from_surface.cpp, prism_polynomial_approximation.cpp, quad_polynomial_approximation.cpp, schur_test_diag_mat.cpp, simple_contact_thermal.cpp, simple_elasticity.cpp, and simple_interface.cpp.

Definition at line 549 of file ForcesAndSourcesCore.hpp.

Member Typedef Documentation

◆ AdjCache

Definition at line 904 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.
  • OPSPACE no field is defined for such operator. Is usually used to modify base

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 

operator doWork function is executed on FE rows

OPCOL 

operator doWork function is executed on FE columns

OPROWCOL 

operator doWork is executed on FE rows &columns

OPSPACE 

operator do Work is execute on space data

OPLAST 
Deprecated:
would be removed

Definition at line 566 of file ForcesAndSourcesCore.hpp.

566  {
567  OPROW = 1 << 0, ///< operator doWork function is executed on FE rows
568  OPCOL = 1 << 1, ///< operator doWork function is executed on FE columns
569  OPROWCOL = 1 << 2, ///< operator doWork is executed on FE rows &columns
570  OPSPACE = 1 << 3, ///< operator do Work is execute on space data
571  OPLAST = 1 << 3 ///< @deprecated would be removed
572  };

Constructor & Destructor Documentation

◆ UserDataOperator() [1/3]

UserDataOperator::UserDataOperator ( const FieldSpace  space,
const char  type = OPSPACE,
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
EshelbianPlasticity.cpp.

Definition at line 1977 of file ForcesAndSourcesCore.cpp.

1980  : DataOperator(symm), opType(type), sPace(space), ptrFE(nullptr) {}

◆ UserDataOperator() [2/3]

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

Definition at line 1982 of file ForcesAndSourcesCore.cpp.

◆ UserDataOperator() [3/3]

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

Definition at line 1987 of file ForcesAndSourcesCore.cpp.

1990  : DataOperator(symm), opType(type), rowFieldName(row_field_name),
1991  colFieldName(col_field_name), sPace(LASTSPACE), ptrFE(nullptr) {}

Member Function Documentation

◆ addOpType()

void UserDataOperator::addOpType ( const OpType  type)
inline

Add operator type.

Definition at line 1052 of file ForcesAndSourcesCore.hpp.

1052  {
1053  opType |= type;
1054 }

◆ getCoordsAtGaussPts()

MatrixDouble & UserDataOperator::getCoordsAtGaussPts ( )
inline

Gauss points and weight, matrix (nb. of points x 3)

Column 0-2 integration points coordinate x, y and z, respectively. At rows are integration points.

Examples
MagneticElement.hpp, and UnsaturatedFlow.hpp.

Definition at line 1265 of file ForcesAndSourcesCore.hpp.

1265  {
1266  return static_cast<ForcesAndSourcesCore *>(ptrFE)->coordsAtGaussPts;
1267 }

◆ getDataCtx()

const PetscData::Switches & UserDataOperator::getDataCtx ( ) const
inline

Definition at line 1069 of file ForcesAndSourcesCore.hpp.

1069  {
1070  return getFEMethod()->data_ctx;
1071 }

◆ getFEDim()

int UserDataOperator::getFEDim ( ) const
inline

Get dimension of finite element.

Returns
int

Definition at line 1008 of file ForcesAndSourcesCore.hpp.

1008  {
1010 };

◆ getFEEntityHandle()

EntityHandle UserDataOperator::getFEEntityHandle ( ) const
inline

Return finite element entity handle.

Returns
Finite element entity handle
Examples
HookeInternalStressElement.hpp, and UnsaturatedFlow.hpp.

Definition at line 1004 of file ForcesAndSourcesCore.hpp.

1004  {
1005  return getNumeredEntFiniteElementPtr()->getEnt();
1006 }

◆ getFEMethod()

const FEMethod * UserDataOperator::getFEMethod ( ) const
inline

Return raw pointer to Finite Element Method object.

Examples
analytical_poisson_field_split.cpp, and UnsaturatedFlow.hpp.

Definition at line 1042 of file ForcesAndSourcesCore.hpp.

1042  {
1043  return ptrFE;
1044 }

◆ getFEName()

std::string UserDataOperator::getFEName ( ) const
inline

Get name of the element.

Definition at line 1064 of file ForcesAndSourcesCore.hpp.

1064  {
1065  return getFEMethod()->getFEName();
1066 }

◆ getFEType()

EntityType UserDataOperator::getFEType ( ) const
inline

Get dimension of finite element.

Returns
int

Definition at line 1012 of file ForcesAndSourcesCore.hpp.

1012  {
1014 };

◆ getFTensor0IntegrationWeight()

auto UserDataOperator::getFTensor0IntegrationWeight ( )
inline

Get integration weights.

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

Definition at line 1240 of file ForcesAndSourcesCore.hpp.

1240  {
1242  &(getGaussPts()(getGaussPts().size1() - 1, 0)));
1243 }

◆ getFTensor1CoordsAtGaussPts()

auto UserDataOperator::getFTensor1CoordsAtGaussPts ( )
inline

Get coordinates at integration points assuming linear geometry.

auto t_coords = getFTensor1CoordsAtGaussPts();
for(int gg = 0;gg!=nb_int_ptrs;gg++) {
// do something
++t_coords;
}
Examples
prism_polynomial_approximation.cpp, quad_polynomial_approximation.cpp, and UnsaturatedFlow.hpp.

Definition at line 1269 of file ForcesAndSourcesCore.hpp.

1269  {
1271  &getCoordsAtGaussPts()(0, 0), &getCoordsAtGaussPts()(0, 1),
1272  &getCoordsAtGaussPts()(0, 2));
1273 }

◆ getGaussPts()

MatrixDouble & UserDataOperator::getGaussPts ( )
inline

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
analytical_poisson_field_split.cpp, MagneticElement.hpp, and UnsaturatedFlow.hpp.

Definition at line 1236 of file ForcesAndSourcesCore.hpp.

1236  {
1237  return static_cast<ForcesAndSourcesCore *>(ptrFE)->gaussPts;
1238 }

◆ getKSPA()

Mat UserDataOperator::getKSPA ( ) const
inline

Definition at line 1096 of file ForcesAndSourcesCore.hpp.

1096  {
1097 #ifndef NDEBUG
1098  if (getFEMethod()->ksp_A == PETSC_NULL)
1099  THROW_MESSAGE("KSP not set A vector");
1100 #endif
1101  return getFEMethod()->ksp_A;
1102 }

◆ getKSPB()

Mat UserDataOperator::getKSPB ( ) const
inline
Examples
EshelbianPlasticity.cpp.

Definition at line 1104 of file ForcesAndSourcesCore.hpp.

1104  {
1105 #ifndef NDEBUG
1106  if (getFEMethod()->ksp_B == PETSC_NULL)
1107  THROW_MESSAGE("KSP not set B vector");
1108 #endif
1109  return getFEMethod()->ksp_B;
1110 }

◆ getKSPCtx()

const KspMethod::KSPContext UserDataOperator::getKSPCtx ( ) const
inline

Definition at line 1074 of file ForcesAndSourcesCore.hpp.

1074  {
1075  return getFEMethod()->ksp_ctx;
1076 }

◆ getKSPf()

Vec UserDataOperator::getKSPf ( ) const
inline

Definition at line 1088 of file ForcesAndSourcesCore.hpp.

1088  {
1089 #ifndef NDEBUG
1090  if (getFEMethod()->ksp_f == PETSC_NULL)
1091  THROW_MESSAGE("KSP not set F vector");
1092 #endif
1093  return getFEMethod()->ksp_f;
1094 }

◆ getLoopSize()

int UserDataOperator::getLoopSize ( ) const
inline

get size of elements in the loop

Returns
loop size

Definition at line 1060 of file ForcesAndSourcesCore.hpp.

1060  {
1061  return getFEMethod()->getLoopSize();
1062 }

◆ getMeasure() [1/2]

double& MoFEM::ForcesAndSourcesCore::UserDataOperator::getMeasure ( )
inline

get measure of element

Returns
volume

◆ getMeasure() [2/2]

double & UserDataOperator::getMeasure ( ) const
inline

get measure of element

Returns
volume
Examples
quad_polynomial_approximation.cpp.

Definition at line 1275 of file ForcesAndSourcesCore.hpp.

1275  {
1276  return static_cast<ForcesAndSourcesCore *>(ptrFE)->elementMeasure;
1277 }

◆ getNinTheLoop()

int UserDataOperator::getNinTheLoop ( ) const
inline

get number of finite element in the loop

Returns
number of finite element

Definition at line 1056 of file ForcesAndSourcesCore.hpp.

1056  {
1057  return getFEMethod()->getNinTheLoop();
1058 }

◆ getNumberOfNodesOnElement()

int UserDataOperator::getNumberOfNodesOnElement ( ) const
inline

Get the number of nodes on finite element.

Returns
int

Definition at line 1038 of file ForcesAndSourcesCore.hpp.

1038  {
1039  return ptrFE->getNumberOfNodes();
1040 }

◆ getNumeredEntFiniteElementPtr()

boost::shared_ptr< const NumeredEntFiniteElement > UserDataOperator::getNumeredEntFiniteElementPtr ( ) const
inline

Return raw pointer to NumeredEntFiniteElement.

Examples
ElasticityMixedFormulation.hpp, and UnsaturatedFlow.hpp.

Definition at line 1000 of file ForcesAndSourcesCore.hpp.

1000  {
1002 };

◆ getOpType()

int UserDataOperator::getOpType ( ) const
inline

Get operator types.

Returns
Return operator type

Definition at line 1046 of file ForcesAndSourcesCore.hpp.

1046 { return opType; }

◆ getProblemColIndices()

MoFEMErrorCode 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 1669 of file ForcesAndSourcesCore.cpp.

1671  {
1673  if (ptrFE == NULL)
1674  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
1675 
1676  switch (type) {
1677  case MBVERTEX:
1679  break;
1680  default:
1682  }
1684 }

◆ getProblemRowIndices()

MoFEMErrorCode 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 1652 of file ForcesAndSourcesCore.cpp.

1654  {
1656  if (ptrFE == NULL)
1657  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
1658 
1659  switch (type) {
1660  case MBVERTEX:
1662  break;
1663  default:
1665  }
1667 }

◆ getPtrFE()

ForcesAndSourcesCore * UserDataOperator::getPtrFE ( ) const
inline

Definition at line 1251 of file ForcesAndSourcesCore.hpp.

1251  {
1252  return ptrFE;
1253 }

◆ getRefinePtrFE()

ForcesAndSourcesCore * UserDataOperator::getRefinePtrFE ( ) const
inline

Definition at line 1261 of file ForcesAndSourcesCore.hpp.

1261  {
1262  return ptrFE->refinePtrFE;
1263 }

◆ getSideEntity()

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

Get the side entity.

Note
For vertex is expectation. 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)
// Do somthing
} else {
// Do somthing
}
}
Parameters
side_number
type

Definition at line 1030 of file ForcesAndSourcesCore.hpp.

1031  {
1032  if (auto side_ptr = getSideNumberPtr(side_number, type).lock())
1033  return side_ptr->ent;
1034  else
1035  return 0;
1036 }

◆ getSideNumberPtr()

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

Get the side number pointer.

Note
For vertex is expectation. 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 1017 of file ForcesAndSourcesCore.hpp.

1018  {
1019  auto &side_table_by_side_and_type =
1020  ptrFE->numeredEntFiniteElementPtr->getSideNumberTable().get<1>();
1021  auto side_it =
1022  side_table_by_side_and_type.find(boost::make_tuple(type, side_number));
1023  if (side_it != side_table_by_side_and_type.end())
1024  return *side_it;
1025  else
1026  return boost::weak_ptr<SideNumber>();
1027 }

◆ getSidePtrFE()

ForcesAndSourcesCore * UserDataOperator::getSidePtrFE ( ) const
inline

Definition at line 1256 of file ForcesAndSourcesCore.hpp.

1256  {
1257  return ptrFE->sidePtrFE;
1258 }

◆ getSNESA()

Mat UserDataOperator::getSNESA ( ) const
inline

Definition at line 1128 of file ForcesAndSourcesCore.hpp.

1128  {
1129 #ifndef NDEBUG
1130  if (getFEMethod()->snes_A == PETSC_NULL)
1131  THROW_MESSAGE("SNES not set A vector");
1132 #endif
1133  return getFEMethod()->snes_A;
1134 }

◆ getSNESB()

Mat UserDataOperator::getSNESB ( ) const
inline

Definition at line 1136 of file ForcesAndSourcesCore.hpp.

1136  {
1137 #ifndef NDEBUG
1138  if (getFEMethod()->snes_B == PETSC_NULL)
1139  THROW_MESSAGE("SNES not set A matrix");
1140 #endif
1141  return getFEMethod()->snes_B;
1142 }

◆ getSNESCtx()

const SnesMethod::SNESContext UserDataOperator::getSNESCtx ( ) const
inline

Definition at line 1079 of file ForcesAndSourcesCore.hpp.

1079  {
1080  return getFEMethod()->snes_ctx;
1081 }

◆ getSNESf()

Vec UserDataOperator::getSNESf ( ) const
inline

Definition at line 1112 of file ForcesAndSourcesCore.hpp.

1112  {
1113 #ifndef NDEBUG
1114  if (getFEMethod()->snes_f == PETSC_NULL)
1115  THROW_MESSAGE("SNES not set F vector");
1116 #endif
1117  return getFEMethod()->snes_f;
1118 }

◆ getSNESx()

Vec UserDataOperator::getSNESx ( ) const
inline

Definition at line 1120 of file ForcesAndSourcesCore.hpp.

1120  {
1121 #ifndef NDEBUG
1122  if (getFEMethod()->snes_x == PETSC_NULL)
1123  THROW_MESSAGE("SNESnot set X vector");
1124 #endif
1125  return getFEMethod()->snes_x;
1126 }

◆ getTSA()

Mat UserDataOperator::getTSA ( ) const
inline

Definition at line 1176 of file ForcesAndSourcesCore.hpp.

1176  {
1177 #ifndef NDEBUG
1178  if (getFEMethod()->ts_A == PETSC_NULL)
1179  THROW_MESSAGE("TS not set A matrix");
1180 #endif
1181  return getFEMethod()->ts_A;
1182 }

◆ getTSa()

double UserDataOperator::getTSa ( ) const
inline

Definition at line 1216 of file ForcesAndSourcesCore.hpp.

1216  {
1217 #ifndef NDEBUG
1219  .none() ||
1220  (getFEMethod()->data_ctx & (PetscData::CtxSetX_T)).none())
1221  THROW_MESSAGE("TS not set B matrix");
1222 #endif
1223  return getFEMethod()->ts_a;
1224 }

◆ getTSaa()

double UserDataOperator::getTSaa ( ) const
inline

Definition at line 1226 of file ForcesAndSourcesCore.hpp.

1226  {
1227 #ifndef NDEBUG
1229  .none() ||
1230  (getFEMethod()->data_ctx & (PetscData::CtxSetX_TT)).none())
1231  THROW_MESSAGE("TS not set B matrix");
1232 #endif
1233  return getFEMethod()->ts_aa;
1234 }

◆ getTSB()

Mat UserDataOperator::getTSB ( ) const
inline

Definition at line 1184 of file ForcesAndSourcesCore.hpp.

1184  {
1185 #ifndef NDEBUG
1186  if (getFEMethod()->ts_B == PETSC_NULL)
1187  THROW_MESSAGE("TS not set B matrix");
1188 #endif
1189  return getFEMethod()->ts_B;
1190 }

◆ getTSCtx()

const TSMethod::TSContext UserDataOperator::getTSCtx ( ) const
inline

Definition at line 1084 of file ForcesAndSourcesCore.hpp.

1084  {
1085  return getFEMethod()->ts_ctx;
1086 }

◆ getTSf()

Vec UserDataOperator::getTSf ( ) const
inline

Definition at line 1168 of file ForcesAndSourcesCore.hpp.

1168  {
1169 #ifndef NDEBUG
1170  if (getFEMethod()->ts_F == PETSC_NULL)
1171  THROW_MESSAGE("TS not set F vector");
1172 #endif
1173  return getFEMethod()->ts_F;
1174 }

◆ getTSstep()

int UserDataOperator::getTSstep ( ) const
inline

Definition at line 1192 of file ForcesAndSourcesCore.hpp.

1192  {
1193 #ifndef NDEBUG
1194  if ((getFEMethod()->data_ctx & PetscData::PetscData::CtxSetTime).none())
1195  THROW_MESSAGE("TS not set time");
1196 #endif
1197  return getFEMethod()->ts_step;
1198 }

◆ getTStime()

double UserDataOperator::getTStime ( ) const
inline

Definition at line 1200 of file ForcesAndSourcesCore.hpp.

1200  {
1201 #ifndef NDEBUG
1202  if ((getFEMethod()->data_ctx & PetscData::PetscData::CtxSetTime).none())
1203  THROW_MESSAGE("TS not set time");
1204 #endif
1205  return getFEMethod()->ts_t;
1206 }

◆ getTStimeStep()

double UserDataOperator::getTStimeStep ( ) const
inline

Definition at line 1208 of file ForcesAndSourcesCore.hpp.

1208  {
1209 #ifndef NDEBUG
1210  if ((getFEMethod()->data_ctx & PetscData::PetscData::CtxSetTime).none())
1211  THROW_MESSAGE("TS not set time");
1212 #endif
1213  return getFEMethod()->ts_dt;
1214 }

◆ getTSu()

Vec UserDataOperator::getTSu ( ) const
inline

Definition at line 1144 of file ForcesAndSourcesCore.hpp.

1144  {
1145 #ifndef NDEBUG
1146  if (getFEMethod()->ts_u == PETSC_NULL)
1147  THROW_MESSAGE("TS not set U vector");
1148 #endif
1149  return getFEMethod()->ts_u;
1150 }

◆ getTSu_t()

Vec UserDataOperator::getTSu_t ( ) const
inline

Definition at line 1152 of file ForcesAndSourcesCore.hpp.

1152  {
1153 #ifndef NDEBUG
1154  if (getFEMethod()->ts_u_t == PETSC_NULL)
1155  THROW_MESSAGE("TS not set U_t vector");
1156 #endif
1157  return getFEMethod()->ts_u_t;
1158 }

◆ getTSu_tt()

Vec UserDataOperator::getTSu_tt ( ) const
inline

Definition at line 1160 of file ForcesAndSourcesCore.hpp.

1160  {
1161 #ifndef NDEBUG
1162  if (getFEMethod()->ts_u_tt == PETSC_NULL)
1163  THROW_MESSAGE("TS not set U_tt vector");
1164 #endif
1165  return getFEMethod()->ts_u_tt;
1166 }

◆ loopChildren()

MoFEMErrorCode UserDataOperator::loopChildren ( const string &  fe_name,
ForcesAndSourcesCore child_fe,
const int  verb = QUIET,
const LogManager::SeverityLevel  sev = Sev::noisy 
)

User calls this function to loop over parent elements. This function calls finite element with its operator to do calculations.

Parameters
fe_name
child_fe
verb
sev
Returns
MoFEMErrorCode

Definition at line 1874 of file ForcesAndSourcesCore.cpp.

1876  {
1878 
1879  auto fe_miit = ptrFE->mField.get_finite_elements()
1880  ->get<FiniteElement_name_mi_tag>()
1881  .find(fe_name);
1882  if (fe_miit != ptrFE->mField.get_finite_elements()
1883  ->get<FiniteElement_name_mi_tag>()
1884  .end()) {
1885 
1886  const auto *problem_ptr = getFEMethod()->problemPtr;
1887  auto &ref_ents = *getPtrFE()->mField.get_ref_ents();
1888  auto &numered_fe =
1889  problem_ptr->numeredFiniteElementsPtr->get<Unique_mi_tag>();
1890 
1891  const auto parent_ent = getNumeredEntFiniteElementPtr()->getEnt();
1892  const auto parent_type = getNumeredEntFiniteElementPtr()->getEntType();
1893  auto range =
1894  ref_ents.get<Composite_ParentEnt_And_EntType_mi_tag>().equal_range(
1895  boost::make_tuple(parent_type, parent_ent));
1896 
1897  if (auto size = std::distance(range.first, range.second)) {
1898 
1899  std::vector<EntityHandle> childs_vec;
1900  childs_vec.reserve(size);
1901  for (; range.first != range.second; ++range.first)
1902  childs_vec.emplace_back((*range.first)->getEnt());
1903 
1904  Range childs;
1905 
1906  if ((childs_vec.back() - childs_vec.front() + 1) == size)
1907  childs = Range(childs_vec.front(), childs_vec.back());
1908  else
1909  childs.insert_list(childs_vec.begin(), childs_vec.end());
1910 
1911  child_fe->feName = fe_name;
1912 
1913  CHKERR child_fe->setRefineFEPtr(ptrFE);
1914  CHKERR child_fe->copyBasicMethod(*getFEMethod());
1915  CHKERR child_fe->copyPetscData(*getFEMethod());
1916  CHKERR child_fe->copyKsp(*getFEMethod());
1917  CHKERR child_fe->copySnes(*getFEMethod());
1918  CHKERR child_fe->copyTs(*getFEMethod());
1919 
1920  child_fe->cacheWeakPtr = getFEMethod()->cacheWeakPtr;
1921 
1922  CHKERR child_fe->preProcess();
1923 
1924  int nn = 0;
1925  child_fe->loopSize = size;
1926 
1927  for (auto p = childs.pair_begin(); p != childs.pair_end(); ++p) {
1928 
1929  auto miit =
1930  numered_fe.lower_bound(EntFiniteElement::getLocalUniqueIdCalculate(
1931  p->first, (*fe_miit)->getFEUId()));
1932  auto hi_miit =
1933  numered_fe.upper_bound(EntFiniteElement::getLocalUniqueIdCalculate(
1934  p->second, (*fe_miit)->getFEUId()));
1935 
1936  for (; miit != hi_miit; ++miit) {
1937 
1938  if (verb >= VERBOSE)
1939  MOFEM_LOG("SELF", sev)
1940  << "Child finite element " << "(" << nn << "): " << **miit;
1941 
1942  child_fe->nInTheLoop = nn++;
1943  child_fe->numeredEntFiniteElementPtr = *miit;
1944  CHKERR (*child_fe)();
1945  }
1946  }
1947  }
1948 
1949  CHKERR child_fe->postProcess();
1950  }
1951 
1953 }

◆ loopParent()

MoFEMErrorCode UserDataOperator::loopParent ( const string &  fe_name,
ForcesAndSourcesCore parent_fe,
const int  verb = QUIET,
const LogManager::SeverityLevel  sev = Sev::noisy 
)

User calls this function to loop over parent elements. This function calls finite element with its operator to do calculations.

Parameters
fe_name
parent_fe
verb
sev
Returns
MoFEMErrorCode

Definition at line 1825 of file ForcesAndSourcesCore.cpp.

1827  {
1829 
1830  auto &fes =
1831  ptrFE->mField.get_finite_elements()->get<FiniteElement_name_mi_tag>();
1832  auto fe_miit = fes.find(fe_name);
1833  if (fe_miit != fes.end()) {
1834 
1835  const auto *problem_ptr = getFEMethod()->problemPtr;
1836  auto &numered_fe =
1837  problem_ptr->numeredFiniteElementsPtr->get<Unique_mi_tag>();
1838 
1839  parent_fe->feName = fe_name;
1840  CHKERR parent_fe->setRefineFEPtr(ptrFE);
1841  CHKERR parent_fe->copyBasicMethod(*getFEMethod());
1842  CHKERR parent_fe->copyPetscData(*getFEMethod());
1843  CHKERR parent_fe->copyKsp(*getFEMethod());
1844  CHKERR parent_fe->copySnes(*getFEMethod());
1845  CHKERR parent_fe->copyTs(*getFEMethod());
1846 
1847  parent_fe->cacheWeakPtr = getFEMethod()->cacheWeakPtr;
1848 
1849  const auto parent_ent = getNumeredEntFiniteElementPtr()->getParentEnt();
1850  auto miit = numered_fe.find(EntFiniteElement::getLocalUniqueIdCalculate(
1851  parent_ent, (*fe_miit)->getFEUId()));
1852  if (miit != numered_fe.end()) {
1853  if (verb >= VERBOSE)
1854  MOFEM_LOG("SELF", sev) << "Parent finite element: " << **miit;
1855  parent_fe->loopSize = 1;
1856  parent_fe->nInTheLoop = 0;
1857  parent_fe->numeredEntFiniteElementPtr = *miit;
1858  CHKERR parent_fe->preProcess();
1859  CHKERR (*parent_fe)();
1860  CHKERR parent_fe->postProcess();
1861  } else {
1862  if (verb >= VERBOSE)
1863  MOFEM_LOG("SELF", sev) << "Parent finite element: no parent";
1864  parent_fe->loopSize = 0;
1865  parent_fe->nInTheLoop = 0;
1866  CHKERR parent_fe->preProcess();
1867  CHKERR parent_fe->postProcess();
1868  }
1869  }
1870 
1872 }

◆ loopSide()

MoFEMErrorCode UserDataOperator::loopSide ( const string &  fe_name,
ForcesAndSourcesCore side_fe,
const size_t  dim,
const EntityHandle  ent_for_side = 0,
const int  verb = QUIET,
const LogManager::SeverityLevel  sev = Sev::noisy,
AdjCache adj_cache = nullptr 
)

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

Parameters
fe_namename of the side element
side_fepointer to the side element instance
dimdimension the of side element
ent_for_sideentity handle for which adjacent volume or face will be accessed
verb
sev
adj_cache
Returns
MoFEMErrorCode

Definition at line 1700 of file ForcesAndSourcesCore.cpp.

1703  {
1705 
1706  const auto *problem_ptr = getFEMethod()->problemPtr;
1707  auto &numered_fe =
1708  problem_ptr->numeredFiniteElementsPtr->get<Unique_mi_tag>();
1709 
1710  auto fe_miit = ptrFE->mField.get_finite_elements()
1711  ->get<FiniteElement_name_mi_tag>()
1712  .find(fe_name);
1713  if (fe_miit != ptrFE->mField.get_finite_elements()
1714  ->get<FiniteElement_name_mi_tag>()
1715  .end()) {
1716 
1717  const auto ent = ent_for_side ? ent_for_side : getFEEntityHandle();
1718 
1719  side_fe->feName = fe_name;
1720 
1721  CHKERR side_fe->setSideFEPtr(ptrFE);
1722  CHKERR side_fe->copyBasicMethod(*getFEMethod());
1723  CHKERR side_fe->copyPetscData(*getFEMethod());
1724  CHKERR side_fe->copyKsp(*getFEMethod());
1725  CHKERR side_fe->copySnes(*getFEMethod());
1726  CHKERR side_fe->copyTs(*getFEMethod());
1727 
1728  side_fe->cacheWeakPtr = getFEMethod()->cacheWeakPtr;
1729 
1730  CHKERR side_fe->preProcess();
1731 
1732  std::vector<boost::weak_ptr<NumeredEntFiniteElement>> fe_vec;
1733  auto get_numered_fe_ptr = [&](auto &fe_uid, Range &&adjacent_ents)
1734  -> std::vector<boost::weak_ptr<NumeredEntFiniteElement>> & {
1735  fe_vec.reserve(adjacent_ents.size());
1736  for (auto fe_ent : adjacent_ents) {
1737  auto miit = numered_fe.find(
1739  if (miit != numered_fe.end()) {
1740  fe_vec.emplace_back(*miit);
1741  }
1742  }
1743  return fe_vec;
1744  };
1745 
1746  auto get_bit_entity_adjacency = [&]() {
1747  Range adjacent_ents;
1748  CHKERR ptrFE->mField.getInterface<BitRefManager>()->getAdjacenciesAny(
1749  ent, side_dim, adjacent_ents);
1750  return adjacent_ents;
1751  };
1752 
1753  auto get_adj = [&](auto &fe_uid)
1754  -> std::vector<boost::weak_ptr<NumeredEntFiniteElement>> & {
1755  if (adj_cache) {
1756  try {
1757  return (*adj_cache).at(ent);
1758  } catch (const std::out_of_range &) {
1759  return (*adj_cache)[ent] =
1760  get_numered_fe_ptr(fe_uid, get_bit_entity_adjacency());
1761  }
1762  } else {
1763  return get_numered_fe_ptr(fe_uid, get_bit_entity_adjacency());
1764  }
1765  };
1766 
1767  auto adj = get_adj((*fe_miit)->getFEUId());
1768 
1769  if (verb >= VERBOSE && !adj.empty())
1770  MOFEM_LOG("SELF", sev) << "Number of side finite elements " << adj.size();
1771 
1772  int nn = 0;
1773  side_fe->loopSize = adj.size();
1774  for (auto fe_weak_ptr : adj) {
1775  if (auto fe_ptr = fe_weak_ptr.lock()) {
1776  if (verb >= VERBOSE)
1777  MOFEM_LOG("SELF", sev)
1778  << "Side finite element " << "(" << nn << "): " << *fe_ptr;
1779  side_fe->nInTheLoop = nn;
1780  side_fe->numeredEntFiniteElementPtr = fe_ptr;
1781  CHKERR (*side_fe)();
1782  }
1783  ++nn;
1784  }
1785 
1786  CHKERR side_fe->postProcess();
1787  }
1788 
1790 }

◆ loopThis()

MoFEMErrorCode UserDataOperator::loopThis ( const string &  fe_name,
ForcesAndSourcesCore this_fe,
const int  verb = QUIET,
const LogManager::SeverityLevel  sev = Sev::noisy 
)

User calls this function to loop over the same element using a different set of integration points. This function calls finite element with its operator to do calculations.

Parameters
fe_name
this_fe
verb
sev
Returns
MoFEMErrorCode

Definition at line 1792 of file ForcesAndSourcesCore.cpp.

1794  {
1796 
1797  if (verb >= VERBOSE)
1798  MOFEM_LOG("SELF", sev) << "This finite element: "
1800 
1801  this_fe->feName = fe_name;
1802 
1803  CHKERR this_fe->setRefineFEPtr(ptrFE);
1804  CHKERR this_fe->copyBasicMethod(*getFEMethod());
1805  CHKERR this_fe->copyPetscData(*getFEMethod());
1806  CHKERR this_fe->copyKsp(*getFEMethod());
1807  CHKERR this_fe->copySnes(*getFEMethod());
1808  CHKERR this_fe->copyTs(*getFEMethod());
1809 
1810  this_fe->cacheWeakPtr = getFEMethod()->cacheWeakPtr;
1811 
1812  CHKERR this_fe->preProcess();
1813 
1814  this_fe->nInTheLoop = getNinTheLoop();
1815  this_fe->loopSize = getLoopSize();
1816  this_fe->numeredEntFiniteElementPtr = getNumeredEntFiniteElementPtr();
1817 
1818  CHKERR (*this_fe)();
1819 
1820  CHKERR this_fe->postProcess();
1821 
1823 }

◆ setOpType()

void UserDataOperator::setOpType ( const OpType  type)
inline

Set operator type.

Parameters
Operatortype

Definition at line 1048 of file ForcesAndSourcesCore.hpp.

1048  {
1049  opType = type;
1050 }

◆ setPtrFE()

MoFEMErrorCode UserDataOperator::setPtrFE ( ForcesAndSourcesCore ptr)
protectedvirtual

Friends And Related Function Documentation

◆ ContactPrismElementForcesAndSourcesCore

Definition at line 993 of file ForcesAndSourcesCore.hpp.

◆ EdgeElementForcesAndSourcesCore

friend class EdgeElementForcesAndSourcesCore
friend

Definition at line 991 of file ForcesAndSourcesCore.hpp.

◆ FaceElementForcesAndSourcesCore

friend class FaceElementForcesAndSourcesCore
friend

Definition at line 992 of file ForcesAndSourcesCore.hpp.

◆ ForcesAndSourcesCore

friend class ForcesAndSourcesCore
friend

Definition at line 990 of file ForcesAndSourcesCore.hpp.

Member Data Documentation

◆ colFieldName

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

Definition at line 578 of file ForcesAndSourcesCore.hpp.

◆ opType

char MoFEM::ForcesAndSourcesCore::UserDataOperator::opType

Definition at line 576 of file ForcesAndSourcesCore.hpp.

◆ OpTypeNames

const char *const UserDataOperator::OpTypeNames
static
Initial value:
= {
"OPROW", " OPCOL", "OPROWCOL", "OPSPACE", "OPLAST"}

Definition at line 574 of file ForcesAndSourcesCore.hpp.

◆ ptrFE

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

Definition at line 985 of file ForcesAndSourcesCore.hpp.

◆ rowFieldName

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

Definition at line 577 of file ForcesAndSourcesCore.hpp.

◆ sPace

FieldSpace MoFEM::ForcesAndSourcesCore::UserDataOperator::sPace

Definition at line 579 of file ForcesAndSourcesCore.hpp.


The documentation for this struct was generated from the following files:
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:460
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Definition: UnknownInterface.hpp:93
MoFEM::dimension_from_handle
auto dimension_from_handle(const EntityHandle h)
get entity dimension form handle
Definition: Templates.hpp:1914
MoFEM::ForcesAndSourcesCore::UserDataOperator::getSideEntity
EntityHandle getSideEntity(const int side_number, const EntityType type)
Get the side entity.
Definition: ForcesAndSourcesCore.hpp:1030
MoFEM::BasicMethod::getNinTheLoop
int getNinTheLoop() const
get number of evaluated element in the loop
Definition: LoopMethods.hpp:207
MoFEM::ForcesAndSourcesCore::UserDataOperator::getSideNumberPtr
boost::weak_ptr< SideNumber > getSideNumberPtr(const int side_number, const EntityType type)
Get the side number pointer.
Definition: ForcesAndSourcesCore.hpp:1017
MoFEM::TSMethod::ts_a
PetscReal ts_a
shift for U_t (see PETSc Time Solver)
Definition: LoopMethods.hpp:160
MoFEM::KspMethod::ksp_B
Mat & ksp_B
Definition: LoopMethods.hpp:91
FTensor::Tensor1
Definition: Tensor1_value.hpp:8
EntityHandle
MoFEM::ForcesAndSourcesCore::UserDataOperator::OPSPACE
@ OPSPACE
operator do Work is execute on space data
Definition: ForcesAndSourcesCore.hpp:570
MoFEM::PetscData::CtxSetB
static constexpr Switches CtxSetB
Definition: LoopMethods.hpp:38
MoFEM::TSMethod::ts_F
Vec & ts_F
residual vector
Definition: LoopMethods.hpp:168
MoFEM::PetscData::CtxSetA
static constexpr Switches CtxSetA
Definition: LoopMethods.hpp:37
MoFEM::SnesMethod::snes_x
Vec & snes_x
state vector
Definition: LoopMethods.hpp:121
MoFEM::SnesMethod::snes_ctx
SNESContext snes_ctx
Definition: LoopMethods.hpp:118
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
MoFEM::SnesMethod::snes_B
Mat & snes_B
preconditioner of jacobian matrix
Definition: LoopMethods.hpp:124
MoFEM::ForcesAndSourcesCore::sidePtrFE
ForcesAndSourcesCore * sidePtrFE
Element to integrate on the sides.
Definition: ForcesAndSourcesCore.hpp:507
MoFEM::CoreInterface::get_finite_elements
virtual const FiniteElement_multiIndex * get_finite_elements() const =0
Get the finite elements object.
MoFEM::ForcesAndSourcesCore::UserDataOperator::rowFieldName
std::string rowFieldName
Definition: ForcesAndSourcesCore.hpp:577
MoFEM::SnesMethod::snes_A
Mat & snes_A
jacobian matrix
Definition: LoopMethods.hpp:123
MoFEM::ForcesAndSourcesCore::UserDataOperator::OPROWCOL
@ OPROWCOL
operator doWork is executed on FE rows &columns
Definition: ForcesAndSourcesCore.hpp:569
MoFEM::DataOperator::doWork
virtual MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
Operator for bi-linear form, usually to calculate values on left hand side.
Definition: DataOperators.hpp:40
MoFEM::ForcesAndSourcesCore::getProblemTypeRowIndices
MoFEMErrorCode getProblemTypeRowIndices(const std::string &field_name, EntityType type, int side_number, VectorInt &indices) const
Definition: ForcesAndSourcesCore.cpp:572
THROW_MESSAGE
#define THROW_MESSAGE(msg)
Throw MoFEM exception.
Definition: definitions.h:574
MoFEM::ForcesAndSourcesCore::getProblemTypeColIndices
MoFEMErrorCode getProblemTypeColIndices(const std::string &field_name, EntityType type, int side_number, VectorInt &indices) const
Definition: ForcesAndSourcesCore.cpp:586
MoFEM::TSMethod::ts_A
Mat & ts_A
Definition: LoopMethods.hpp:170
MoFEM::DataOperator::DataOperator
DataOperator(const bool symm=true)
Definition: DataOperators.cpp:21
MoFEM::ForcesAndSourcesCore::UserDataOperator::getFTensor0IntegrationWeight
auto getFTensor0IntegrationWeight()
Get integration weights.
Definition: ForcesAndSourcesCore.hpp:1240
MoFEM::ForcesAndSourcesCore::UserDataOperator::getGaussPts
MatrixDouble & getGaussPts()
matrix of integration (Gauss) points for Volume Element
Definition: ForcesAndSourcesCore.hpp:1236
MoFEM::ForcesAndSourcesCore::UserDataOperator::ForcesAndSourcesCore
friend class ForcesAndSourcesCore
Definition: ForcesAndSourcesCore.hpp:990
VERBOSE
@ VERBOSE
Definition: definitions.h:222
MoFEM::ForcesAndSourcesCore::elementMeasure
double elementMeasure
Definition: ForcesAndSourcesCore.hpp:545
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
MoFEM::ForcesAndSourcesCore::UserDataOperator::colFieldName
std::string colFieldName
Definition: ForcesAndSourcesCore.hpp:578
MoFEM::ForcesAndSourcesCore::refinePtrFE
ForcesAndSourcesCore * refinePtrFE
Element to integrate parent or child.
Definition: ForcesAndSourcesCore.hpp:524
MoFEM::BasicMethod::getLoopSize
int getLoopSize() const
get loop size
Definition: LoopMethods.hpp:211
MoFEM::ForcesAndSourcesCore::coordsAtGaussPts
MatrixDouble coordsAtGaussPts
coordinated at gauss points
Definition: ForcesAndSourcesCore.hpp:544
MoFEM::PetscData::data_ctx
Switches data_ctx
Definition: LoopMethods.hpp:44
MoFEM::FEMethod::getNumberOfNodes
auto getNumberOfNodes() const
Definition: LoopMethods.hpp:456
convert.type
type
Definition: convert.py:64
MoFEM::ForcesAndSourcesCore::UserDataOperator::OPCOL
@ OPCOL
operator doWork function is executed on FE columns
Definition: ForcesAndSourcesCore.hpp:568
MoFEM::ForcesAndSourcesCore::UserDataOperator::opType
char opType
Definition: ForcesAndSourcesCore.hpp:576
MoFEM::EntFiniteElement::getLocalUniqueIdCalculate
UId getLocalUniqueIdCalculate() const
Generate UId for finite element entity.
Definition: FEMultiIndices.hpp:528
MoFEM::TSMethod::ts_u_t
Vec & ts_u_t
time derivative of state vector
Definition: LoopMethods.hpp:166
MoFEM::ForcesAndSourcesCore::UserDataOperator::getNinTheLoop
int getNinTheLoop() const
get number of finite element in the loop
Definition: ForcesAndSourcesCore.hpp:1056
MoFEM::ForcesAndSourcesCore::UserDataOperator::getLoopSize
int getLoopSize() const
get size of elements in the loop
Definition: ForcesAndSourcesCore.hpp:1060
MoFEM::PetscData::CtxSetX_T
static constexpr Switches CtxSetX_T
Definition: LoopMethods.hpp:40
MoFEM::ForcesAndSourcesCore::UserDataOperator::getPtrFE
ForcesAndSourcesCore * getPtrFE() const
Definition: ForcesAndSourcesCore.hpp:1251
MoFEM::ForcesAndSourcesCore::UserDataOperator::getNumeredEntFiniteElementPtr
boost::shared_ptr< const NumeredEntFiniteElement > getNumeredEntFiniteElementPtr() const
Return raw pointer to NumeredEntFiniteElement.
Definition: ForcesAndSourcesCore.hpp:1000
LASTSPACE
@ LASTSPACE
FieldSpace in [ 0, LASTSPACE )
Definition: definitions.h:89
MoFEM::TSMethod::ts_ctx
TSContext ts_ctx
Definition: LoopMethods.hpp:157
MoFEM::FEMethod::getFEName
auto getFEName() const
get finite element name
Definition: LoopMethods.hpp:391
MoFEM::ForcesAndSourcesCore::UserDataOperator::getCoordsAtGaussPts
MatrixDouble & getCoordsAtGaussPts()
Gauss points and weight, matrix (nb. of points x 3)
Definition: ForcesAndSourcesCore.hpp:1265
MoFEM::type_from_handle
auto type_from_handle(const EntityHandle h)
get type from entity handle
Definition: Templates.hpp:1898
MoFEM::Problem::numeredFiniteElementsPtr
boost::shared_ptr< NumeredEntFiniteElement_multiIndex > numeredFiniteElementsPtr
store finite elements
Definition: ProblemsMultiIndices.hpp:77
EntData
EntitiesFieldData::EntData EntData
Definition: child_and_parent.cpp:37
field_name
constexpr auto field_name
Definition: poisson_2d_homogeneous.cpp:13
MoFEM::ForcesAndSourcesCore::UserDataOperator::getFEMethod
const FEMethod * getFEMethod() const
Return raw pointer to Finite Element Method object.
Definition: ForcesAndSourcesCore.hpp:1042
convert.n
n
Definition: convert.py:82
MoFEM::PetscData::CtxSetX_TT
static constexpr Switches CtxSetX_TT
Definition: LoopMethods.hpp:41
MoFEM::TSMethod::ts_step
PetscInt ts_step
time step number
Definition: LoopMethods.hpp:159
MoFEM::KspMethod::ksp_A
Mat & ksp_A
Definition: LoopMethods.hpp:90
Range
MoFEM::SnesMethod::snes_f
Vec & snes_f
residual
Definition: LoopMethods.hpp:122
MOFEM_LOG
#define MOFEM_LOG(channel, severity)
Log.
Definition: LogManager.hpp:308
MoFEM::TSMethod::ts_t
PetscReal ts_t
time
Definition: LoopMethods.hpp:162
FTensor::Tensor0
Definition: Tensor0.hpp:16
MoFEM::ForcesAndSourcesCore::UserDataOperator::getFEEntityHandle
EntityHandle getFEEntityHandle() const
Return finite element entity handle.
Definition: ForcesAndSourcesCore.hpp:1004
MoFEM::TSMethod::ts_aa
PetscReal ts_aa
shift for U_tt shift for U_tt
Definition: LoopMethods.hpp:161
MoFEM::TSMethod::ts_u
Vec & ts_u
state vector
Definition: LoopMethods.hpp:165
MoFEM::ForcesAndSourcesCore::getProblemNodesColIndices
MoFEMErrorCode getProblemNodesColIndices(const std::string &field_name, VectorInt &nodes_indices) const
Definition: ForcesAndSourcesCore.cpp:579
MoFEM::CoreInterface::get_ref_ents
virtual const RefEntity_multiIndex * get_ref_ents() const =0
Get the ref ents object.
MoFEM::BasicMethod::cacheWeakPtr
boost::weak_ptr< CacheTuple > cacheWeakPtr
Definition: LoopMethods.hpp:351
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
MoFEM::BasicMethod::problemPtr
const Problem * problemPtr
raw pointer to problem
Definition: LoopMethods.hpp:250
MoFEM::ForcesAndSourcesCore::gaussPts
MatrixDouble gaussPts
Matrix of integration points.
Definition: ForcesAndSourcesCore.hpp:109
MoFEM::ForcesAndSourcesCore::UserDataOperator::sPace
FieldSpace sPace
Definition: ForcesAndSourcesCore.hpp:579
MoFEM::KspMethod::ksp_f
Vec & ksp_f
Definition: LoopMethods.hpp:89
MoFEM::FEMethod::numeredEntFiniteElementPtr
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
Definition: LoopMethods.hpp:383
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:453
MoFEM::TSMethod::ts_dt
PetscReal ts_dt
time step size
Definition: LoopMethods.hpp:163
MoFEM::TSMethod::ts_B
Mat & ts_B
Preconditioner for ts_A.
Definition: LoopMethods.hpp:172
MoFEM::ForcesAndSourcesCore::UserDataOperator::OPLAST
@ OPLAST
Definition: ForcesAndSourcesCore.hpp:571
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
MoFEM::ForcesAndSourcesCore::getProblemNodesRowIndices
MoFEMErrorCode getProblemNodesRowIndices(const std::string &field_name, VectorInt &nodes_indices) const
Definition: ForcesAndSourcesCore.cpp:565
MoFEM::ForcesAndSourcesCore::UserDataOperator::ptrFE
ForcesAndSourcesCore * ptrFE
Definition: ForcesAndSourcesCore.hpp:985
MoFEM::TSMethod::ts_u_tt
Vec & ts_u_tt
second time derivative of state vector
Definition: LoopMethods.hpp:167
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
MoFEM::ForcesAndSourcesCore::UserDataOperator::getFTensor1CoordsAtGaussPts
auto getFTensor1CoordsAtGaussPts()
Get coordinates at integration points assuming linear geometry.
Definition: ForcesAndSourcesCore.hpp:1269
MoFEM::ForcesAndSourcesCore::mField
Interface & mField
Definition: ForcesAndSourcesCore.hpp:24
MoFEM::ForcesAndSourcesCore::UserDataOperator::OPROW
@ OPROW
operator doWork function is executed on FE rows
Definition: ForcesAndSourcesCore.hpp:567
MoFEM::KspMethod::ksp_ctx
KSPContext ksp_ctx
Context.
Definition: LoopMethods.hpp:86