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