v0.13.1
Simple.hpp
Go to the documentation of this file.
1/** \file Simple.hpp
2 * \brief Header file for simple interface
3 * \ingroup mofem_simple_interface
4 *
5 * Make simplified interface, to speedup problem setup and analysts.
6 * See discussion here
7 * <a
8 * href=https://groups.google.com/d/msg/mofem-group/Vkc00aia4dU/o9RF3ZmPAAAJ>link
9 * to google groups</a>
10 *
11 */
12
13#ifndef __SIMPLE_HPP__
14#define __SIMPLE_HPP__
15
16#include "UnknownInterface.hpp"
17
18namespace MoFEM {
19
20template <int DIM> struct ParentFiniteElementAdjacencyFunction;
21
22/**
23 * \brief Simple interface for fast problem set-up
24 * \ingroup mofem_simple_interface
25 */
26struct Simple : public UnknownInterface {
27
28 MoFEMErrorCode query_interface(boost::typeindex::type_index type_index,
29 UnknownInterface **iface) const;
30
31 Simple(const MoFEM::Core &core);
32 virtual ~Simple() = default;
33
34 /**
35 * \brief get options
36 * @return error code
37 */
39
40 /**
41 * \brief Load mesh file
42 * @param options file load options
43 * @param mesh_file_name file name if not set default or set by command line
44 * is used.
45 * @return error code
46 */
47 MoFEMErrorCode loadFile(const std::string options,
48 const std::string mesh_file_name);
49 /**
50 * \brief Load mesh file with parallel options if number of cores > 1
51 * @param mesh_file_name file name if not set default or set by command line
52 * is used.
53 * @return error code
54 */
55 MoFEMErrorCode loadFile(const std::string mesh_file_name = "");
56 /**
57 * \brief Add field on domain
58 * @param name name of the filed
59 * @param space space (L2,H1,Hdiv,Hcurl)
60 * @param base approximation base, see FieldApproximationBase
61 * @param nb_of_coefficients number of field coefficients
62 * @param tag_type type of the tag MB_TAG_DENSE or MB_TAG_SPARSE
63 * (DENSE is faster and uses less memory, SPARSE is more flexible if you
64 * define field on subdomains)
65 * @param bh if MF_EXCL throws error if field exits, MF_ZERO
66 * no error if field exist
67 * @param verb verbosity level
68 * @return error code
69 */
71 addDomainField(const std::string &name, const FieldSpace space,
72 const FieldApproximationBase base,
73 const FieldCoefficientsNumber nb_of_coefficients,
74 const TagType tag_type = MB_TAG_SPARSE,
75 const enum MoFEMTypes bh = MF_ZERO, int verb = -1);
76
77 /**
78 * \brief Add field on boundary
79 * @param name name of the filed
80 * @param space space (L2,H1,Hdiv,Hcurl)
81 * @param base approximation base, see FieldApproximationBase
82 * @param nb_of_coefficients number of field coefficients
83 * @param tag_type type of the tag MB_TAG_DENSE or MB_TAG_SPARSE
84 * (DENSE is faster and uses less memory, SPARSE is more flexible if you
85 * define field on subdomains)
86 * @param bh if MF_EXCL throws error if field exits, MF_ZERO
87 * no error if field exist
88 * @param verb verbosity level
89 * @return error code
90 */
92 addBoundaryField(const std::string &name, const FieldSpace space,
93 const FieldApproximationBase base,
94 const FieldCoefficientsNumber nb_of_coefficients,
95 const TagType tag_type = MB_TAG_SPARSE,
96 const enum MoFEMTypes bh = MF_ZERO, int verb = -1);
97
98 /**
99 * \brief Add field on skeleton
100 * @param name name of the filed
101 * @param space space (L2,H1,Hdiv,Hcurl)
102 * @param base approximation base, see FieldApproximationBase
103 * @param nb_of_coefficients number of field coefficients
104 * @param tag_type type of the tag MB_TAG_DENSE or MB_TAG_SPARSE
105 * (DENSE is faster and uses less memory, SPARSE is more flexible if you
106 * define field on subdomains)
107 * @param bh if MF_EXCL throws error if field exits, MF_ZERO
108 * no error if field exist
109 * @param verb verbosity level
110 * @return error code
111 */
113 addSkeletonField(const std::string &name, const FieldSpace space,
114 const FieldApproximationBase base,
115 const FieldCoefficientsNumber nb_of_coefficients,
116 const TagType tag_type = MB_TAG_SPARSE,
117 const enum MoFEMTypes bh = MF_ZERO, int verb = -1);
118
119 /**
120 * \brief Add field on domain
121 * @param name name of the filed
122 * @param space space (L2,H1,Hdiv,Hcurl)
123 * @param base approximation base, see FieldApproximationBase
124 * @param nb_of_coefficients number of field coefficients
125 * @param tag_type type of the tag MB_TAG_DENSE or MB_TAG_SPARSE
126 * (DENSE is faster and uses less memory, SPARSE is more flexible if you
127 * define field on subdomains)
128 * @param bh if MF_EXCL throws error if field exits, MF_ZERO
129 * no error if field exist
130 * @param verb verbosity level
131 * @return error code
132 */
133 MoFEMErrorCode addDataField(const std::string &name, const FieldSpace space,
134 const FieldApproximationBase base,
135 const FieldCoefficientsNumber nb_of_coefficients,
136 const TagType tag_type = MB_TAG_SPARSE,
137 const enum MoFEMTypes bh = MF_ZERO,
138 int verb = -1);
139
140 /**
141 * @brief Remove field form domain
142 *
143 * @param name
144 * @return MoFEMErrorCode
145 */
146 MoFEMErrorCode removeDomainField(const std::string &name);
147
148 /**
149 * @brief Remove field form boundary
150 *
151 * @param name
152 * @return MoFEMErrorCode
153 */
154 MoFEMErrorCode removeBoundaryField(const std::string &name);
155
156 /**
157 * @brief Remove field form skeleton
158 *
159 * @param name
160 * @return MoFEMErrorCode
161 */
162 MoFEMErrorCode removeSkeletonField(const std::string &name);
163
164 /**
165 * \brief Define finite elements
166 * @return Error code
167 */
169
170 /**
171 * \brief define problem
172 * @return error code
173 */
174 MoFEMErrorCode defineProblem(const PetscBool is_partitioned = PETSC_TRUE);
175
176 /**
177 * \brief Set field order
178 * @param std::field_name field name
179 * @param order order
180 * @param range of entities to which order is set (If null it sat to all
181 * entities)
182 * @return error code
183 */
184 MoFEMErrorCode setFieldOrder(const std::string field_name, const int order,
185 const Range *ents = NULL);
186
187 /**
188 * \brief Build fields
189 * @return error code
190 */
192
193 /**
194 * \brief Build finite elements
195 * @return error code
196 */
198
199 /**
200 * @brief Set the skeleton adjacency object
201 *
202 * @param dim
203 * @param fe_name
204 * @return MoFEMErrorCode
205 */
206 MoFEMErrorCode setSkeletonAdjacency(int dim = -1, std::string fe_name = "");
207
208 /**
209 * \brief Build problem
210 * @return error code
211 */
213
214 /**
215 * \brief Setup problem
216 * @return error code
217 */
218 MoFEMErrorCode setUp(const PetscBool is_partitioned = PETSC_TRUE);
219
220 /**
221 * @brief Rebuild internal MoFEM data structures
222 *
223 * Call this function after you add field or remove it.
224 *
225 * \note If you add field, or remove it, finite element and problem needs to
226 * be rebuild. However DM can remain the same.
227 *
228 * @return MoFEMErrorCode
229 */
230 MoFEMErrorCode reSetUp(bool only_dm = false);
231
232 /**
233 * \brief Get DM
234 * @param dm discrete manager
235 * @return error code
236 */
237 MoFEMErrorCode getDM(DM *dm);
238
239 /**
240 * @brief Return smart DM object
241 *
242 * \code
243 * {
244 * auto dm = simple_interface->getDM();
245 *
246 * // ...
247 *
248 * // dm is automatically destroyed when out of the scope
249 * }
250 * \endcode
251 *
252 * @return SmartPetscObj<DM>
253 */
254 inline SmartPetscObj<DM> getDM() { return dM; }
255
256 /**
257 * @brief Get the problem dimension
258 *
259 * Problem dimension is determined by highest dimension of entities on the
260 * mesh.
261 *
262 * @return int
263 */
264 inline int getDim() const { return dIm; }
265
266 /**
267 * @brief Set the problem dimension
268 *
269 * @return int
270 */
271 void setDim(int dim) { dIm = dim; };
272
273 /**
274 * @brief Get the MeshSet object
275 *
276 * @return EntityHandle&
277 */
278 inline EntityHandle &getMeshSet() { return meshSet; }
279
280 /**
281 * @brief Get the BoundaryMeshSet object
282 *
283 * @return EntityHandle&
284 */
285 inline EntityHandle &getBoundaryMeshSet() { return boundaryMeshset; }
286
287 /**
288 * @brief Get the SkeletonMeshSet object
289 *
290 * @return EntityHandle&
291 */
292 inline EntityHandle &getSkeletonMeshSet() { return skeletonMeshset; }
293
294 /**
295 * @brief Get the BitRefLevel
296 *
297 * @return BitRefLevel
298 */
299 inline BitRefLevel &getBitRefLevel() { return bitLevel; }
300
301 /**
302 * @brief Get the BitRefLevel
303 *
304 * @return BitRefLevel
305 */
307
308 /**
309 * @brief Get the Domain FE Name
310 *
311 * @return const std::string
312 */
313 inline const std::string getDomainFEName() const { return domainFE; }
314
315 /**
316 * @brief Get the Boundary FE Name
317 *
318 * @return const std::string
319 */
320 inline const std::string getBoundaryFEName() const { return boundaryFE; }
321
322 /**
323 * @brief Get the Skeleton FE Name
324 *
325 * @return const std::string
326 */
327 inline const std::string getSkeletonFEName() const { return skeletonFE; }
328
329 /**
330 * @brief Get the Problem Name
331 *
332 * @return const std::string
333 */
334 inline const std::string getProblemName() const { return nameOfProblem; }
335
336 /**
337 * @brief Get the Domain FE Name
338 *
339 * @return std::string&
340 */
341 inline std::string &getDomainFEName() { return domainFE; }
342
343 /**
344 * @brief Get the Boundary FE Name
345 *
346 * @return std::string&
347 */
348 inline std::string &getBoundaryFEName() { return boundaryFE; }
349
350 /**
351 * @brief Get the Skeleton FE Name
352 *
353 * @return std::string&
354 */
355 inline std::string &getSkeletonFEName() { return skeletonFE; }
356
357 /**
358 * @brief Get the Problem Name
359 *
360 * @return std::string&
361 */
362 inline std::string &getProblemName() { return nameOfProblem; }
363
364 /**
365 * @brief Get the Other Finite Elements
366 *
367 * User can create finite elements using directly core interface and
368 * and them to simple problem by this function
369 *
370 * @return std::vector<std::string>&
371 */
372 inline std::vector<std::string> &getOtherFiniteElements() { return otherFEs; }
373
374 /**
375 * @brief Delete dm
376 *
377 * @return MoFEMErrorCode
378 */
380
381 /**
382 * @brief Delete finite elements
383 *
384 * @return MoFEMErrorCode
385 */
387
388 /**
389 * @brief Get the addSkeletonFE
390 *
391 * If variable set to true, skeleton element is created regardless field on
392 * skelton is added or not.
393 *
394 * @return true
395 * @return false
396 */
397 bool &getAddSkeletonFE() { return addSkeletonFE; }
398
399 /**
400 * @brief Get the addSkeletonFE
401 *
402 * If variable set to true, boundary element is created regardless field on
403 * skelton is added or not.
404 *
405 * @return true
406 * @return false
407 */
408 bool &getAddBoundaryFE() { return addBoundaryFE; }
409
410 /**
411 * @brief Get the addParentAdjacencies
412 *
413 * If set true add parent adjacencies
414 *
415 * @return true
416 * @return false
417 */
419
420 /**
421 * @brief bit ref level for parent
422 *
423 * @return auto&
424 */
425 auto &getBitAdjParent() { return bitAdjParent; }
426
427 /**
428 * @brief bit ref level for parent marent
429 *
430 * @return auto&
431 */
433
434 /**
435 * @brief bit ref level for parent
436 *
437 * @return auto&
438 */
439 auto &getBitAdjEnt() { return bitAdjEnt; }
440
441 /**
442 * @brief bit ref level for parent marent
443 *
444 * @return auto&
445 */
446 auto &getBitAdjEntMask() { return bitAdjEntMask; }
447
448 /**
449 * @brief add empty block to problem
450 *
451 * MatrixManager assumes that all blocks, i.e. all fields combinations are non
452 * zero. This is not always the case, to optimise code and reduce memory usage
453 * user can specifi which blocks are empty.
454 *
455 * @param row_field row filed name
456 * @param col_field col field name
457 * @return MoFEMErrorCode
458 */
459 MoFEMErrorCode addFieldToEmptyFieldBlocks(const std::string row_field,
460 const std::string col_field) const;
461
462
463private:
465
466 BitRefLevel bitLevel; ///< BitRefLevel of the probelm
467 BitRefLevel bitLevelMask; ///< BitRefLevel of the probelm
468
475
476 EntityHandle meshSet; ///< domain meshset
477 EntityHandle boundaryMeshset; ///< meshset with boundary
478 EntityHandle skeletonMeshset; ///< skeleton meshset with boundary
479
480 bool addSkeletonFE; ///< Add skeleton FE
481 bool addBoundaryFE; ///< Add boundary FE
482 bool addParentAdjacencies; ///< If set to true parent adjacencies are build
483
484 BitRefLevel bitAdjParent; ///< bit ref level for parent
485 BitRefLevel bitAdjParentMask; ///< bit ref level for parent marent
486 BitRefLevel bitAdjEnt; ///< bit ref level for parent
487 BitRefLevel bitAdjEntMask; ///< bit ref level for parent marent
488
489 std::vector<std::string> domainFields; ///< domain fields
490 std::vector<std::string> boundaryFields; ///< boundary fields
491 std::vector<std::string> skeletonFields; ///< fields on the skeleton
492 std::vector<std::string> dataFields; ///< Data fields
493 std::vector<std::string> noFieldFields; ///< NOFIELD field name
494 std::vector<std::string> noFieldDataFields; ///< NOFIELD field name
495
496 std::multimap<std::string, std::pair<int, Range>>
497 fieldsOrder; ///< fields order
498
499 std::string nameOfProblem; ///< problem name
500 std::string domainFE; ///< domain finite element
501 std::string boundaryFE; ///< boundary finite element
502 std::string skeletonFE; ///< skeleton finite element
503
504 std::vector<std::string> otherFEs; ///< Other finite elements
505
506 char meshFileName[255]; ///< mesh file name
507 int dIm; ///< dimension of problem
508
510 dM; ///< Discrete manager (interface to PETSc/MoFEM functions)
511
515
516 template <int DIM = -1>
517 MoFEMErrorCode setSkeletonAdjacency(std::string fe_name);
518
519 template <int DIM = -1> MoFEMErrorCode setParentAdjacency();
520
521 boost::shared_ptr<ParentFiniteElementAdjacencyFunction<3>>
523 boost::shared_ptr<ParentFiniteElementAdjacencyFunction<2>>
525 boost::shared_ptr<ParentFiniteElementAdjacencyFunction<1>>
527};
528
529} // namespace MoFEM
530
531#endif // __SIMPLE_HPP__
532
533/**
534 * \defgroup mofem_simple_interface Simple interface
535 * \brief Implementation of simple interface for fast problem set-up.
536 *
537 * \ingroup mofem
538 **/
MoFEM interface.
MoFEMTypes
Those types control how functions respond on arguments, f.e. error handling.
Definition: definitions.h:97
@ MF_ZERO
Definition: definitions.h:98
FieldApproximationBase
approximation base
Definition: definitions.h:58
FieldSpace
approximation spaces
Definition: definitions.h:82
const int dim
char mesh_file_name[255]
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
int FieldCoefficientsNumber
Number of field coefficients.
Definition: Types.hpp:27
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: MoFEM.hpp:24
constexpr auto field_name
Core (interface) class.
Definition: Core.hpp:82
Simple interface for fast problem set-up.
Definition: Simple.hpp:26
boost::shared_ptr< ParentFiniteElementAdjacencyFunction< 2 > > parentAdjFunctionDim2
Definition: Simple.hpp:524
EntityHandle boundaryMeshset
meshset with boundary
Definition: Simple.hpp:477
MoFEMErrorCode buildProblem()
Build problem.
Definition: Simple.cpp:735
std::string & getDomainFEName()
Get the Domain FE Name.
Definition: Simple.hpp:341
bool & getAddBoundaryFE()
Get the addSkeletonFE.
Definition: Simple.hpp:408
const std::string getBoundaryFEName() const
Get the Boundary FE Name.
Definition: Simple.hpp:320
std::vector< std::string > dataFields
Data fields.
Definition: Simple.hpp:492
MoFEMErrorCode setParentAdjacency()
Definition: Simple.cpp:177
std::vector< std::string > domainFields
domain fields
Definition: Simple.hpp:489
MoFEMErrorCode removeDomainField(const std::string &name)
Remove field form domain.
Definition: Simple.cpp:447
MoFEMErrorCode defineFiniteElements()
Define finite elements.
Definition: Simple.cpp:493
MoFEM::Core & cOre
Definition: Simple.hpp:464
BitRefLevel bitAdjEnt
bit ref level for parent
Definition: Simple.hpp:486
char meshFileName[255]
mesh file name
Definition: Simple.hpp:506
EntityHandle meshSet
domain meshset
Definition: Simple.hpp:476
auto & getBitAdjParent()
bit ref level for parent
Definition: Simple.hpp:425
int getDim() const
Get the problem dimension.
Definition: Simple.hpp:264
std::vector< std::string > otherFEs
Other finite elements.
Definition: Simple.hpp:504
MoFEMErrorCode buildFiniteElements()
Build finite elements.
Definition: Simple.cpp:709
MoFEMErrorCode addDataField(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_ZERO, int verb=-1)
Add field on domain.
Definition: Simple.cpp:430
MoFEMErrorCode exchangeGhostCells()
Definition: Simple.cpp:931
std::multimap< std::string, std::pair< int, Range > > fieldsOrder
fields order
Definition: Simple.hpp:497
PetscLogEvent MOFEM_EVENT_SimpleBuildFields
Definition: Simple.hpp:471
const std::string getSkeletonFEName() const
Get the Skeleton FE Name.
Definition: Simple.hpp:327
PetscLogEvent MOFEM_EVENT_SimpleLoadMesh
Definition: Simple.hpp:470
MoFEMErrorCode addDomainField(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_ZERO, int verb=-1)
Add field on domain.
Definition: Simple.cpp:374
MoFEMErrorCode reSetUp(bool only_dm=false)
Rebuild internal MoFEM data structures.
Definition: Simple.cpp:772
std::vector< std::string > noFieldFields
NOFIELD field name.
Definition: Simple.hpp:493
bool addSkeletonFE
Add skeleton FE.
Definition: Simple.hpp:480
BitRefLevel bitAdjParent
bit ref level for parent
Definition: Simple.hpp:484
bool addBoundaryFE
Add boundary FE.
Definition: Simple.hpp:481
std::vector< std::string > & getOtherFiniteElements()
Get the Other Finite Elements.
Definition: Simple.hpp:372
MoFEMErrorCode addFieldToEmptyFieldBlocks(const std::string row_field, const std::string col_field) const
add empty block to problem
Definition: Simple.cpp:841
std::string skeletonFE
skeleton finite element
Definition: Simple.hpp:502
std::vector< std::string > noFieldDataFields
NOFIELD field name.
Definition: Simple.hpp:494
EntityHandle & getSkeletonMeshSet()
Get the SkeletonMeshSet object.
Definition: Simple.hpp:292
std::string boundaryFE
boundary finite element
Definition: Simple.hpp:501
std::string nameOfProblem
problem name
Definition: Simple.hpp:499
void setDim(int dim)
Set the problem dimension.
Definition: Simple.hpp:271
boost::shared_ptr< ParentFiniteElementAdjacencyFunction< 3 > > parentAdjFunctionDim3
Definition: Simple.hpp:522
MoFEMErrorCode getOptions()
get options
Definition: Simple.cpp:290
virtual ~Simple()=default
std::string domainFE
domain finite element
Definition: Simple.hpp:500
std::vector< std::string > boundaryFields
boundary fields
Definition: Simple.hpp:490
std::string & getBoundaryFEName()
Get the Boundary FE Name.
Definition: Simple.hpp:348
BitRefLevel bitAdjEntMask
bit ref level for parent marent
Definition: Simple.hpp:487
auto & getBitAdjEnt()
bit ref level for parent
Definition: Simple.hpp:439
SmartPetscObj< DM > dM
Discrete manager (interface to PETSc/MoFEM functions)
Definition: Simple.hpp:510
MoFEMErrorCode buildFields()
Build fields.
Definition: Simple.cpp:597
SmartPetscObj< DM > getDM()
Return smart DM object.
Definition: Simple.hpp:254
MoFEMErrorCode loadFile(const std::string options, const std::string mesh_file_name)
Load mesh file.
Definition: Simple.cpp:304
BitRefLevel bitLevel
BitRefLevel of the probelm.
Definition: Simple.hpp:466
MoFEMErrorCode setFieldOrder(const std::string field_name, const int order, const Range *ents=NULL)
Set field order.
Definition: Simple.cpp:589
MoFEMErrorCode addBoundaryField(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_ZERO, int verb=-1)
Add field on boundary.
Definition: Simple.cpp:392
MoFEMErrorCode defineProblem(const PetscBool is_partitioned=PETSC_TRUE)
define problem
Definition: Simple.cpp:565
BitRefLevel & getBitRefLevelMask()
Get the BitRefLevel.
Definition: Simple.hpp:306
std::vector< std::string > skeletonFields
fields on the skeleton
Definition: Simple.hpp:491
MoFEMErrorCode addSkeletonField(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_ZERO, int verb=-1)
Add field on skeleton.
Definition: Simple.cpp:410
EntityHandle & getBoundaryMeshSet()
Get the BoundaryMeshSet object.
Definition: Simple.hpp:285
EntityHandle & getMeshSet()
Get the MeshSet object.
Definition: Simple.hpp:278
MoFEMErrorCode removeSkeletonField(const std::string &name)
Remove field form skeleton.
Definition: Simple.cpp:478
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
Definition: Simple.cpp:9
auto & getBitAdjEntMask()
bit ref level for parent marent
Definition: Simple.hpp:446
bool & getAddSkeletonFE()
Get the addSkeletonFE.
Definition: Simple.hpp:397
MoFEMErrorCode setUp(const PetscBool is_partitioned=PETSC_TRUE)
Setup problem.
Definition: Simple.cpp:749
MoFEMErrorCode removeBoundaryField(const std::string &name)
Remove field form boundary.
Definition: Simple.cpp:463
std::string & getSkeletonFEName()
Get the Skeleton FE Name.
Definition: Simple.hpp:355
const std::string getProblemName() const
Get the Problem Name.
Definition: Simple.hpp:334
const std::string getDomainFEName() const
Get the Domain FE Name.
Definition: Simple.hpp:313
int dIm
dimension of problem
Definition: Simple.hpp:507
MoFEMErrorCode setSkeletonAdjacency(int dim=-1, std::string fe_name="")
Set the skeleton adjacency object.
Definition: Simple.cpp:253
MoFEMErrorCode deleteDM()
Delete dm.
Definition: Simple.cpp:818
std::string & getProblemName()
Get the Problem Name.
Definition: Simple.hpp:362
PetscLogEvent MOFEM_EVENT_SimpleSetUP
Definition: Simple.hpp:469
PetscLogEvent MOFEM_EVENT_SimpleKSPSolve
Definition: Simple.hpp:474
MoFEMErrorCode createBoundaryMeshset()
Definition: Simple.cpp:850
PetscLogEvent MOFEM_EVENT_SimpleBuildFiniteElements
Definition: Simple.hpp:472
BitRefLevel bitLevelMask
BitRefLevel of the probelm.
Definition: Simple.hpp:467
auto & getBitAdjParentMask()
bit ref level for parent marent
Definition: Simple.hpp:432
BitRefLevel bitAdjParentMask
bit ref level for parent marent
Definition: Simple.hpp:485
PetscLogEvent MOFEM_EVENT_SimpleBuildProblem
Definition: Simple.hpp:473
boost::shared_ptr< ParentFiniteElementAdjacencyFunction< 1 > > parentAdjFunctionDim1
Definition: Simple.hpp:526
BitRefLevel & getBitRefLevel()
Get the BitRefLevel.
Definition: Simple.hpp:299
bool & getParentAdjacencies()
Get the addParentAdjacencies.
Definition: Simple.hpp:418
MoFEMErrorCode deleteFiniteElements()
Delete finite elements.
Definition: Simple.cpp:829
Simple(const MoFEM::Core &core)
Definition: Simple.cpp:272
EntityHandle skeletonMeshset
skeleton meshset with boundary
Definition: Simple.hpp:478
MoFEMErrorCode createSkeletonMeshset()
Definition: Simple.cpp:897
bool addParentAdjacencies
If set to true parent adjacencies are build.
Definition: Simple.hpp:482
base class for all interface classes