v0.14.0
Classes | Typedefs | Functions | Variables
hcurl_divergence_operator_2d_embaded_in_3d.cpp File Reference
#include <MoFEM.hpp>

Go to the source code of this file.

Classes

struct  OpDivergence
 
struct  OpFlux
 

Typedefs

using FaceEle = MoFEM::FaceElementForcesAndSourcesCore
 
using EdgeEle = MoFEM::EdgeElementForcesAndSourcesCore
 
using FaceEleOp = FaceEle::UserDataOperator
 
using EdgeEleOp = EdgeEle::UserDataOperator
 

Functions

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

Variables

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

Typedef Documentation

◆ EdgeEle

◆ EdgeEleOp

◆ FaceEle

◆ FaceEleOp

Function Documentation

◆ main()

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

Definition at line 76 of file hcurl_divergence_operator_2d_embaded_in_3d.cpp.

76  {
77 
78  MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
79 
80  try {
81 
82  moab::Core mb_instance;
83  moab::Interface &moab = mb_instance;
84 
85  PetscBool flg_file = PETSC_TRUE;
86  char mesh_file_name[255];
87  CHKERR PetscOptionsGetString(PETSC_NULL, "", "-my_file", mesh_file_name,
88  255, &flg_file);
89  if (flg_file != PETSC_TRUE)
90  SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA,
91  "*** ERROR -my_file (MESH FILE NEEDED)");
92 
93  // Read mesh to MOAB
94  CHKERR moab.load_file(mesh_file_name, 0, "");
95 
96  // Create MoFEM instance
97  MoFEM::Core core(moab);
98  MoFEM::Interface &m_field = core;
99 
100  // set entities bit level
101  BitRefLevel bit_level0 = BitRefLevel().set(0);
102  CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(
103  0, 2, bit_level0);
104 
105  // Declare elements
106  enum bases { AINSWORTH, DEMKOWICZ, LASBASETOP };
107  const char *list_bases[] = {"ainsworth", "demkowicz"};
108  PetscBool flg;
109  PetscInt choice_base_value = AINSWORTH;
110  CHKERR PetscOptionsGetEList(PETSC_NULL, NULL, "-base", list_bases,
111  LASBASETOP, &choice_base_value, &flg);
112  if (flg != PETSC_TRUE)
113  SETERRQ(PETSC_COMM_SELF, MOFEM_IMPOSSIBLE_CASE, "base not set");
115  if (choice_base_value == AINSWORTH)
117  else if (choice_base_value == DEMKOWICZ)
118  base = DEMKOWICZ_JACOBI_BASE;
119  int order = 5;
120  CHKERR PetscOptionsGetInt(PETSC_NULL, "", "-order", &order, PETSC_NULL);
121 
122  CHKERR m_field.add_field("FIELD1", HCURL, base, 1);
123  CHKERR m_field.add_finite_element("FACE_FE");
124  // Define rows/cols and element data
125  CHKERR m_field.modify_finite_element_add_field_row("FACE_FE", "FIELD1");
126  CHKERR m_field.modify_finite_element_add_field_col("FACE_FE", "FIELD1");
127  CHKERR m_field.modify_finite_element_add_field_data("FACE_FE", "FIELD1");
128 
129  CHKERR m_field.add_finite_element("EDGE_FE");
130  // Define rows/cols and element data
131  CHKERR m_field.modify_finite_element_add_field_row("EDGE_FE", "FIELD1");
132  CHKERR m_field.modify_finite_element_add_field_col("EDGE_FE", "FIELD1");
133  CHKERR m_field.modify_finite_element_add_field_data("EDGE_FE", "FIELD1");
134 
135  // Problem
136  CHKERR m_field.add_problem("TEST_PROBLEM");
137  // set finite elements for problem
138  CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM", "FACE_FE");
139  CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM", "EDGE_FE");
140  // set refinement level for problem
141  CHKERR m_field.modify_problem_ref_level_add_bit("TEST_PROBLEM", bit_level0);
142 
143  // Add entities
144 
145  CHKERR m_field.add_ents_to_field_by_type(0, MBTRI, "FIELD1");
146  CHKERR m_field.add_ents_to_field_by_type(0, MBQUAD, "FIELD1");
147  // Set order
148  CHKERR m_field.set_field_order(0, MBTRI, "FIELD1", order);
149  CHKERR m_field.set_field_order(0, MBQUAD, "FIELD1", order);
150  CHKERR m_field.set_field_order(0, MBEDGE, "FIELD1", order);
151 
152  // Add entities to elements
153  CHKERR m_field.add_ents_to_finite_element_by_type(0, MBTRI, "FACE_FE");
154  CHKERR m_field.add_ents_to_finite_element_by_type(0, MBQUAD, "FACE_FE");
155 
156  auto set_edge_elements_entities_on_mesh_skin = [&]() {
158  Range faces;
159  CHKERR moab.get_entities_by_dimension(0, 2, faces, false);
160  Skinner skin(&m_field.get_moab());
161  Range faces_skin;
162  CHKERR skin.find_skin(0, faces, false, faces_skin);
163  CHKERR m_field.add_ents_to_finite_element_by_type(faces_skin, MBEDGE,
164  "EDGE_FE");
166  };
167  CHKERR set_edge_elements_entities_on_mesh_skin();
168 
169  // Build database
170  CHKERR m_field.build_fields();
171  // build finite elemnts
172  CHKERR m_field.build_finite_elements();
173  // build adjacencies
174  CHKERR m_field.build_adjacencies(bit_level0);
175 
176  // build problem
177  ProblemsManager *prb_mng_ptr;
178  CHKERR m_field.getInterface(prb_mng_ptr);
179  CHKERR prb_mng_ptr->buildProblem("TEST_PROBLEM", true);
180  // Partition
181  CHKERR prb_mng_ptr->partitionSimpleProblem("TEST_PROBLEM");
182  CHKERR prb_mng_ptr->partitionFiniteElements("TEST_PROBLEM");
183  CHKERR prb_mng_ptr->partitionGhostDofs("TEST_PROBLEM");
184 
185  // integration rule
186  auto rule = [&](int, int, int p) { return 2 * p; };
187 
188  auto calculate_divergence = [&]() {
189  double div = 0;
190  FaceEle fe_face(m_field);
191  fe_face.getRuleHook = rule;
192  CHKERR AddHOOps<2, 2, 3>::add(fe_face.getOpPtrVector(), {HDIV});
193  fe_face.getOpPtrVector().push_back(new OpDivergence(div));
194  CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "FACE_FE", fe_face);
195  return div;
196  };
197 
198  auto calculate_flux = [&]() {
199  double flux = 0;
200  EdgeEle fe_edge(m_field);
201  fe_edge.getRuleHook = rule;
202  fe_edge.getOpPtrVector().push_back(new OpFlux(flux));
203  CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "EDGE_FE", fe_edge);
204  return flux;
205  };
206 
207  const double div = calculate_divergence();
208  const double flux = calculate_flux();
209 
210  PetscPrintf(PETSC_COMM_WORLD, "Div = %4.3e Flux = %3.4e Error = %4.3e\n",
211  div, flux, div + flux);
212 
213  constexpr double tol = 1e-8;
214  if (std::abs(div + flux) > tol)
215  SETERRQ2(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
216  "Test failed (div != flux) %3.4e != %3.4e", div, flux);
217  }
218  CATCH_ERRORS;
219 
221 }

Variable Documentation

◆ help

char help[] = "...\n\n"
static
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Definition: UnknownInterface.hpp:93
MoFEM::EdgeElementForcesAndSourcesCore
Edge finite element.
Definition: EdgeElementForcesAndSourcesCore.hpp:30
MoFEM::CoreInterface::loop_finite_elements
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.
MoFEM::ProblemsManager::buildProblem
MoFEMErrorCode buildProblem(const std::string name, const bool square_matrix, int verb=VERBOSE)
build problem data structures
Definition: ProblemsManager.cpp:87
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
MoFEM::ProblemsManager
Problem manager is used to build and partition problems.
Definition: ProblemsManager.hpp:21
MoFEM::CoreInterface::modify_finite_element_add_field_row
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
MoFEM::CoreTmp< 0 >::Finalize
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:112
MOFEM_IMPOSSIBLE_CASE
@ MOFEM_IMPOSSIBLE_CASE
Definition: definitions.h:35
MoFEM::CoreInterface::add_ents_to_field_by_type
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.
order
constexpr int order
Definition: dg_projection.cpp:18
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:2010
MoFEM::CoreInterface::add_ents_to_finite_element_by_type
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
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
MoFEM::CoreInterface::add_finite_element
virtual MoFEMErrorCode add_finite_element(const std::string &fe_name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
add finite element
MoFEM::CoreInterface::get_moab
virtual moab::Interface & get_moab()=0
MoFEM::CoreInterface::modify_finite_element_add_field_col
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
MoFEM::CoreInterface::build_finite_elements
virtual MoFEMErrorCode build_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elements.
OpFlux
Definition: hcurl_divergence_operator_2d.cpp:56
MoFEM::FaceElementForcesAndSourcesCore
Face finite element.
Definition: FaceElementForcesAndSourcesCore.hpp:23
MoFEM::CoreInterface::modify_finite_element_add_field_data
virtual MoFEMErrorCode modify_finite_element_add_field_data(const std::string &fe_name, const std::string name_field)=0
set finite element field data
MoFEM::ProblemsManager::partitionFiniteElements
MoFEMErrorCode partitionFiniteElements(const std::string name, bool part_from_moab=false, int low_proc=-1, int hi_proc=-1, int verb=VERBOSE)
partition finite elements
Definition: ProblemsManager.cpp:2167
MoFEM::CoreInterface::modify_problem_ref_level_add_bit
virtual MoFEMErrorCode modify_problem_ref_level_add_bit(const std::string &name_problem, const BitRefLevel &bit)=0
add ref level to problem
mesh_file_name
char mesh_file_name[255]
Definition: mesh_smoothing.cpp:22
MoFEM::AddHOOps
Add operators pushing bases from local to physical configuration.
Definition: HODataOperators.hpp:413
Range
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
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:385
DEMKOWICZ_JACOBI_BASE
@ DEMKOWICZ_JACOBI_BASE
Definition: definitions.h:66
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1148
MoFEM::ProblemsManager::partitionSimpleProblem
MoFEMErrorCode partitionSimpleProblem(const std::string name, int verb=VERBOSE)
partition problem dofs
Definition: ProblemsManager.cpp:1537
AINSWORTH_LEGENDRE_BASE
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
HCURL
@ HCURL
field with continuous tangents
Definition: definitions.h:86
MoFEM::PetscOptionsGetString
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
Definition: DeprecatedPetsc.hpp:172
MoFEM::PetscOptionsGetEList
PetscErrorCode PetscOptionsGetEList(PetscOptions *, const char pre[], const char name[], const char *const *list, PetscInt next, PetscInt *value, PetscBool *set)
Definition: DeprecatedPetsc.hpp:203
MoFEM::CoreInterface::build_fields
virtual MoFEMErrorCode build_fields(int verb=DEFAULT_VERBOSITY)=0
FieldApproximationBase
FieldApproximationBase
approximation base
Definition: definitions.h:58
MoFEM::BitRefManager
Managing BitRefLevels.
Definition: BitRefManager.hpp:21
MoFEM::CoreInterface::modify_problem_add_finite_element
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
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_ATOM_TEST_INVALID
@ MOFEM_ATOM_TEST_INVALID
Definition: definitions.h:40
MoFEM::CoreInterface::build_adjacencies
virtual MoFEMErrorCode build_adjacencies(const Range &ents, int verb=DEFAULT_VERBOSITY)=0
build adjacencies
MoFEM::CoreInterface::set_field_order
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.
MoFEM::ProblemsManager::partitionGhostDofs
MoFEMErrorCode partitionGhostDofs(const std::string name, int verb=VERBOSE)
determine ghost nodes
Definition: ProblemsManager.cpp:2339
MoFEM::CoreInterface::add_problem
virtual MoFEMErrorCode add_problem(const std::string &name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add problem.
convert.int
int
Definition: convert.py:64
MoFEM::PetscOptionsGetInt
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
Definition: DeprecatedPetsc.hpp:142
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
MoFEM::CoreInterface::add_field
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.
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
MOFEM_INVALID_DATA
@ MOFEM_INVALID_DATA
Definition: definitions.h:36
OpDivergence
Definition: hcurl_divergence_operator_2d.cpp:29
help
static char help[]
Definition: hcurl_divergence_operator_2d_embaded_in_3d.cpp:13
tol
double tol
Definition: mesh_smoothing.cpp:26