v0.14.0
hello_world.cpp
Go to the documentation of this file.
1 /**
2  * \file hello_world.cpp
3  * \ingroup mofem_simple_interface
4  * \example hello_world.cpp
5  *
6  * Prints basic information about users data operator.
7  * See more details in \ref hello_world_tut1
8  *
9  */
10 
11 #include <BasicFiniteElements.hpp>
12 using namespace MoFEM;
13 
14 static char help[] = "...\n\n";
15 
16 #define HelloFunctionBegin \
17  MoFEMFunctionBegin; \
18  MOFEM_LOG_CHANNEL("SYNC"); \
19  MOFEM_LOG_FUNCTION(); \
20  MOFEM_LOG_TAG("WORLD", "HelloWorld");
21 
23  OpRow(const std::string &field_name)
25  MoFEMErrorCode doWork(int side, EntityType type,
28  if (type == MBVERTEX) {
29  // get number of evaluated element in the loop
30  MOFEM_LOG("SYNC", Sev::inform) << "**** " << getNinTheLoop() << " ****";
31  MOFEM_LOG("SYNC", Sev::inform) << "**** Operators ****";
32  }
33  MOFEM_LOG("SYNC", Sev::inform)
34  << "Hello Operator OpRow:"
35  << " field name " << rowFieldName << " side " << side << " type "
36  << CN::EntityTypeName(type) << " nb dofs on entity "
37  << data.getIndices().size();
39  }
40 };
41 
43  OpRowCol(const std::string row_field, const std::string col_field,
44  const bool symm)
45  : ForcesAndSourcesCore::UserDataOperator(row_field, col_field, OPROWCOL,
46  symm) {}
47  MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type,
48  EntityType col_type,
50  EntitiesFieldData::EntData &col_data) {
52  MOFEM_LOG("SYNC", Sev::inform)
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();
61  }
62 };
63 
65  OpVolume(const std::string &field_name)
67  field_name, OPROW) {
68  }
69  MoFEMErrorCode doWork(int side, EntityType type,
72  if (type == MBVERTEX) {
73  MOFEM_LOG("SYNC", Sev::inform) << "Hello Operator OpVolume:"
74  << " volume " << getVolume();
75  }
77  }
78 };
79 
81  OpFace(const std::string &field_name)
83  MoFEMErrorCode doWork(int side, EntityType type,
86  if (type == MBVERTEX) {
87  MOFEM_LOG("SYNC", Sev::inform) << "Hello Operator OpFace:"
88  << " normal " << getNormal();
89  }
91  }
92 };
93 
95  boost::shared_ptr<VolumeElementForcesAndSourcesCoreOnSide> &feSidePtr;
97  const std::string &field_name,
98  boost::shared_ptr<VolumeElementForcesAndSourcesCoreOnSide> &fe_side_ptr)
100  feSidePtr(fe_side_ptr) {}
101  MoFEMErrorCode doWork(int side, EntityType type,
103 
105  if (type == MBVERTEX) {
106  MOFEM_LOG("SYNC", Sev::inform) << "Hello Operator OpSideFace";
107  CHKERR loopSideVolumes("dFE", *feSidePtr);
108  }
110  }
111 };
112 
115  OpVolumeSide(const std::string &field_name)
117  field_name, field_name, OPROW) {}
118  MoFEMErrorCode doWork(int side, EntityType type,
121  if (type == MBVERTEX) {
122  MOFEM_LOG("SYNC", Sev::inform)
123  << "Hello Operator OpVolumeSide:"
124  << " volume " << getVolume() << " normal " << getNormal();
125  }
127  }
128 };
129 
130 int main(int argc, char *argv[]) {
131 
132  // initialize petsc
133  MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
134 
135  try {
136 
137  // Register DM Manager
138  DMType dm_name = "DMMOFEM";
139  CHKERR DMRegister_MoFEM(dm_name);
140 
141  // Create MoAB database
142  moab::Core moab_core;
143  moab::Interface &moab = moab_core;
144 
145  // Create MoFEM database and link it to MoAB
146  MoFEM::Core mofem_core(moab);
147  MoFEM::Interface &m_field = mofem_core;
148 
149  // Simple interface
150  Simple *simple;
151  CHKERR m_field.getInterface(simple);
152 
153  // get options from command line
154  CHKERR simple->getOptions();
155  // load mesh file
156  CHKERR simple->loadFile();
157  // add fields
158  CHKERR simple->addDomainField("U", H1, AINSWORTH_LEGENDRE_BASE, 3);
159  CHKERR simple->addBoundaryField("L", H1, AINSWORTH_LEGENDRE_BASE, 3);
160  CHKERR simple->addSkeletonField("S", H1, AINSWORTH_LEGENDRE_BASE, 3);
161  // set fields order
162  CHKERR simple->setFieldOrder("U", 4);
163  CHKERR simple->setFieldOrder("L", 3);
164  CHKERR simple->setFieldOrder("S", 3);
165 
166  // setup problem
167  CHKERR simple->setUp();
168 
169  auto pipeline_mng = m_field.getInterface<PipelineManager>();
170 
171  //! [set operator to the volume element instance]
172  pipeline_mng->getOpDomainRhsPipeline().push_back(new OpRow("U"));
173  pipeline_mng->getOpDomainRhsPipeline().push_back(new OpVolume("U"));
174  pipeline_mng->getOpDomainLhsPipeline().push_back(
175  new OpRowCol("U", "U", true));
176  //! [set operator to the volume element instance]
177 
178  //! [set operator to the face element instance]
179  pipeline_mng->getOpBoundaryRhsPipeline().push_back(new OpRow("L"));
180  pipeline_mng->getOpBoundaryRhsPipeline().push_back(new OpFace("L"));
181  pipeline_mng->getOpBoundaryLhsPipeline().push_back(
182  new OpRowCol("U", "L", false));
183  //! [set operator to the face element instance]
184 
185  //! [create skeleton side element and push operator to it]
186  boost::shared_ptr<VolumeElementForcesAndSourcesCoreOnSide> side_fe(
188  side_fe->getOpPtrVector().push_back(new OpVolumeSide("U"));
189  //! [create skeleton side element and push operator to it]
190 
191  //! [set operator to the face element on skeleton instance]
192  pipeline_mng->getOpSkeletonRhsPipeline().push_back(new OpRow("S"));
193  pipeline_mng->getOpSkeletonRhsPipeline().push_back(
194  new OpFaceSide("S", side_fe));
195  //! [set operator to the face element on skeleton instance]
196 
197  //! [executing finite elements]
198  CHKERR pipeline_mng->loopFiniteElements();
199  //! [executing finite elements]
200 
201  MOFEM_LOG_SYNCHRONISE(m_field.get_comm());
202  }
203  CATCH_ERRORS;
204 
205  // finish work cleaning memory, getting statistics, etc.
207 
208  return 0;
209 }
MoFEM::VolumeElementForcesAndSourcesCoreOnSide::UserDataOperator
default operator for TET element
Definition: VolumeElementForcesAndSourcesCoreOnSide.hpp:97
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.
Definition: UnknownInterface.hpp:93
MoFEM::EntitiesFieldData::EntData
Data on single entity (This is passed as argument to DataOperator::doWork)
Definition: EntitiesFieldData.hpp:127
OpFaceSide::OpFaceSide
OpFaceSide(const std::string &field_name, boost::shared_ptr< VolumeElementForcesAndSourcesCoreOnSide > &fe_side_ptr)
Definition: hello_world.cpp:96
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
H1
@ H1
continuous field
Definition: definitions.h:85
OpFaceSide::doWork
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Definition: hello_world.cpp:101
OpRowCol
Definition: hello_world.cpp:42
HelloFunctionBegin
#define HelloFunctionBegin
Definition: hello_world.cpp:16
MoFEM::CoreInterface::get_comm
virtual MPI_Comm & get_comm() const =0
help
static char help[]
Definition: hello_world.cpp:14
OpFace::doWork
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Definition: hello_world.cpp:83
OpVolume::OpVolume
OpVolume(const std::string &field_name)
Definition: hello_world.cpp:65
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
OpRow::OpRow
OpRow(const std::string &field_name)
Definition: hello_world.cpp:23
MoFEM::PipelineManager::getOpDomainRhsPipeline
boost::ptr_deque< UserDataOperator > & getOpDomainRhsPipeline()
Get the Op Domain Rhs Pipeline object.
Definition: PipelineManager.hpp:773
MoFEM::PipelineManager
PipelineManager interface.
Definition: PipelineManager.hpp:24
BasicFiniteElements.hpp
MoFEM::CoreTmp< 0 >::Finalize
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:112
MoFEM::Simple
Simple interface for fast problem set-up.
Definition: Simple.hpp:27
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:1975
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:535
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
OpVolumeSide::OpVolumeSide
OpVolumeSide(const std::string &field_name)
Definition: hello_world.cpp:115
MoFEM::ForcesAndSourcesCore::UserDataOperator
Definition: ForcesAndSourcesCore.hpp:549
MoFEM::FaceElementForcesAndSourcesCore::UserDataOperator
default operator for TRI element
Definition: FaceElementForcesAndSourcesCore.hpp:94
simple
void simple(double P1[], double P2[], double P3[], double c[], const int N)
Definition: acoustic.cpp:69
OpVolumeSide
Definition: hello_world.cpp:113
OpRow
Definition: hello_world.cpp:22
convert.type
type
Definition: convert.py:64
MOFEM_LOG_SYNCHRONISE
#define MOFEM_LOG_SYNCHRONISE(comm)
Synchronise "SYNC" channel.
Definition: LogManager.hpp:345
MoFEM::EntitiesFieldData::EntData::getIndices
const VectorInt & getIndices() const
Get global indices of dofs on entity.
Definition: EntitiesFieldData.hpp:1201
MoFEM::DMRegister_MoFEM
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition: DMMoFEM.cpp:47
OpFace
Definition: boundary_marker.cpp:19
MoFEM::FaceElementForcesAndSourcesCore
Face finite element.
Definition: FaceElementForcesAndSourcesCore.hpp:23
OpRowCol::doWork
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: hello_world.cpp:47
MoFEM::VolumeElementForcesAndSourcesCore
Volume finite element base.
Definition: VolumeElementForcesAndSourcesCore.hpp:26
MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator
Definition: VolumeElementForcesAndSourcesCore.hpp:108
field_name
constexpr auto field_name
Definition: poisson_2d_homogeneous.cpp:13
OpRowCol::OpRowCol
OpRowCol(const std::string row_field, const std::string col_field, const bool symm)
Definition: hello_world.cpp:43
MoFEM::ForcesAndSourcesCore
structure to get information form mofem into EntitiesFieldData
Definition: ForcesAndSourcesCore.hpp:22
MoFEM::CoreTmp< 0 >::Initialize
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
Definition: Core.cpp:72
MOFEM_LOG
#define MOFEM_LOG(channel, severity)
Log.
Definition: LogManager.hpp:308
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:372
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1094
OpVolume::doWork
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Definition: hello_world.cpp:69
main
int main(int argc, char *argv[])
Definition: hello_world.cpp:130
OpVolumeSide::doWork
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Definition: hello_world.cpp:118
OpFaceSide
Definition: hello_world.cpp:94
OpRow::doWork
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Definition: hello_world.cpp:25
AINSWORTH_LEGENDRE_BASE
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
OpFace::OpFace
OpFace(const std::string &field_name)
Definition: hello_world.cpp:81
MoFEM::VolumeElementForcesAndSourcesCoreOnSide
Base volume element used to integrate on skeleton.
Definition: VolumeElementForcesAndSourcesCoreOnSide.hpp:22
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
OpVolume
Definition: simple_interface.cpp:18
OpFaceSide::feSidePtr
boost::shared_ptr< VolumeElementForcesAndSourcesCoreOnSide > & feSidePtr
Definition: hello_world.cpp:95