v0.14.0
Loading...
Searching...
No Matches
hcurl_divergence_operator_2d.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 * Note 0: This is low-level implementation.
8 * Note 1: Generic implementation for Quad/Tri mesh of arbitrary order.
9 *
10 */
11
12
13
14#include <MoFEM.hpp>
15
16using namespace MoFEM;
17
18static char help[] = "...\n\n";
19
22
25
26constexpr int SPACE_DIM = 2;
28
29struct OpDivergence : public FaceEleOp {
30
31 double &dIv;
32 OpDivergence(double &div) : FaceEleOp("FIELD1", OPROW), dIv(div) {}
33
37 const int nb_dofs = data.getIndices().size();
38 if (nb_dofs == 0)
40 const int nb_gauss_pts = data.getN().size1();
42 auto t_diff_base_fun = data.getFTensor2DiffN<3, 2>();
43 const auto area = getMeasure();
44 for (int gg = 0; gg != nb_gauss_pts; gg++) {
45 const double val = area * t_w;
46 for (int bb = 0; bb != nb_dofs; bb++) {
47 dIv += val * t_diff_base_fun(i, i);
48 ++t_diff_base_fun;
49 }
50 ++t_w;
51 }
53 }
54};
55
56struct OpFlux : public EdgeEleOp {
57
58 double &fLux;
59 OpFlux(double &flux) : EdgeEleOp("FIELD1", OPROW), fLux(flux) {}
60
64 const int nb_dofs = data.getIndices().size();
65 if (nb_dofs == 0)
67 const int nb_gauss_pts = data.getN().size1();
68 auto t_normal = getFTensor1Normal();
69 auto t_base_fun = data.getFTensor1N<3>();
71 for (int gg = 0; gg != nb_gauss_pts; gg++) {
72 for (int bb = 0; bb != nb_dofs; bb++) {
73 fLux += t_w * t_normal(i) * t_base_fun(i);
74 ++t_base_fun;
75 }
76 ++t_w;
77 }
78
80 }
81};
82
83int main(int argc, char *argv[]) {
84
85 MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
86
87 try {
88
89 moab::Core mb_instance;
90 moab::Interface &moab = mb_instance;
91
92 PetscBool flg_file = PETSC_TRUE;
93 char mesh_file_name[255];
94 CHKERR PetscOptionsGetString(PETSC_NULL, "", "-my_file", mesh_file_name,
95 255, &flg_file);
96 if (flg_file != PETSC_TRUE)
97 SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA,
98 "*** ERROR -my_file (MESH FILE NEEDED)");
99
100 // Read mesh to MOAB
101 CHKERR moab.load_file(mesh_file_name, 0, "");
102
103 // Create MoFEM instance
104 MoFEM::Core core(moab);
105 MoFEM::Interface &m_field = core;
106
107 // set entities bit level
108 BitRefLevel bit_level0 = BitRefLevel().set(0);
109 CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(
110 0, 2, bit_level0);
111
112 // Declare elements
113 enum bases { AINSWORTH, DEMKOWICZ, LASBASETOP };
114 const char *list_bases[] = {"ainsworth", "demkowicz"};
115 PetscBool flg;
116 PetscInt choice_base_value = AINSWORTH;
117 CHKERR PetscOptionsGetEList(PETSC_NULL, NULL, "-base", list_bases,
118 LASBASETOP, &choice_base_value, &flg);
119 if (flg != PETSC_TRUE)
120 SETERRQ(PETSC_COMM_SELF, MOFEM_IMPOSSIBLE_CASE, "base not set");
122 if (choice_base_value == AINSWORTH)
124 else if (choice_base_value == DEMKOWICZ)
126 int order = 5;
127 CHKERR PetscOptionsGetInt(PETSC_NULL, "", "-order", &order, PETSC_NULL);
128
129 CHKERR m_field.add_field("FIELD1", HCURL, base, 1);
130 CHKERR m_field.add_finite_element("FACE_FE");
131 // Define rows/cols and element data
132 CHKERR m_field.modify_finite_element_add_field_row("FACE_FE", "FIELD1");
133 CHKERR m_field.modify_finite_element_add_field_col("FACE_FE", "FIELD1");
134 CHKERR m_field.modify_finite_element_add_field_data("FACE_FE", "FIELD1");
135
136 CHKERR m_field.add_finite_element("EDGE_FE");
137 // Define rows/cols and element data
138 CHKERR m_field.modify_finite_element_add_field_row("EDGE_FE", "FIELD1");
139 CHKERR m_field.modify_finite_element_add_field_col("EDGE_FE", "FIELD1");
140 CHKERR m_field.modify_finite_element_add_field_data("EDGE_FE", "FIELD1");
141
142 // Problem
143 CHKERR m_field.add_problem("TEST_PROBLEM");
144 // set finite elements for problem
145 CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM", "FACE_FE");
146 CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM", "EDGE_FE");
147 // set refinement level for problem
148 CHKERR m_field.modify_problem_ref_level_add_bit("TEST_PROBLEM", bit_level0);
149
150 // Add entities
151
152 CHKERR m_field.add_ents_to_field_by_dim(0, SPACE_DIM, "FIELD1");
153 // Set order
154 for (auto t : {MBEDGE, MBTRI, MBQUAD})
155 CHKERR m_field.set_field_order(0, t, "FIELD1", order);
156
157 // Add entities to elements
159 auto set_edge_elements_entities_on_mesh_skin = [&]() {
161 Range faces;
162 CHKERR moab.get_entities_by_dimension(0, 2, faces, false);
163 Skinner skin(&m_field.get_moab());
164 Range faces_skin;
165 CHKERR skin.find_skin(0, faces, false, faces_skin);
167 faces_skin, SPACE_DIM - 1, "EDGE_FE");
169 };
170 CHKERR set_edge_elements_entities_on_mesh_skin();
171
172 // Build database
173 CHKERR m_field.build_fields();
174 // build finite elemnts
176 // build adjacencies
177 CHKERR m_field.build_adjacencies(bit_level0);
178
179 // build problem
180 ProblemsManager *prb_mng_ptr;
181 CHKERR m_field.getInterface(prb_mng_ptr);
182 CHKERR prb_mng_ptr->buildProblem("TEST_PROBLEM", true);
183 // Partition
184 CHKERR prb_mng_ptr->partitionSimpleProblem("TEST_PROBLEM");
185 CHKERR prb_mng_ptr->partitionFiniteElements("TEST_PROBLEM");
186 CHKERR prb_mng_ptr->partitionGhostDofs("TEST_PROBLEM");
187
188 // integration rule
189 auto rule = [&](int, int, int p) { return 2 * p; };
190
191 auto calculate_divergence = [&]() {
192 double div = 0;
193 FaceEle fe_face(m_field);
194 fe_face.getRuleHook = rule;
196 fe_face.getOpPtrVector().push_back(new OpDivergence(div));
197 CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "FACE_FE", fe_face);
198 return div;
199 };
200
201 auto calculate_flux = [&]() {
202 double flux = 0;
203 EdgeEle fe_edge(m_field);
204 fe_edge.getRuleHook = rule;
206 fe_edge.getOpPtrVector().push_back(new OpFlux(flux));
207 CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "EDGE_FE", fe_edge);
208 return flux;
209 };
210
211 const double div = calculate_divergence();
212 const double flux = calculate_flux();
213
214 MOFEM_LOG_CHANNEL("WOLD"); // reset channel
215 MOFEM_LOG_C("WORLD", Sev::inform,
216 "Div = %4.3e Flux = %3.4e Error = %4.3e\n", div, flux,
217 div - flux);
218
219 constexpr double tol = 1e-8;
220 if (std::abs(div - flux) > tol)
221 SETERRQ2(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
222 "Test failed (div != flux) %3.4e != %3.4e", div, flux);
223 }
225
227}
static Index< 'p', 3 > p
#define MOFEM_LOG_C(channel, severity, format,...)
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_ents_to_finite_element_by_dim(const EntityHandle entities, const int dim, const std::string &name, const bool recursive=true)=0
add entities to 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
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 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 add_ents_to_field_by_dim(const Range &ents, const int dim, const std::string &name, int verb=DEFAULT_VERBOSITY)=0
Add entities to field meshset.
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.
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
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
static char help[]
FTensor::Index< 'i', SPACE_DIM > i
constexpr int SPACE_DIM
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)
constexpr double t
plate stiffness
Definition plate.cpp:59
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.
auto getFTensor1Normal()
get edge normal NOTE: it should be used only in 2D analysis
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.
auto getFTensor0IntegrationWeight()
Get integration weights.
double getMeasure() const
get measure of element
@ OPROW
operator doWork function is executed on FE rows
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.