14static char help[] =
"...\n\n";
16#define HelloFunctionBegin \
18 MOFEM_LOG_CHANNEL("SYNC"); \
19 MOFEM_LOG_FUNCTION(); \
20 MOFEM_LOG_TAG("WORLD", "HelloWorld");
22struct OpRow :
public ForcesAndSourcesCore::UserDataOperator {
28 if (type == MBVERTEX) {
30 MOFEM_LOG(
"SYNC", Sev::inform) <<
"**** " << getNinTheLoop() <<
" ****";
31 MOFEM_LOG(
"SYNC", Sev::inform) <<
"**** Operators ****";
34 <<
"Hello Operator OpRow:"
35 <<
" field name " << rowFieldName <<
" side " << side <<
" type "
36 << CN::EntityTypeName(type) <<
" nb dofs on entity "
42struct OpRowCol :
public ForcesAndSourcesCore::UserDataOperator {
43 OpRowCol(
const std::string row_field,
const std::string col_field,
53 <<
"Hello Operator OpRowCol:"
54 <<
" row field name " << rowFieldName <<
" row side " << row_side
55 <<
" row type " << CN::EntityTypeName(row_type)
56 <<
" nb dofs on row entity " << row_data.
getIndices().size() <<
" : "
57 <<
" col field name " << colFieldName <<
" col side " << col_side
58 <<
" col type " << CN::EntityTypeName(col_type)
59 <<
" nb dofs on col entity " << col_data.
getIndices().size();
64struct OpVolume :
public VolumeElementForcesAndSourcesCore::UserDataOperator {
72 if (type == MBVERTEX) {
73 MOFEM_LOG(
"SYNC", Sev::inform) <<
"Hello Operator OpVolume:"
74 <<
" volume " << getVolume();
80struct OpFace :
public FaceElementForcesAndSourcesCore::UserDataOperator {
86 if (type == MBVERTEX) {
87 MOFEM_LOG(
"SYNC", Sev::inform) <<
"Hello Operator OpFace:"
94struct OpFaceSide :
public FaceElementForcesAndSourcesCore::UserDataOperator {
95 boost::shared_ptr<VolumeElementForcesAndSourcesCoreOnSide> &
feSidePtr;
98 boost::shared_ptr<VolumeElementForcesAndSourcesCoreOnSide> &fe_side_ptr)
105 if (type == MBVERTEX) {
106 MOFEM_LOG(
"SYNC", Sev::inform) <<
"Hello Operator OpSideFace";
114 :
public VolumeElementForcesAndSourcesCoreOnSide::UserDataOperator {
121 if (type == MBVERTEX) {
123 <<
"Hello Operator OpVolumeSide:"
124 <<
" volume " << getVolume() <<
" normal " << getNormal();
130int main(
int argc,
char *argv[]) {
138 DMType dm_name =
"DMMOFEM";
142 moab::Core moab_core;
143 moab::Interface &moab = moab_core;
173 pipeline_mng->getOpDomainRhsPipeline().push_back(
new OpVolume(
"U"));
174 pipeline_mng->getOpDomainLhsPipeline().push_back(
179 pipeline_mng->getOpBoundaryRhsPipeline().push_back(
new OpRow(
"L"));
180 pipeline_mng->getOpBoundaryRhsPipeline().push_back(
new OpFace(
"L"));
181 pipeline_mng->getOpBoundaryLhsPipeline().push_back(
186 boost::shared_ptr<VolumeElementForcesAndSourcesCoreOnSide> side_fe(
188 side_fe->getOpPtrVector().push_back(
new OpVolumeSide(
"U"));
192 pipeline_mng->getOpSkeletonRhsPipeline().push_back(
new OpRow(
"S"));
193 pipeline_mng->getOpSkeletonRhsPipeline().push_back(
198 CHKERR pipeline_mng->loopFiniteElements();
ForcesAndSourcesCore::UserDataOperator UserDataOperator
#define MOFEM_LOG_SYNCHRONISE(comm)
Synchronise "SYNC" channel.
void simple(double P1[], double P2[], double P3[], double c[], const int N)
#define CATCH_ERRORS
Catch errors.
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base nme:nme847.
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
boost::ptr_deque< UserDataOperator > & getOpDomainRhsPipeline()
Get the Op Domain Rhs Pipeline object.
#define MOFEM_LOG(channel, severity)
Log.
#define HelloFunctionBegin
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
implementation of Data Operators for Forces and Sources
constexpr auto field_name
virtual MPI_Comm & get_comm() const =0
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Deprecated interface functions.
Data on single entity (This is passed as argument to DataOperator::doWork)
const VectorInt & getIndices() const
Get global indices of dofs on entity.
VectorDouble & getNormal()
get triangle normal
@ OPROW
operator doWork function is executed on FE rows
structure to get information form mofem into EntitiesFieldData
PipelineManager interface.
Simple interface for fast problem set-up.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Volume finite element base.
Base volume element used to integrate on skeleton.
OpFace(const std::string &field_name)
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
boost::shared_ptr< VolumeElementForcesAndSourcesCoreOnSide > & feSidePtr
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
OpFaceSide(const std::string &field_name, boost::shared_ptr< VolumeElementForcesAndSourcesCoreOnSide > &fe_side_ptr)
MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
OpRowCol(const std::string row_field, const std::string col_field, const bool symm)
OpRow(const std::string &field_name)
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
OpVolume(const std::string &field_name)
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
OpVolumeSide(const std::string &field_name)
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)