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
10
11#include <MoFEM.hpp>
12
13using namespace MoFEM;
14
15static char help[] = "...\n\n";
16
18
20
23
24struct OpDivergence : public FaceEleOp {
25
26 double &dIv;
27 OpDivergence(double &div) : FaceEleOp("FIELD1", OPROW), dIv(div) {}
28
32 const int nb_dofs = data.getIndices().size();
33 if (nb_dofs == 0)
35 const int nb_gauss_pts = data.getN().size1();
36 auto t_diff_base_fun = data.getFTensor2DiffN<3, 3>();
37 for (int gg = 0; gg != nb_gauss_pts; gg++) {
38 const double val = getArea() * getGaussPts()(2, gg);
39 for (int bb = 0; bb != nb_dofs; bb++) {
40 dIv += val * (t_diff_base_fun(0, 0) + t_diff_base_fun(1, 1) +
41 t_diff_base_fun(2, 2));
42 ++t_diff_base_fun;
43 }
44 }
46 }
47};
48
49struct OpFlux : public EdgeEleOp {
50
51 double &fLux;
52 OpFlux(double &flux) : EdgeEleOp("FIELD1", OPROW), fLux(flux) {}
53
55
59 const int nb_dofs = data.getIndices().size();
60 if (nb_dofs == 0)
62 const int nb_gauss_pts = data.getN().size1();
63
64 auto t_base_fun = data.getFTensor1N<3>();
66 for (int gg = 0; gg != nb_gauss_pts; gg++) {
67 const double val = getGaussPts()(1, gg);
68 for (int bb = 0; bb != nb_dofs; bb++) {
69 fLux += val * t_base_fun(0);
70 ++t_base_fun;
71 }
72 }
73
75 }
76};
77
78int main(int argc, char *argv[]) {
79
80 MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
81
82 try {
83
84 moab::Core mb_instance;
85 moab::Interface &moab = mb_instance;
86
87 PetscBool flg_file = PETSC_TRUE;
88 char mesh_file_name[255];
89 CHKERR PetscOptionsGetString(PETSC_NULL, "", "-my_file", mesh_file_name,
90 255, &flg_file);
91 if (flg_file != PETSC_TRUE)
92 SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA,
93 "*** ERROR -my_file (MESH FILE NEEDED)");
94
95 // Read mesh to MOAB
96 CHKERR moab.load_file(mesh_file_name, 0, "");
97
98 // Create MoFEM instance
99 MoFEM::Core core(moab);
100 MoFEM::Interface &m_field = core;
101
102 // set entities bit level
103 BitRefLevel bit_level0 = BitRefLevel().set(0);
104 CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(
105 0, 2, bit_level0);
106
107 // Declare elements
108 enum bases { AINSWORTH, DEMKOWICZ, LASBASETOP };
109 const char *list_bases[] = {"ainsworth", "demkowicz"};
110 PetscBool flg;
111 PetscInt choice_base_value = AINSWORTH;
112 CHKERR PetscOptionsGetEList(PETSC_NULL, NULL, "-base", list_bases,
113 LASBASETOP, &choice_base_value, &flg);
114 if (flg != PETSC_TRUE)
115 SETERRQ(PETSC_COMM_SELF, MOFEM_IMPOSSIBLE_CASE, "base not set");
117 if (choice_base_value == AINSWORTH)
119 else if (choice_base_value == DEMKOWICZ)
121 int order = 5;
122 CHKERR PetscOptionsGetInt(PETSC_NULL, "", "-order", &order, PETSC_NULL);
123
124 CHKERR m_field.add_field("FIELD1", HCURL, base, 1);
125 CHKERR m_field.add_finite_element("FACE_FE");
126 // Define rows/cols and element data
127 CHKERR m_field.modify_finite_element_add_field_row("FACE_FE", "FIELD1");
128 CHKERR m_field.modify_finite_element_add_field_col("FACE_FE", "FIELD1");
129 CHKERR m_field.modify_finite_element_add_field_data("FACE_FE", "FIELD1");
130
131 CHKERR m_field.add_finite_element("EDGE_FE");
132 // Define rows/cols and element data
133 CHKERR m_field.modify_finite_element_add_field_row("EDGE_FE", "FIELD1");
134 CHKERR m_field.modify_finite_element_add_field_col("EDGE_FE", "FIELD1");
135 CHKERR m_field.modify_finite_element_add_field_data("EDGE_FE", "FIELD1");
136
137 // Problem
138 CHKERR m_field.add_problem("TEST_PROBLEM");
139 // set finite elements for problem
140 CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM", "FACE_FE");
141 CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM", "EDGE_FE");
142 // set refinement level for problem
143 CHKERR m_field.modify_problem_ref_level_add_bit("TEST_PROBLEM", bit_level0);
144
145 // Add entities
146
147 CHKERR m_field.add_ents_to_field_by_type(0, MBTRI, "FIELD1");
148 CHKERR m_field.add_ents_to_field_by_type(0, MBQUAD, "FIELD1");
149 // Set order
150 CHKERR m_field.set_field_order(0, MBTRI, "FIELD1", order);
151 CHKERR m_field.set_field_order(0, MBQUAD, "FIELD1", order);
152 CHKERR m_field.set_field_order(0, MBEDGE, "FIELD1", order);
153
154 // Add entities to elements
155 CHKERR m_field.add_ents_to_finite_element_by_type(0, MBTRI, "FACE_FE");
156 CHKERR m_field.add_ents_to_finite_element_by_type(0, MBQUAD, "FACE_FE");
157
158 auto set_edge_elements_entities_on_mesh_skin = [&]() {
160 Range faces;
161 CHKERR moab.get_entities_by_dimension(0, 2, faces, false);
162 Skinner skin(&m_field.get_moab());
163 Range faces_skin;
164 CHKERR skin.find_skin(0, faces, false, faces_skin);
165 CHKERR m_field.add_ents_to_finite_element_by_type(faces_skin, MBEDGE,
166 "EDGE_FE");
168 };
169 CHKERR set_edge_elements_entities_on_mesh_skin();
170
171 // Build database
172 CHKERR m_field.build_fields();
173 // build finite elemnts
175 // build adjacencies
176 CHKERR m_field.build_adjacencies(bit_level0);
177
178 // build problem
179 ProblemsManager *prb_mng_ptr;
180 CHKERR m_field.getInterface(prb_mng_ptr);
181 CHKERR prb_mng_ptr->buildProblem("TEST_PROBLEM", true);
182 // Partition
183 CHKERR prb_mng_ptr->partitionSimpleProblem("TEST_PROBLEM");
184 CHKERR prb_mng_ptr->partitionFiniteElements("TEST_PROBLEM");
185 CHKERR prb_mng_ptr->partitionGhostDofs("TEST_PROBLEM");
186
187 // integration rule
188 auto rule = [&](int, int, int p) { return 2 * p; };
189
190 auto calculate_divergence = [&]() {
191 double div = 0;
192 FaceEle fe_face(m_field);
193 fe_face.getRuleHook = rule;
194
195 auto jac_ptr = boost::make_shared<MatrixDouble>();
196 auto inv_jac_ptr = boost::make_shared<MatrixDouble>();
197 auto det_ptr = boost::make_shared<VectorDouble>();
198
199 fe_face.getOpPtrVector().push_back(
201 fe_face.getOpPtrVector().push_back(new OpMakeHdivFromHcurl());
202 fe_face.getOpPtrVector().push_back(
204 jac_ptr));
205
206 fe_face.getOpPtrVector().push_back(
207 new OpInvertMatrix<3>(jac_ptr, det_ptr, inv_jac_ptr));
208
209 fe_face.getOpPtrVector().push_back(
211 fe_face.getOpPtrVector().push_back(new OpSetHOWeightsOnFace());
212 fe_face.getOpPtrVector().push_back(new OpDivergence(div));
213 CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "FACE_FE", fe_face);
214 return div;
215 };
216
217 auto calculate_flux = [&]() {
218 double flux = 0;
219 EdgeEle fe_edge(m_field);
220 fe_edge.getRuleHook = rule;
221 fe_edge.getOpPtrVector().push_back(new OpFlux(flux));
222 CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "EDGE_FE", fe_edge);
223 return flux;
224 };
225
226 const double div = calculate_divergence();
227 const double flux = calculate_flux();
228
229 PetscPrintf(PETSC_COMM_WORLD, "Div = %4.3e Flux = %3.4e Error = %4.3e\n",
230 div, flux, div + flux);
231
232 constexpr double tol = 1e-8;
233 if (std::abs(div + flux) > tol)
234 SETERRQ2(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
235 "Test failed (div != flux) %3.4e != %3.4e", div, flux);
236 }
238
240}
static Index< 'p', 3 > p
int main()
Definition: adol-c_atom.cpp:46
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:372
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()
Definition: definitions.h:447
@ 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 ...
Definition: definitions.h:346
@ 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()
Definition: definitions.h:416
#define CHKERR
Inline error check.
Definition: definitions.h:535
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.
Definition: Exceptions.hpp:56
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
OpSetContravariantPiolaTransformOnFace2DImpl< 3 > OpSetContravariantPiolaTransformOnFace2DEmbeddedIn3DSpace
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
OpSetInvJacHcurlFaceImpl< 3 > OpSetInvJacHcurlFaceEmbeddedIn3DSpace
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)
OpCalculateHOJacForFaceImpl< 3 > OpCalculateHOJacForFaceEmbeddedIn3DSpace
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.
Make Hdiv space from Hcurl space in 2d.
Modify integration weights on face to take in account higher-order geometry.
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.