v0.13.2
Loading...
Searching...
No Matches
continuity_check_on_skeleton_3d.cpp
Go to the documentation of this file.
1/** \file continuity_check_on_skeleton_3d.cpp
2 * \example continuity_check_on_skeleton_3d.cpp
3 *
4 * \brief Testing integration on skeleton for 3D
5 *
6 * Checking continuity of hdiv and hcurl spaces on faces, and testing methods
7 * for integration on the skeleton.
8 *
9 */
10
11#include <MoFEM.hpp>
12using namespace MoFEM;
13
14static char help[] = "...\n\n";
15
16int main(int argc, char *argv[]) {
17
18 MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
19
20 auto core_log = logging::core::get();
21 core_log->add_sink(
23 LogManager::setLog("ATOM_TEST");
24
25 try {
26
27 moab::Core mb_instance;
28 moab::Interface &moab = mb_instance;
29 int rank;
30 MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
31
32 PetscBool flg = PETSC_TRUE;
33 char mesh_file_name[255];
34#if PETSC_VERSION_GE(3, 6, 4)
35 CHKERR PetscOptionsGetString(PETSC_NULL, "", "-my_file", mesh_file_name,
36 255, &flg);
37#else
38 CHKERR PetscOptionsGetString(PETSC_NULL, PETSC_NULL, "-my_file",
39 mesh_file_name, 255, &flg);
40#endif
41 if (flg != PETSC_TRUE) {
42 SETERRQ(PETSC_COMM_SELF, 1, "*** ERROR -my_file (MESH FILE NEEDED)");
43 }
44
45 const char *option;
46 option = "";
47 CHKERR moab.load_file(mesh_file_name, 0, option);
48
49 // Create MoFEM database
50 MoFEM::Core core(moab);
51 // Access to database through interface
52 MoFEM::Interface &m_field = core;
53
54 // set entireties bit level
55 BitRefLevel bit_level0 = BitRefLevel().set(0);
56 CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(
57 0, 3, bit_level0);
58
59 // Fields
60 auto get_base = []() -> FieldApproximationBase {
61 enum bases { AINSWORTH, DEMKOWICZ, LASTBASEOP };
62 const char *list_bases[] = {"ainsworth", "demkowicz"};
63 PetscBool flg;
64 PetscInt choice_base_value = AINSWORTH;
65 CHKERR PetscOptionsGetEList(PETSC_NULL, NULL, "-base", list_bases,
66 LASTBASEOP, &choice_base_value, &flg);
67 if (flg == PETSC_TRUE) {
69 if (choice_base_value == AINSWORTH)
71 else if (choice_base_value == DEMKOWICZ)
73 return base;
74 }
75 return LASTBASE;
76 };
77
78 auto get_space = []() -> FieldSpace {
79 enum spaces { HDIV, HCURL, LASTSPACEOP };
80 const char *list_spaces[] = {"hdiv", "hcurl"};
81 PetscBool flg;
82 PetscInt choice_space_value = HDIV;
83 CHKERR PetscOptionsGetEList(PETSC_NULL, NULL, "-space", list_spaces,
84 LASTSPACEOP, &choice_space_value, &flg);
85 if (flg == PETSC_TRUE) {
87 if (choice_space_value == HDIV)
88 space = FieldSpace::HDIV;
89 else if (choice_space_value == HCURL)
90 space = FieldSpace::HCURL;
91 return space;
92 }
93 return LASTSPACE;
94 };
95
96 CHKERR m_field.add_field("F2", get_space(), get_base(), 1);
97
98 // meshset consisting all entities in mesh
99 EntityHandle root_set = moab.get_root_set();
100 // add entities to field
101 CHKERR m_field.add_ents_to_field_by_dim(root_set, 3, "F2");
102
103 // set app. order
104 int order = 3;
105 CHKERR m_field.set_field_order(root_set, MBTET, "F2", order);
106 CHKERR m_field.set_field_order(root_set, MBHEX, "F2", order);
107 CHKERR m_field.set_field_order(root_set, MBTRI, "F2", order);
108 CHKERR m_field.set_field_order(root_set, MBQUAD, "F2", order);
109 CHKERR m_field.set_field_order(root_set, MBEDGE, "F2", order);
110
111 CHKERR m_field.build_fields();
112
113 // add elements
114 CHKERR m_field.add_finite_element("V1");
115 CHKERR m_field.add_finite_element("S2");
122
123 CHKERR m_field.add_ents_to_finite_element_by_dim(root_set, 3, "V1");
124 Range faces;
125 CHKERR m_field.getInterface<BitRefManager>()->getEntitiesByDimAndRefLevel(
126 bit_level0, BitRefLevel().set(), 2, faces);
127 CHKERR m_field.add_ents_to_finite_element_by_dim(faces, 2, "S2");
128
130 CHKERR m_field.build_adjacencies(bit_level0);
131
132 // Problems
133 CHKERR m_field.add_problem("P1");
134
135 // set refinement level for problem
136 CHKERR m_field.modify_problem_ref_level_add_bit("P1", bit_level0);
137 CHKERR m_field.modify_problem_add_finite_element("P1", "V1");
138 CHKERR m_field.modify_problem_add_finite_element("P1", "S2");
139
140 // build problems
141 ProblemsManager *prb_mng_ptr;
142 CHKERR m_field.getInterface(prb_mng_ptr);
143 CHKERR prb_mng_ptr->buildProblem("P1", true);
144 CHKERR prb_mng_ptr->partitionProblem("P1");
145 CHKERR prb_mng_ptr->partitionFiniteElements("P1");
146 CHKERR prb_mng_ptr->partitionGhostDofs("P1");
147
148 struct CommonData {
149 VectorDouble dotNormalFace;
150 VectorDouble dotNormalEleLeft;
151 VectorDouble dotNormalEleRight;
152 };
153
154 struct SkeletonFE
156
158 struct OpVolSide
160
161 CommonData &elemData;
162 OpVolSide(CommonData &elem_data)
164 "F2", UserDataOperator::OPROW),
165 elemData(elem_data) {}
166 MoFEMErrorCode doWork(int side, EntityType type,
169
170 if (CN::Dimension(type) == 3) {
171 MatrixDouble diff =
172 getCoordsAtGaussPts() - getFaceCoordsAtGaussPts();
173
174 MOFEM_LOG("ATOM_TEST", Sev::noisy)
175 << "getCoordsAtGaussPts() " << getCoordsAtGaussPts();
176 MOFEM_LOG("ATOM_TEST", Sev::noisy)
177 << "getFaceCoordsAtGaussPts() " << getFaceCoordsAtGaussPts();
178
179 constexpr double eps = 1e-12;
180 if (norm_inf(diff) > eps) {
181 MOFEM_LOG("ATOM_TEST", Sev::error) << "diff " << diff;
182 SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
183 "Coordinates at integration pts are different");
184 }
185 }
186
187 const size_t nb_dofs = data.getFieldData().size();
188 if (nb_dofs) {
189 const size_t nb_integration_pts = getGaussPts().size2();
190
192 auto t_to_do_dot = getFTensor1Normal();
193 if (data.getSpace() == HCURL) {
194 auto s1 = getFTensor1Tangent1();
195 auto s2 = getFTensor1Tangent2();
196 t_to_do_dot(i) = s1(i) + s2(i);
197 }
198
199 VectorDouble *ptr_dot_elem_data = nullptr;
200 if (getFEMethod()->nInTheLoop == 0)
201 ptr_dot_elem_data = &elemData.dotNormalEleLeft;
202 else
203 ptr_dot_elem_data = &elemData.dotNormalEleRight;
204 auto &dot_elem_data = *ptr_dot_elem_data;
205 if (dot_elem_data.size() == 0) {
206 dot_elem_data.resize(nb_integration_pts, false);
207 dot_elem_data.clear();
208 }
209
210 auto t_hdiv_base = data.getFTensor1N<3>();
211 for (size_t gg = 0; gg != nb_integration_pts; ++gg) {
212 auto t_data = data.getFTensor0FieldData();
213 for (size_t bb = 0; bb != nb_dofs; ++bb) {
214 dot_elem_data(gg) += (t_to_do_dot(i) * t_hdiv_base(i)) * t_data;
215 ++t_hdiv_base;
216 ++t_data;
217 }
218 }
219 }
221 }
222 };
223
225 SkeletonFE(MoFEM::Interface &m_field, CommonData &elem_data)
227 "F2", UserDataOperator::OPROW),
228 volSideFe(m_field), elemData(elem_data) {
229
230 auto jac_ptr = boost::make_shared<MatrixDouble>();
231 auto inv_jac_ptr = boost::make_shared<MatrixDouble>();
232 auto det_ptr = boost::make_shared<VectorDouble>();
233
234 volSideFe.getOpPtrVector().push_back(new OpCalculateHOJac<3>(jac_ptr));
235 volSideFe.getOpPtrVector().push_back(
236 new OpInvertMatrix<3>(jac_ptr, det_ptr, inv_jac_ptr));
237 volSideFe.getOpPtrVector().push_back(new OpSetHOWeights(det_ptr));
238 volSideFe.getOpPtrVector().push_back(
239 new OpSetHOCovariantPiolaTransform(HCURL, inv_jac_ptr));
240 volSideFe.getOpPtrVector().push_back(
241 new OpSetHOContravariantPiolaTransform(HDIV, det_ptr, jac_ptr));
242 volSideFe.getOpPtrVector().push_back(
243 new SkeletonFE::OpVolSide(elemData));
244 }
245
246 MoFEMErrorCode doWork(int side, EntityType type,
248
250
251 const size_t nb_integration_pts = getGaussPts().size2();
252
253 if (side == 0 && type == MBEDGE) {
254 elemData.dotNormalFace.resize(nb_integration_pts, false);
255 elemData.dotNormalFace.clear();
256 }
257
258 const size_t nb_dofs = data.getFieldData().size();
259 if (nb_dofs) {
261 auto t_to_do_dot = getFTensor1Normal();
262 if (data.getSpace() == HCURL) {
263 auto s1 = getFTensor1Tangent1();
264 auto s2 = getFTensor1Tangent2();
265 t_to_do_dot(i) = s1(i) + s2(i);
266 }
267 auto t_hdiv_base = data.getFTensor1N<3>();
268 for (size_t gg = 0; gg != nb_integration_pts; ++gg) {
269 auto t_data = data.getFTensor0FieldData();
270 for (size_t bb = 0; bb != nb_dofs; ++bb) {
271 elemData.dotNormalFace(gg) +=
272 (t_to_do_dot(i) * t_hdiv_base(i)) * t_data;
273 ++t_hdiv_base;
274 ++t_data;
275 }
276 }
277 }
278
279 if (CN::Dimension(type) == 2) {
280
281 elemData.dotNormalEleLeft.resize(0, false);
282 elemData.dotNormalEleRight.resize(0, false);
283 CHKERR loopSideVolumes("V1", volSideFe);
284
285 auto check_continuity_of_base = [&](auto &vol_dot_data) {
287 if (vol_dot_data.size() != elemData.dotNormalFace.size())
288 SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
289 "Inconsistent number of integration points");
290 const double eps = 1e-12;
291 for (size_t gg = 0; gg != nb_integration_pts; ++gg) {
292 const double error =
293 std::abs(vol_dot_data(gg) - elemData.dotNormalFace(gg));
294 MOFEM_LOG("ATOM_TEST", Sev::noisy) << "Error: " << error;
295 if (error > eps)
296 SETERRQ2(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
297 "Inconsistency %3.4e != %3.4e", vol_dot_data(gg),
298 elemData.dotNormalFace(gg));
299 }
301 };
302 if (elemData.dotNormalEleLeft.size() != 0)
303 CHKERR check_continuity_of_base(elemData.dotNormalEleLeft);
304 else if (elemData.dotNormalEleRight.size() != 0)
305 CHKERR check_continuity_of_base(elemData.dotNormalEleRight);
306 }
307
309 }
310 };
311
312 CommonData elem_data;
313 FaceElementForcesAndSourcesCore face_fe(m_field);
314 face_fe.getRuleHook = [](int, int, int) { return 1; };
315 face_fe.getOpPtrVector().push_back(
317 face_fe.getOpPtrVector().push_back(
319 face_fe.getOpPtrVector().push_back(new SkeletonFE(m_field, elem_data));
320
321 auto field_ents_ptr = m_field.get_field_ents();
322
323 auto cache_ptr = boost::make_shared<CacheTuple>();
324 CHKERR m_field.cache_problem_entities("P1", cache_ptr);
325
326 for (auto &ent_ptr : (*field_ents_ptr)) {
327 MOFEM_LOG("ATOM_TEST", Sev::verbose) << *ent_ptr;
328 for(auto &v : ent_ptr->getEntFieldData())
329 v = 1;
330 CHKERR m_field.loop_finite_elements("P1", "S2", face_fe, 0, 1, nullptr,
331 MF_ZERO, cache_ptr);
332 for (auto &v : ent_ptr->getEntFieldData())
333 v = 0;
334 }
335 }
337
339
340 return 0;
341}
int main()
Definition: adol-c_atom.cpp:46
static const double eps
static char help[]
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:372
@ MF_ZERO
Definition: definitions.h:98
FieldApproximationBase
approximation base
Definition: definitions.h:58
@ LASTBASE
Definition: definitions.h:69
@ 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
FieldSpace
approximation spaces
Definition: definitions.h:82
@ LASTSPACE
FieldSpace in [ 0, LASTSPACE )
Definition: definitions.h:89
@ HCURL
field with continuous tangents
Definition: definitions.h:86
@ HDIV
field with continuous normal traction
Definition: definitions.h:87
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
@ MOFEM_ATOM_TEST_INVALID
Definition: definitions.h:40
#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
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:440
virtual const FieldEntity_multiIndex * get_field_ents() const =0
Get the field ents object.
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.
static LoggerType & setLog(const std::string channel)
Set ans resset chanel logger.
Definition: LogManager.cpp:391
#define MOFEM_LOG(channel, severity)
Log.
Definition: LogManager.hpp:308
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 buildProblem(const std::string name, const bool square_matrix, int verb=VERBOSE)
build problem data structures
MoFEMErrorCode partitionProblem(const std::string name, int verb=VERBOSE)
partition problem dofs (collective)
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
FTensor::Index< 'i', SPACE_DIM > i
const double v
phase velocity of light in medium (cm/ns)
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
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)
Managing BitRefLevels.
virtual MoFEMErrorCode cache_problem_entities(const std::string prb_name, CacheTupleWeakPtr cache_ptr)=0
Cache variables.
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)
const VectorDouble & getFieldData() const
get dofs values
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0FieldData()
Resturn scalar files as a FTensor of rank 0.
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1N(FieldApproximationBase base)
Get base functions for Hdiv/Hcurl spaces.
FieldSpace & getSpace()
Get field space.
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.
static boost::shared_ptr< SinkType > createSink(boost::shared_ptr< std::ostream > stream_ptr, std::string comm_filter)
Create a sink object.
Definition: LogManager.cpp:300
static boost::shared_ptr< std::ostream > getStrmSelf()
Get the strm self object.
Definition: LogManager.cpp:342
transform Hdiv base fluxes from reference element to physical triangle
transform Hcurl base fluxes from reference element to physical triangle
Apply contravariant (Piola) transfer to Hdiv space for HO geometry.
Apply covariant (Piola) transfer to Hcurl space for HO geometry.
Set inverse jacobian to base functions.
Problem manager is used to build and partition problems.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.
Base volume element used to integrate on skeleton.
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.