v0.15.0
Loading...
Searching...
No Matches
Classes | Macros | Functions | Variables
hello_world.cpp File Reference
#include <MoFEM.hpp>

Go to the source code of this file.

Classes

struct  OpRow
 
struct  OpRowCol
 
struct  OpVolume
 
struct  OpFace
 
struct  OpFaceSide
 
struct  OpVolumeSide
 

Macros

#define HelloFunctionBegin
 

Functions

int main (int argc, char *argv[])
 

Variables

static char help [] = "...\n\n"
 

Macro Definition Documentation

◆ HelloFunctionBegin

#define HelloFunctionBegin
Value:
MOFEM_LOG_CHANNEL("SYNC"); \
MOFEM_LOG_FUNCTION(); \
MOFEM_LOG_TAG("WORLD", "HelloWorld");
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Examples
hello_world.cpp.

Definition at line 16 of file hello_world.cpp.

21 : public ForcesAndSourcesCore::UserDataOperator {
22 OpRow(const std::string &field_name)
24 MoFEMErrorCode doWork(int side, EntityType type,
27 if (type == MBVERTEX) {
28 // get number of evaluated element in the loop
29 MOFEM_LOG("SYNC", Sev::inform) << "**** " << getNinTheLoop() << " ****";
30 MOFEM_LOG("SYNC", Sev::inform) << "**** Operators ****";
31 }
32 MOFEM_LOG("SYNC", Sev::inform)
33 << "Hello Operator OpRow:"
34 << " field name " << rowFieldName << " side " << side << " type "
35 << CN::EntityTypeName(type) << " nb dofs on entity "
36 << data.getIndices().size();
38 }
39};
40
41struct OpRowCol : public ForcesAndSourcesCore::UserDataOperator {
42 OpRowCol(const std::string row_field, const std::string col_field,
43 const bool symm)
44 : ForcesAndSourcesCore::UserDataOperator(row_field, col_field, OPROWCOL,
45 symm) {}
46 MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type,
47 EntityType col_type,
51 MOFEM_LOG("SYNC", Sev::inform)
52 << "Hello Operator OpRowCol:"
53 << " row field name " << rowFieldName << " row side " << row_side
54 << " row type " << CN::EntityTypeName(row_type)
55 << " nb dofs on row entity " << row_data.getIndices().size() << " : "
56 << " col field name " << colFieldName << " col side " << col_side
57 << " col type " << CN::EntityTypeName(col_type)
58 << " nb dofs on col entity " << col_data.getIndices().size();
60 }
61};
62
63struct OpVolume : public VolumeElementForcesAndSourcesCore::UserDataOperator {
64 OpVolume(const std::string &field_name)
66 field_name, OPROW) {
67 }
68 MoFEMErrorCode doWork(int side, EntityType type,
71 if (type == MBVERTEX) {
72 MOFEM_LOG("SYNC", Sev::inform) << "Hello Operator OpVolume:"
73 << " volume " << getVolume();
74 }
76 }
77};
78
79struct OpFace : public FaceElementForcesAndSourcesCore::UserDataOperator {
80 OpFace(const std::string &field_name)
82 MoFEMErrorCode doWork(int side, EntityType type,
85 if (type == MBVERTEX) {
86 MOFEM_LOG("SYNC", Sev::inform) << "Hello Operator OpFace:"
87 << " normal " << getNormal();
88 }
90 }
91};
92
93struct OpFaceSide : public FaceElementForcesAndSourcesCore::UserDataOperator {
94 boost::shared_ptr<VolumeElementForcesAndSourcesCoreOnSide> &feSidePtr;
96 const std::string &field_name,
97 boost::shared_ptr<VolumeElementForcesAndSourcesCoreOnSide> &fe_side_ptr)
99 feSidePtr(fe_side_ptr) {}
100 MoFEMErrorCode doWork(int side, EntityType type,
102
104 if (type == MBVERTEX) {
105 MOFEM_LOG("SYNC", Sev::inform) << "Hello Operator OpSideFace";
106 CHKERR loopSideVolumes("dFE", *feSidePtr);
107 }
109 }
110};
111
112struct OpVolumeSide
113 : public VolumeElementForcesAndSourcesCoreOnSide::UserDataOperator {
114 OpVolumeSide(const std::string &field_name)
116 field_name, field_name, OPROW) {}
117 MoFEMErrorCode doWork(int side, EntityType type,
120 if (type == MBVERTEX) {
121 MOFEM_LOG("SYNC", Sev::inform)
122 << "Hello Operator OpVolumeSide:"
123 << " volume " << getVolume() << " normal " << getNormal();
124 }
126 }
127};
128
129int main(int argc, char *argv[]) {
130
131 // initialize petsc
132 MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
133
134 try {
135
136 // Register DM Manager
137 DMType dm_name = "DMMOFEM";
138 CHKERR DMRegister_MoFEM(dm_name);
139
140 // Create MoAB database
141 moab::Core moab_core;
142 moab::Interface &moab = moab_core;
143
144 // Create MoFEM database and link it to MoAB
145 MoFEM::Core mofem_core(moab);
146 MoFEM::Interface &m_field = mofem_core;
147
148 // Simple interface
149 Simple *simple;
150 CHKERR m_field.getInterface(simple);
151
152 // get options from command line
153 CHKERR simple->getOptions();
154 // load mesh file
155 CHKERR simple->loadFile();
156 // add fields
157 CHKERR simple->addDomainField("U", H1, AINSWORTH_LEGENDRE_BASE, 3);
158 CHKERR simple->addBoundaryField("L", H1, AINSWORTH_LEGENDRE_BASE, 3);
159 CHKERR simple->addSkeletonField("S", H1, AINSWORTH_LEGENDRE_BASE, 3);
160 // set fields order
161 CHKERR simple->setFieldOrder("U", 4);
162 CHKERR simple->setFieldOrder("L", 3);
163 CHKERR simple->setFieldOrder("S", 3);
164
165 // setup problem
166 CHKERR simple->setUp();
167
168 auto pipeline_mng = m_field.getInterface<PipelineManager>();
169
170 //! [set operator to the volume element instance]
171 pipeline_mng->getOpDomainRhsPipeline().push_back(new OpRow("U"));
172 pipeline_mng->getOpDomainRhsPipeline().push_back(new OpVolume("U"));
173 pipeline_mng->getOpDomainLhsPipeline().push_back(
174 new OpRowCol("U", "U", true));
175 //! [set operator to the volume element instance]
176
177 //! [set operator to the face element instance]
178 pipeline_mng->getOpBoundaryRhsPipeline().push_back(new OpRow("L"));
179 pipeline_mng->getOpBoundaryRhsPipeline().push_back(new OpFace("L"));
180 pipeline_mng->getOpBoundaryLhsPipeline().push_back(
181 new OpRowCol("U", "L", false));
182 //! [set operator to the face element instance]
183
184 //! [create skeleton side element and push operator to it]
185 boost::shared_ptr<VolumeElementForcesAndSourcesCoreOnSide> side_fe(
187 side_fe->getOpPtrVector().push_back(new OpVolumeSide("U"));
188 //! [create skeleton side element and push operator to it]
189
190 //! [set operator to the face element on skeleton instance]
191 pipeline_mng->getOpSkeletonRhsPipeline().push_back(new OpRow("S"));
192 pipeline_mng->getOpSkeletonRhsPipeline().push_back(
193 new OpFaceSide("S", side_fe));
194 //! [set operator to the face element on skeleton instance]
195
196 //! [executing finite elements]
197 CHKERR pipeline_mng->loopFiniteElements();
198 //! [executing finite elements]
199
201 }
203
204 // finish work cleaning memory, getting statistics, etc.
206
207 return 0;
208}
ForcesAndSourcesCore::UserDataOperator UserDataOperator
#define MOFEM_LOG_SYNCHRONISE(comm)
Synchronise "SYNC" channel.
void simple(double P1[], double P2[], double P3[], double c[], const int N)
Definition acoustic.cpp:69
int main()
#define CATCH_ERRORS
Catch errors.
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition definitions.h:60
@ H1
continuous field
Definition definitions.h:85
#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.
Definition DMMoFEM.cpp:43
boost::ptr_deque< UserDataOperator > & getOpDomainRhsPipeline()
Get the Op Domain Rhs Pipeline object.
#define MOFEM_LOG(channel, severity)
Log.
#define HelloFunctionBegin
static char help[]
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
constexpr auto field_name
virtual MPI_Comm & get_comm() const =0
Core (interface) class.
Definition Core.hpp:82
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
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition Core.cpp:118
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.
@ 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.
Definition Simple.hpp:27
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
boost::shared_ptr< VolumeElementForcesAndSourcesCoreOnSide > & feSidePtr
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)

Function Documentation

◆ main()

int main ( int  argc,
char *  argv[] 
)

[set operator to the volume element instance]

[set operator to the volume element instance]

[set operator to the face element instance]

[set operator to the face element instance]

[create skeleton side element and push operator to it]

[create skeleton side element and push operator to it]

[set operator to the face element on skeleton instance]

[set operator to the face element on skeleton instance]

[executing finite elements]

[executing finite elements]

Definition at line 130 of file hello_world.cpp.

130 {
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
202 }
204
205 // finish work cleaning memory, getting statistics, etc.
207
208 return 0;
209}

Variable Documentation

◆ help

char help[] = "...\n\n"
static

Definition at line 14 of file hello_world.cpp.