v0.14.0
Loading...
Searching...
No Matches
forces_and_sources_testing_users_base.cpp
Go to the documentation of this file.
1/** \file forces_and_sources_testing_users_base.cpp
2 * \example forces_and_sources_testing_users_base.cpp
3 *
4 * Primarily this is used for testing if the code can handle user base. It is
5 * also, an example of how to build and use user approximation base. This is a
6 * test, so we used RT base by Demkowicz recipe.
7 *
8 * Note that triple defines approximation element; element entity type,
9 * approximation space and approximation base. Entity type determines the
10 * integration method; approximation space determines the adjacency of the
11 * matrix and approximation base determines together with space the regularity
12 * of approximation.
13 *
14 */
15
16
17
18#include <MoFEM.hpp>
19#include <Hdiv.hpp>
20
21namespace bio = boost::iostreams;
22using bio::stream;
23using bio::tee_device;
24
25using namespace MoFEM;
26
27// used for profiling
28static PetscBool quiet = PETSC_FALSE;
29static PetscBool base_cache = PETSC_FALSE;
30
31static char help[] = "...\n\n";
32
33/**
34 * @brief Class used to calculate base functions at integration points
35 *
36 */
38
41
42 /**
43 * @brief Return interface to this class when one ask for for tetrahedron,
44 * otherisw return interface class for generic class.
45 *
46 * @param iface interface class
47 * @return MoFEMErrorCode
48 */
49 MoFEMErrorCode query_interface(boost::typeindex::type_index type_index,
50 UnknownInterface **iface) const {
52 *iface = const_cast<SomeUserPolynomialBase *>(this);
54 }
55
56 /**
57 * @brief Calculate base functions at intergeneration points
58 *
59 * @param pts
60 * @param ctx_ptr
61 * @return MoFEMErrorCode
62 */
64 boost::shared_ptr<BaseFunctionCtx> ctx_ptr) {
66
68
69 int nb_gauss_pts = pts.size2();
70 if (!nb_gauss_pts) {
72 }
73
74 if (pts.size1() < 3) {
75 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
76 "Wrong dimension of pts, should be at least 3 rows with "
77 "coordinates");
78 }
79
80 switch (cTx->sPace) {
81 case HDIV:
83 break;
84 default:
85 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Not yet implemented");
86 }
87
89 }
90
91private:
93
95
98
99 const FieldApproximationBase base = cTx->bAse;
100 // This should be used only in case USER_BASE is selected
101 if (cTx->bAse != USER_BASE) {
102 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
103 "Wrong base, should be USER_BASE");
104 }
105
106 // This is example, simply use Demkowicz HDiv base to generate base
107 // functions
108
109 EntitiesFieldData &data = cTx->dAta;
110 int nb_gauss_pts = pts.size2();
111
112 // calculate shape functions, i.e. barycentric coordinates
113 shapeFun.resize(nb_gauss_pts, 4, false);
114 CHKERR ShapeMBTET(&*shapeFun.data().begin(), &pts(0, 0), &pts(1, 0),
115 &pts(2, 0), nb_gauss_pts);
116 // derivatives of shape functions
117 double diff_shape_fun[12];
118 CHKERR ShapeDiffMBTET(diff_shape_fun);
119
120 int volume_order = data.dataOnEntities[MBTET][0].getOrder();
121
122 int p_f[4];
123 double *phi_f[4];
124 double *diff_phi_f[4];
125
126 // Calculate base function on tet faces
127 for (int ff = 0; ff != 4; ff++) {
128 int face_order = data.dataOnEntities[MBTRI][ff].getOrder();
129 int order = volume_order > face_order ? volume_order : face_order;
130 data.dataOnEntities[MBTRI][ff].getN(base).resize(
131 nb_gauss_pts, 3 * NBFACETRI_DEMKOWICZ_HDIV(order), false);
132 data.dataOnEntities[MBTRI][ff].getDiffN(base).resize(
133 nb_gauss_pts, 9 * NBFACETRI_DEMKOWICZ_HDIV(order), false);
134 p_f[ff] = order;
135 phi_f[ff] = &*data.dataOnEntities[MBTRI][ff].getN(base).data().begin();
136 diff_phi_f[ff] =
137 &*data.dataOnEntities[MBTRI][ff].getDiffN(base).data().begin();
139 continue;
141 &data.facesNodes(ff, 0), order, &*shapeFun.data().begin(),
142 diff_shape_fun, phi_f[ff], diff_phi_f[ff], nb_gauss_pts, 4);
143 }
144
145 // Calculate base functions in tet interior
146 if (NBVOLUMETET_DEMKOWICZ_HDIV(volume_order) > 0) {
147 data.dataOnEntities[MBTET][0].getN(base).resize(
148 nb_gauss_pts, 3 * NBVOLUMETET_DEMKOWICZ_HDIV(volume_order), false);
149 data.dataOnEntities[MBTET][0].getDiffN(base).resize(
150 nb_gauss_pts, 9 * NBVOLUMETET_DEMKOWICZ_HDIV(volume_order), false);
151 double *phi_v = &*data.dataOnEntities[MBTET][0].getN(base).data().begin();
152 double *diff_phi_v =
153 &*data.dataOnEntities[MBTET][0].getDiffN(base).data().begin();
155 volume_order, &*shapeFun.data().begin(), diff_shape_fun, p_f, phi_f,
156 diff_phi_f, phi_v, diff_phi_v, nb_gauss_pts);
157 }
158
159 // Set size of face base correctly
160 for (int ff = 0; ff != 4; ff++) {
161 int face_order = data.dataOnEntities[MBTRI][ff].getOrder();
162 data.dataOnEntities[MBTRI][ff].getN(base).resize(
163 nb_gauss_pts, 3 * NBFACETRI_DEMKOWICZ_HDIV(face_order), true);
164 data.dataOnEntities[MBTRI][ff].getDiffN(base).resize(
165 nb_gauss_pts, 9 * NBFACETRI_DEMKOWICZ_HDIV(face_order), true);
166 }
167
169 }
170};
171
172int main(int argc, char *argv[]) {
173
174 // Initialise MoFEM, MPI and petsc
175 MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
176
177 try {
178
179 // create moab
180 moab::Core mb_instance;
181 // get interface to moab databse
182 moab::Interface &moab = mb_instance;
183
184 // get file
185 PetscBool flg = PETSC_TRUE;
186 char mesh_file_name[255];
187#if PETSC_VERSION_GE(3, 6, 4)
188 CHKERR PetscOptionsGetString(PETSC_NULL, "", "-my_file", mesh_file_name,
189 255, &flg);
190#else
191 CHKERR PetscOptionsGetString(PETSC_NULL, PETSC_NULL, "-my_file",
192 mesh_file_name, 255, &flg);
193#endif
194 if (flg != PETSC_TRUE) {
195 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
196 "*** ERROR -my_file (MESH FILE NEEDED)");
197 }
198
199 // create MoFEM database
200 MoFEM::Core core(moab);
201 // get interface to moab database
202 MoFEM::Interface &m_field = core;
203
204 // load mesh file
205 const char *option;
206 option = "";
207 CHKERR moab.load_file(mesh_file_name, 0, option);
208
209 // set bit refinement level
210 CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(
211 0, 3, BitRefLevel().set(0));
212
213 // Create fields, field "FIELD_CGG" has user base, it means that recipe how
214 // to construct approximation is provided by user. Is set that user provided
215 // base is in h-div space.
216 CHKERR m_field.add_field("FILED_CGG", HDIV, USER_BASE, 1);
217 CHKERR m_field.add_field("FILED_RT", HDIV, DEMKOWICZ_JACOBI_BASE, 1);
218
219 // get access to "FIELD_CGG" data structure
220 auto field_ptr = m_field.get_field_structure("FILED_CGG");
221 // get table associating number of dofs to entities depending on
222 // approximation order set on those entities.
223 auto field_order_table =
224 const_cast<Field *>(field_ptr)->getFieldOrderTable();
225
226 // function set zero number of dofs
227 auto get_cgg_bubble_order_zero = [](int p) { return 0; };
228 // function set non-zero number of dofs on tetrahedrons
229 auto get_cgg_bubble_order_face = [](int p) {
230 return NBFACETRI_DEMKOWICZ_HDIV(p);
231 };
232 auto get_cgg_bubble_order_tet = [](int p) {
234 };
235 field_order_table[MBVERTEX] = get_cgg_bubble_order_zero;
236 field_order_table[MBEDGE] = get_cgg_bubble_order_zero;
237 field_order_table[MBTRI] = get_cgg_bubble_order_face;
238 field_order_table[MBTET] = get_cgg_bubble_order_tet;
239 const_cast<Field *>(field_ptr)->rebuildDofsOrderMap();
240
241 auto &dof_order_map = field_ptr->getDofOrderMap(MBTET);
242 for(auto d = 0; d!=10; ++d) {
243 MOFEM_LOG("WORLD", Sev::noisy) << "dof " << dof_order_map[d];
244 }
245
246 CHKERR m_field.add_finite_element("FE");
247
248 // define rows/cols and element data
249 CHKERR m_field.modify_finite_element_add_field_row("FE", "FILED_CGG");
250 CHKERR m_field.modify_finite_element_add_field_col("FE", "FILED_CGG");
251 CHKERR m_field.modify_finite_element_add_field_data("FE", "FILED_CGG");
252 CHKERR m_field.modify_finite_element_add_field_row("FE", "FILED_RT");
253 CHKERR m_field.modify_finite_element_add_field_col("FE", "FILED_RT");
254 CHKERR m_field.modify_finite_element_add_field_data("FE", "FILED_RT");
255
256 // add problem
257 CHKERR m_field.add_problem("PROBLEM");
258
259 // set finite elements for problem
260 CHKERR m_field.modify_problem_add_finite_element("PROBLEM", "FE");
261 // set refinement level for problem
263 BitRefLevel().set(0));
264
265 // meshset consisting all entities in mesh
266 EntityHandle root_set = moab.get_root_set();
267 // add entities to field
268 CHKERR m_field.add_ents_to_field_by_type(root_set, MBTET, "FILED_CGG");
269 CHKERR m_field.add_ents_to_field_by_type(root_set, MBTET, "FILED_RT");
270 // add entities to finite element
271 CHKERR m_field.add_ents_to_finite_element_by_type(root_set, MBTET, "FE");
272
273 // set app. order
274 int order = 3;
275 CHKERR PetscOptionsGetInt(PETSC_NULL, "", "-order", &order, PETSC_NULL);
276 CHKERR m_field.set_field_order(root_set, MBTRI, "FILED_CGG", order);
277 CHKERR m_field.set_field_order(root_set, MBTET, "FILED_CGG", order);
278 CHKERR m_field.set_field_order(root_set, MBTRI, "FILED_RT", order);
279 CHKERR m_field.set_field_order(root_set, MBTET, "FILED_RT", order);
280
281 /****/
282 // build database
283 // build field
284 CHKERR m_field.build_fields();
285 // build finite elemnts
287 // build adjacencies
288 CHKERR m_field.build_adjacencies(BitRefLevel().set(0));
289
290 // build problem
291 CHKERR m_field.getInterface<ProblemsManager>()->buildProblem("PROBLEM",
292 true);
293 // dofs partitioning
294 CHKERR m_field.getInterface<ProblemsManager>()->partitionSimpleProblem(
295 "PROBLEM");
296 CHKERR m_field.getInterface<ProblemsManager>()->partitionFiniteElements(
297 "PROBLEM");
298 // what are ghost nodes, see Petsc Manual
299 CHKERR m_field.getInterface<ProblemsManager>()->partitionGhostDofs(
300 "PROBLEM");
301
302 typedef tee_device<std::ostream, std::ofstream> TeeDevice;
303 typedef stream<TeeDevice> TeeStream;
304
305 std::ofstream ofs("forces_and_sources_testing_users_base.txt");
306 TeeDevice my_tee(std::cout, ofs);
307 TeeStream my_split(my_tee);
308
309 /**
310 * Simple user data operator which main purpose is to print values
311 * of base functions at integration points.
312 *
313 */
314 struct MyOp1 : public VolumeElementForcesAndSourcesCore::UserDataOperator {
315
316 TeeStream &my_split;
317 MyOp1(const std::string &row_field, const std::string &col_field,
318 TeeStream &_my_split, char type)
319 : VolumeElementForcesAndSourcesCore::UserDataOperator(
320 row_field, col_field, type),
321 my_split(_my_split) {
322 sYmm = false;
323 }
324
325 MoFEMErrorCode doWork(int side, EntityType type,
328 if (data.getIndices().empty()) {
330 }
331 if (!quiet) {
332 my_split << rowFieldName << endl;
333 my_split << "side: " << side << " type: " << type << std::endl;
334 my_split << data << endl;
335 my_split << data.getN() << endl;
336 my_split << endl;
337 }
339 }
340
341 MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type,
342 EntityType col_type,
344 EntitiesFieldData::EntData &col_data) {
346 if (row_data.getIndices().empty())
348 if (col_data.getIndices().empty())
350 if (!quiet) {
351 my_split << rowFieldName << " : " << colFieldName << endl;
352 my_split << "row side: " << row_side << " row_type: " << row_type
353 << std::endl;
354 my_split << "col side: " << col_side << " col_type: " << col_type
355 << std::endl;
356 my_split << row_data.getIndices().size() << " : "
357 << col_data.getIndices().size() << endl;
358 my_split << endl;
359 }
361 }
362 };
363
364 // create finite element instance
366 // set class needed to construct user approximation base
368 boost::shared_ptr<BaseFunction>(new SomeUserPolynomialBase());
369
370 CHKERR PetscOptionsGetBool(PETSC_NULL, "", "-quiet", &quiet, PETSC_NULL);
371
372 // push user data operators
373 fe1.getOpPtrVector().push_back(
374 new MyOp1("FILED_CGG", "FILED_CGG", my_split,
375 ForcesAndSourcesCore::UserDataOperator::OPROW));
376 fe1.getOpPtrVector().push_back(
377 new MyOp1("FILED_CGG", "FILED_RT", my_split,
378 ForcesAndSourcesCore::UserDataOperator::OPROWCOL));
379
380 CHKERR PetscOptionsGetBool(PETSC_NULL, "", "-base_cache", &base_cache,
381 PETSC_NULL);
382
383 if (base_cache) {
384 if (!TetPolynomialBase::switchCacheBaseFace<HDIV>(DEMKOWICZ_JACOBI_BASE,
385 &fe1)) {
386 TetPolynomialBase::switchCacheBaseFace<HDIV>(DEMKOWICZ_JACOBI_BASE,
387 &fe1);
388 }
389 if (!TetPolynomialBase::switchCacheBaseInterior<HDIV>(
390 DEMKOWICZ_JACOBI_BASE, &fe1)) {
391 TetPolynomialBase::switchCacheBaseInterior<HDIV>(DEMKOWICZ_JACOBI_BASE,
392 &fe1);
393 }
394 }
395
396 // iterate over finite elements, and execute user data operators on each
397 // of them
398 CHKERR m_field.loop_finite_elements("PROBLEM", "FE", fe1);
399
400 if (base_cache) {
401 if (TetPolynomialBase::switchCacheBaseFace<HDIV>(DEMKOWICZ_JACOBI_BASE,
402 &fe1)) {
403 }
404 if (TetPolynomialBase::switchCacheBaseInterior<HDIV>(DEMKOWICZ_JACOBI_BASE,
405 &fe1)) {
406 };
407 }
408 }
410
412
413 return 0;
414}
Implementation of H-curl base function.
int main()
#define CATCH_ERRORS
Catch errors.
FieldApproximationBase
approximation base
Definition definitions.h:58
@ USER_BASE
user implemented approximation base
Definition definitions.h:68
@ DEMKOWICZ_JACOBI_BASE
Definition definitions.h:66
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
@ HDIV
field with continuous normal traction
Definition definitions.h:87
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ MOFEM_DATA_INCONSISTENCY
Definition definitions.h:31
@ MOFEM_NOT_IMPLEMENTED
Definition definitions.h:32
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
constexpr int order
PetscErrorCode ShapeDiffMBTET(double *diffN)
calculate derivatives of shape functions
Definition fem_tools.c:319
PetscErrorCode ShapeMBTET(double *N, const double *G_X, const double *G_Y, const double *G_Z, int DIM)
calculate shape functions
Definition fem_tools.c:306
tee_device< std::ostream, std::ofstream > TeeDevice
static PetscBool base_cache
static PetscBool quiet
virtual MoFEMErrorCode add_finite_element(const std::string &fe_name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
add finite element
virtual MoFEMErrorCode build_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elements.
virtual MoFEMErrorCode modify_finite_element_add_field_col(const std::string &fe_name, const std::string name_row)=0
set field col which finite element use
virtual MoFEMErrorCode add_ents_to_finite_element_by_type(const EntityHandle entities, const EntityType type, const std::string &name, const bool recursive=true)=0
add entities to finite element
virtual MoFEMErrorCode modify_finite_element_add_field_row(const std::string &fe_name, const std::string name_row)=0
set field row which finite element use
virtual MoFEMErrorCode modify_finite_element_add_field_data(const std::string &fe_name, const std::string name_field)=0
set finite element field data
virtual const Field * get_field_structure(const std::string &name, enum MoFEMTypes bh=MF_EXIST) const =0
get field structure
virtual MoFEMErrorCode build_fields(int verb=DEFAULT_VERBOSITY)=0
virtual MoFEMErrorCode set_field_order(const EntityHandle meshset, const EntityType type, const std::string &name, const ApproximationOrder order, int verb=DEFAULT_VERBOSITY)=0
Set order approximation of the entities in the field.
virtual MoFEMErrorCode add_ents_to_field_by_type(const Range &ents, const EntityType type, const std::string &name, int verb=DEFAULT_VERBOSITY)=0
Add entities to field meshset.
#define MOFEM_LOG(channel, severity)
Log.
virtual MoFEMErrorCode loop_finite_elements(const std::string problem_name, const std::string &fe_name, FEMethod &method, boost::shared_ptr< NumeredEntFiniteElement_multiIndex > fe_ptr=nullptr, MoFEMTypes bh=MF_EXIST, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr(), int verb=DEFAULT_VERBOSITY)=0
Make a loop over finite elements.
virtual MoFEMErrorCode add_problem(const std::string &name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add problem.
virtual MoFEMErrorCode modify_problem_ref_level_add_bit(const std::string &name_problem, const BitRefLevel &bit)=0
add ref level to problem
virtual MoFEMErrorCode modify_problem_add_finite_element(const std::string name_problem, const std::string &fe_name)=0
add finite element to problem, this add entities assigned to finite element to a particular problem
#define NBVOLUMETET_DEMKOWICZ_HDIV(P)
#define NBFACETRI_DEMKOWICZ_HDIV(P)
char mesh_file_name[255]
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition Types.hpp:40
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
MoFEMErrorCode Hdiv_Demkowicz_Face_MBTET_ON_FACE(int *faces_nodes, int p, double *N, double *diffN, double *phi_f, double *diff_phi_f, int gdim, int nb)
Definition Hdiv.cpp:634
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
PetscErrorCode PetscOptionsGetBool(PetscOptions *, const char pre[], const char name[], PetscBool *bval, PetscBool *set)
MoFEMErrorCode Hdiv_Demkowicz_Interior_MBTET(int p, double *N, double *diffN, int p_face[], double *phi_f[4], double *diff_phi_f[4], double *phi_v, double *diff_phi_v, int gdim)
Definition Hdiv.cpp:780
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
Base class if inherited used to calculate base functions.
Managing BitRefLevels.
virtual MoFEMErrorCode build_adjacencies(const Range &ents, int verb=DEFAULT_VERBOSITY)=0
build adjacencies
virtual MoFEMErrorCode add_field(const std::string name, const FieldSpace space, const FieldApproximationBase base, const FieldCoefficientsNumber nb_of_coefficients, const TagType tag_type=MB_TAG_SPARSE, const enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add field.
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:112
Deprecated interface functions.
Class used to pass element data to calculate base functions on tet,triangle,edge.
const FieldApproximationBase bAse
Data on single entity (This is passed as argument to DataOperator::doWork)
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
const VectorInt & getIndices() const
Get global indices of dofs on entity.
data structure for finite element entity
MatrixInt facesNodes
nodes on finite element faces
std::array< boost::ptr_vector< EntData >, MBMAXTYPE > dataOnEntities
Provide data structure for (tensor) field approximation.
const std::array< ApproximationOrder, MAX_DOFS_ON_ENTITY > & getDofOrderMap(const EntityType type) const
get hash-map relating dof index on entity with its order
auto & getUserPolynomialBase()
Get the User Polynomial Base object.
boost::ptr_deque< UserDataOperator > & getOpPtrVector()
Use to push back operator for row operator.
Problem manager is used to build and partition problems.
base class for all interface classes
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Class used to calculate base functions at integration points.
MoFEMErrorCode getValueHdivForCGGBubble(MatrixDouble &pts)
MoFEMErrorCode getValue(MatrixDouble &pts, boost::shared_ptr< BaseFunctionCtx > ctx_ptr)
Calculate base functions at intergeneration points.
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
Return interface to this class when one ask for for tetrahedron, otherisw return interface class for ...
~SomeUserPolynomialBase()=default
SomeUserPolynomialBase()=default