v0.14.0
Loading...
Searching...
No Matches
hcurl_divergence_operator_2d_embaded_in_3d.cpp
Go to the documentation of this file.
1/**
2 * \file hcurl_divergence_operator_2d.cpp
3 * \example hcurl_divergence_operator_2d.cpp
4 *
5 * Testing Hcurl base, transfromed to Hdiv base in 2d using Green theorem
6 *
7 */
8
9#include <MoFEM.hpp>
10
11using namespace MoFEM;
12
13static char help[] = "...\n\n";
14
16
18
21
22struct OpDivergence : public FaceEleOp {
23
24 double &dIv;
25 OpDivergence(double &div) : FaceEleOp("FIELD1", OPROW), dIv(div) {}
26
30 const int nb_dofs = data.getIndices().size();
31 if (nb_dofs == 0)
33 const int nb_gauss_pts = data.getN().size1();
34 auto t_diff_base_fun = data.getFTensor2DiffN<3, 3>();
35 for (int gg = 0; gg != nb_gauss_pts; gg++) {
36 const double val = getArea() * getGaussPts()(2, gg);
37 for (int bb = 0; bb != nb_dofs; bb++) {
38 dIv += val * (t_diff_base_fun(0, 0) + t_diff_base_fun(1, 1) +
39 t_diff_base_fun(2, 2));
40 ++t_diff_base_fun;
41 }
42 }
44 }
45};
46
47struct OpFlux : public EdgeEleOp {
48
49 double &fLux;
50 OpFlux(double &flux) : EdgeEleOp("FIELD1", OPROW), fLux(flux) {}
51
53
57 const int nb_dofs = data.getIndices().size();
58 if (nb_dofs == 0)
60 const int nb_gauss_pts = data.getN().size1();
61
62 auto t_base_fun = data.getFTensor1N<3>();
63 FTensor::Index<'i', 2> i;
64 for (int gg = 0; gg != nb_gauss_pts; gg++) {
65 const double val = getGaussPts()(1, gg);
66 for (int bb = 0; bb != nb_dofs; bb++) {
67 fLux += val * t_base_fun(0);
68 ++t_base_fun;
69 }
70 }
71
73 }
74};
75
76int main(int argc, char *argv[]) {
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)
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
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;
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 }
219
221}
static Index< 'p', 3 > p
int main()
#define CATCH_ERRORS
Catch errors.
FieldApproximationBase
approximation base
Definition definitions.h:58
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition definitions.h:60
@ DEMKOWICZ_JACOBI_BASE
Definition definitions.h:66
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
@ HCURL
field with continuous tangents
Definition definitions.h:86
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ MOFEM_IMPOSSIBLE_CASE
Definition definitions.h:35
@ MOFEM_ATOM_TEST_INVALID
Definition definitions.h:40
@ MOFEM_INVALID_DATA
Definition definitions.h:36
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
constexpr int order
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_data(const std::string &fe_name, const std::string name_filed)=0
set finite element field data
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 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.
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.
MoFEMErrorCode partitionGhostDofs(const std::string name, int verb=VERBOSE)
determine ghost nodes
MoFEMErrorCode partitionSimpleProblem(const std::string name, int verb=VERBOSE)
partition problem dofs
MoFEMErrorCode buildProblem(const std::string name, const bool square_matrix, int verb=VERBOSE)
build problem data structures
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
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
double tol
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
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
PetscErrorCode PetscOptionsGetEList(PetscOptions *, const char pre[], const char name[], const char *const *list, PetscInt next, PetscInt *value, PetscBool *set)
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
Add operators pushing bases from local to physical configuration.
Managing BitRefLevels.
virtual moab::Interface & get_moab()=0
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.
Data on single entity (This is passed as argument to DataOperator::doWork)
FTensor::Tensor2< FTensor::PackPtr< double *, Tensor_Dim0 *Tensor_Dim1 >, Tensor_Dim0, Tensor_Dim1 > getFTensor2DiffN(FieldApproximationBase base)
Get derivatives of base functions for Hdiv space.
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1N(FieldApproximationBase base)
Get base functions for Hdiv/Hcurl spaces.
const VectorInt & getIndices() const
Get global indices of dofs on entity.
@ OPROW
operator doWork function is executed on FE rows
MatrixDouble & getGaussPts()
matrix of integration (Gauss) points for Volume Element
boost::ptr_deque< UserDataOperator > & getOpPtrVector()
Use to push back operator for row operator.
RuleHookFun getRuleHook
Hook to get rule.
Problem manager is used to build and partition problems.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.