v0.10.0
Functions | Variables
continuity_check_on_skeleton_3d.cpp File Reference
#include <MoFEM.hpp>

Go to the source code of this file.

Functions

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

Variables

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

Function Documentation

◆ main()

int main ( int  argc,
char *  argv[] 
)
Examples
continuity_check_on_skeleton_3d.cpp.

Definition at line 30 of file continuity_check_on_skeleton_3d.cpp.

30  {
31 
32  MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
33 
34  try {
35 
36  moab::Core mb_instance;
37  moab::Interface &moab = mb_instance;
38  int rank;
39  MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
40 
41  PetscBool flg = PETSC_TRUE;
42  char mesh_file_name[255];
43 #if PETSC_VERSION_GE(3, 6, 4)
44  CHKERR PetscOptionsGetString(PETSC_NULL, "", "-my_file", mesh_file_name,
45  255, &flg);
46 #else
47  CHKERR PetscOptionsGetString(PETSC_NULL, PETSC_NULL, "-my_file",
48  mesh_file_name, 255, &flg);
49 #endif
50  if (flg != PETSC_TRUE) {
51  SETERRQ(PETSC_COMM_SELF, 1, "*** ERROR -my_file (MESH FILE NEEDED)");
52  }
53 
54  ParallelComm *pcomm = ParallelComm::get_pcomm(&moab, MYPCOMM_INDEX);
55  if (pcomm == NULL)
56  pcomm = new ParallelComm(&moab, PETSC_COMM_WORLD);
57 
58  const char *option;
59  option = "";
60  CHKERR moab.load_file(mesh_file_name, 0, option);
61 
62  // Create MoFEM database
63  MoFEM::Core core(moab);
64  // Access to database through interface
65  MoFEM::Interface &m_field = core;
66 
67  // set entireties bit level
68  BitRefLevel bit_level0 = BitRefLevel().set(0);
69  CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(
70  0, 3, bit_level0);
71 
72  // Fields
73  auto get_base = []() -> FieldApproximationBase {
74  enum bases { AINSWORTH, DEMKOWICZ, LASTBASEOP };
75  const char *list_bases[] = {"ainsworth", "demkowicz"};
76  PetscBool flg;
77  PetscInt choice_base_value = AINSWORTH;
78  CHKERR PetscOptionsGetEList(PETSC_NULL, NULL, "-base", list_bases,
79  LASTBASEOP, &choice_base_value, &flg);
80  if (flg == PETSC_TRUE) {
82  if (choice_base_value == AINSWORTH)
84  else if (choice_base_value == DEMKOWICZ)
85  base = DEMKOWICZ_JACOBI_BASE;
86  return base;
87  }
88  return LASTBASE;
89  };
90 
91  auto get_space = []() -> FieldSpace {
92  enum spaces { HDIV, HCURL, LASTSPACEOP };
93  const char *list_spaces[] = {"hdiv", "hcurl"};
94  PetscBool flg;
95  PetscInt choice_space_value = HDIV;
96  CHKERR PetscOptionsGetEList(PETSC_NULL, NULL, "-space", list_spaces,
97  LASTSPACEOP, &choice_space_value, &flg);
98  if (flg == PETSC_TRUE) {
100  if (choice_space_value == HDIV)
101  space = FieldSpace::HDIV;
102  else if (choice_space_value == HCURL)
103  space = FieldSpace::HCURL;
104  return space;
105  }
106  return LASTSPACE;
107  };
108 
109  CHKERR m_field.add_field("F2", get_space(), get_base(), 1);
110 
111  // meshset consisting all entities in mesh
112  EntityHandle root_set = moab.get_root_set();
113  // add entities to field
114  CHKERR m_field.add_ents_to_field_by_type(root_set, MBTET, "F2");
115 
116  // set app. order
117  int order = 2;
118  CHKERR m_field.set_field_order(root_set, MBTET, "F2", order);
119  CHKERR m_field.set_field_order(root_set, MBTRI, "F2", order);
120  CHKERR m_field.set_field_order(root_set, MBEDGE, "F2", order);
121 
122  CHKERR m_field.build_fields();
123 
124  // add elements
125  CHKERR m_field.add_finite_element("V1");
126  CHKERR m_field.add_finite_element("S2");
127  CHKERR m_field.modify_finite_element_add_field_row("V1", "F2");
128  CHKERR m_field.modify_finite_element_add_field_col("V1", "F2");
129  CHKERR m_field.modify_finite_element_add_field_data("V1", "F2");
130  CHKERR m_field.modify_finite_element_add_field_row("S2", "F2");
131  CHKERR m_field.modify_finite_element_add_field_col("S2", "F2");
132  CHKERR m_field.modify_finite_element_add_field_data("S2", "F2");
133 
134  CHKERR m_field.add_ents_to_finite_element_by_type(root_set, MBTET, "V1");
135  Range faces;
136  CHKERR m_field.getInterface<BitRefManager>()->getEntitiesByTypeAndRefLevel(
137  bit_level0, BitRefLevel().set(), MBTRI, faces);
138  CHKERR m_field.add_ents_to_finite_element_by_type(faces, MBTRI, "S2");
139 
140  CHKERR m_field.build_finite_elements();
141  CHKERR m_field.build_adjacencies(bit_level0);
142 
143  // Problems
144  CHKERR m_field.add_problem("P1");
145 
146  // set refinement level for problem
147  CHKERR m_field.modify_problem_ref_level_add_bit("P1", bit_level0);
148  CHKERR m_field.modify_problem_add_finite_element("P1", "V1");
149  CHKERR m_field.modify_problem_add_finite_element("P1", "S2");
150 
151  // build problems
152  ProblemsManager *prb_mng_ptr;
153  CHKERR m_field.getInterface(prb_mng_ptr);
154  CHKERR prb_mng_ptr->buildProblem("P1", true);
155  CHKERR prb_mng_ptr->partitionProblem("P1");
156  CHKERR prb_mng_ptr->partitionFiniteElements("P1");
157  CHKERR prb_mng_ptr->partitionGhostDofs("P1");
158 
159  struct CommonData {
160  MatrixDouble dotNormalFace;
161  MatrixDouble dotNormalEleLeft;
162  MatrixDouble dotNormalEleRight;
163  };
164 
165  struct SkeletonFE
167 
169  struct OpVolSide
171 
173  OpVolSide(CommonData &elem_data)
175  "F2", UserDataOperator::OPROW),
176  elemData(elem_data) {}
177  MoFEMErrorCode doWork(int side, EntityType type,
180 
181  if (type == MBTRI && side == getFaceSideNumber()) {
182 
183  MatrixDouble diff =
184  getCoordsAtGaussPts() - getFaceCoordsAtGaussPts();
185  const double eps = 1e-12;
186  if (norm_inf(diff) > eps)
187  SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
188  "coordinates at integration pts are different");
189 
190  const size_t nb_dofs = data.getN().size2() / 3;
191  const size_t nb_integration_pts = data.getN().size1();
192 
194  auto t_to_do_dot = getFTensor1Normal();
195  if (data.getSpace() == HCURL) {
196  auto s1 = getFTensor1Tangent1();
197  auto s2 = getFTensor1Tangent2();
198  t_to_do_dot(i) = s1(i) + s2(i);
199  }
200 
201  auto t_hdiv_base = data.getFTensor1N<3>();
202  MatrixDouble *ptr_dot_elem_data = nullptr;
203  if (getFEMethod()->nInTheLoop == 0)
204  ptr_dot_elem_data = &elemData.dotNormalEleLeft;
205  else
206  ptr_dot_elem_data = &elemData.dotNormalEleRight;
207  MatrixDouble &dot_elem_data = *ptr_dot_elem_data;
208  dot_elem_data.resize(nb_integration_pts, nb_dofs, false);
209 
210  for (size_t gg = 0; gg != nb_integration_pts; ++gg) {
211  for (size_t bb = 0; bb != nb_dofs; ++bb) {
212  dot_elem_data(gg, bb) = t_to_do_dot(i) * t_hdiv_base(i);
213  ++t_hdiv_base;
214  }
215  }
216  }
218  }
219  };
220 
222  SkeletonFE(MoFEM::Interface &m_field, CommonData &elem_data)
224  "F2", UserDataOperator::OPROW),
225  volSideFe(m_field), elemData(elem_data) {
226  volSideFe.getOpPtrVector().push_back(
227  new SkeletonFE::OpVolSide(elemData));
228  }
229 
230  MoFEMErrorCode doWork(int side, EntityType type,
232 
234  if (type == MBTRI && side == 0) {
235  const size_t nb_dofs = data.getN().size2() / 3;
236  const size_t nb_integration_pts = data.getN().size1();
237 
239  auto t_to_do_dot = getFTensor1Normal();
240  if (data.getSpace() == HCURL) {
241  auto s1 = getFTensor1Tangent1();
242  auto s2 = getFTensor1Tangent2();
243  t_to_do_dot(i) = s1(i) + s2(i);
244  }
245 
246  auto t_hdiv_base = data.getFTensor1N<3>();
247  elemData.dotNormalFace.resize(nb_integration_pts, nb_dofs, false);
248  elemData.dotNormalEleLeft.resize(nb_integration_pts, 0, false);
249  elemData.dotNormalEleRight.resize(nb_integration_pts, 0, false);
250 
251  for (size_t gg = 0; gg != nb_integration_pts; ++gg) {
252  for (size_t bb = 0; bb != nb_dofs; ++bb) {
253  elemData.dotNormalFace(gg, bb) = t_to_do_dot(i) * t_hdiv_base(i);
254  ++t_hdiv_base;
255  }
256  }
257 
258  CHKERR loopSideVolumes("V1", volSideFe);
259 
260  auto check_continuity_of_base = [&](auto &vol_dot_data) {
262 
263  if (vol_dot_data.size1() != elemData.dotNormalFace.size1())
264  SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
265  "Inconsistent number of integration points");
266 
267  if (vol_dot_data.size2() != elemData.dotNormalFace.size2())
268  SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
269  "Inconsistent number of base functions");
270  const double eps = 1e-12;
271  for (size_t gg = 0; gg != vol_dot_data.size1(); ++gg)
272  for (size_t bb = 0; bb != vol_dot_data.size2(); ++bb) {
273  const double error = std::abs(vol_dot_data(gg, bb) -
274  elemData.dotNormalFace(gg, bb));
275  if (error > eps)
276  SETERRQ2(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
277  "Inconsistency %3.4e != %3.4e", vol_dot_data(gg, bb),
278  elemData.dotNormalFace(gg, bb));
279  }
281  };
282 
283  if (elemData.dotNormalEleLeft.size2() != 0)
284  CHKERR check_continuity_of_base(elemData.dotNormalEleLeft);
285  else if (elemData.dotNormalEleRight.size2() != 0)
286  CHKERR check_continuity_of_base(elemData.dotNormalEleRight);
287  }
289  }
290  };
291 
292  CommonData elem_data;
293  FaceElementForcesAndSourcesCore face_fe(m_field);
294  face_fe.getOpPtrVector().push_back(new SkeletonFE(m_field, elem_data));
295  CHKERR m_field.loop_finite_elements("P1", "S2", face_fe);
296  }
297  CATCH_ERRORS;
298 
300 
301  return 0;
302 }

Variable Documentation

◆ help

char help[] = "...\n\n"
static
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:516
MoFEM::PetscOptionsGetEList
PetscErrorCode PetscOptionsGetEList(PetscOptions *, const char pre[], const char name[], const char *const *list, PetscInt next, PetscInt *value, PetscBool *set)
Definition: DeprecatedPetsc.hpp:214
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:26
EntityHandle
SkeletonFE::SkeletonFE
SkeletonFE(MoFEM::Interface &m_field, CommonData &elem_data)
Definition: continuity_check_on_skeleton_with_simple_2d_for_h1.cpp:103
MoFEM::CoreInterface::build_adjacencies
virtual MoFEMErrorCode build_adjacencies(const Range &ents, int verb=DEFAULT_VERBOSITY)=0
build adjacencies
MYPCOMM_INDEX
#define MYPCOMM_INDEX
default communicator number PCOMM
Definition: definitions.h:292
MoFEM::ProblemsManager::buildProblem
MoFEMErrorCode buildProblem(const std::string name, const bool square_matrix, int verb=VERBOSE)
build problem data structures
Definition: ProblemsManager.cpp:428
MoFEM::DataForcesAndSourcesCore::EntData::getSpace
FieldSpace & getSpace()
Get field space.
Definition: DataStructures.hpp:1285
MoFEM::VolumeElementForcesAndSourcesCoreOnSideBase::UserDataOperator
default operator for TET element
Definition: VolumeElementForcesAndSourcesCoreOnSide.hpp:86
LASTBASE
@ LASTBASE
Definition: definitions.h:161
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:67
MoFEM::Types::MatrixDouble
ublas::matrix< double, ublas::row_major, DoubleAllocator > MatrixDouble
Definition: Types.hpp:76
MoFEM::VolumeElementForcesAndSourcesCoreOnSideSwitch
Volume side finite element with switches.
Definition: VolumeElementForcesAndSourcesCoreOnSide.hpp:198
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::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.
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::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:1943
MoFEM::PetscOptionsGetString
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
Definition: DeprecatedPetsc.hpp:183
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.
eps
constexpr double eps
Definition: forces_and_sources_getting_mult_H1_H1_atom_test.cpp:30
FieldSpace
FieldSpace
approximation spaces
Definition: definitions.h:174
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
MoFEM::ForcesAndSourcesCore::UserDataOperator
Data operator to do calculations at integration points.
Definition: ForcesAndSourcesCore.hpp:84
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:604
MoFEM::DataForcesAndSourcesCore::EntData::getFTensor1N
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1N(FieldApproximationBase base)
Get base functions for Hdiv/Hcurl spaces.
Definition: DataStructures.cpp:587
UserDataOperator
ForcesAndSourcesCore::UserDataOperator UserDataOperator
Definition: HookeElement.hpp:88
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::build_finite_elements
virtual MoFEMErrorCode build_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elements.
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
convert.type
type
Definition: convert.py:66
FTensor::Index< 'i', 3 >
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:60
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_filed)=0
set finite element field data
StdRDOperators::order
const int order
Definition: rd_stdOperators.hpp:28
LASTSPACE
@ LASTSPACE
FieldSpace in [ 0, LASTSPACE )
Definition: definitions.h:181
SkeletonFE::elemData
CommonData & elemData
Definition: continuity_check_on_skeleton_with_simple_2d_for_h1.cpp:100
MoFEM::ProblemsManager
Problem manager is used to build and partition problems.
Definition: ProblemsManager.hpp:34
FaceElementForcesAndSourcesCore
MoFEM::Types::BitRefLevel
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition: Types.hpp:51
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:2416
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
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:77
MoFEM::CoreTmp< 0 >::Finalize
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:100
mesh_file_name
char mesh_file_name[255]
Definition: mesh_smoothing.cpp:35
CommonData
Definition: continuity_check_on_skeleton_with_simple_2d_for_h1.cpp:40
MoFEM::DeprecatedCoreInterface::loop_finite_elements
DEPRECATED MoFEMErrorCode loop_finite_elements(const Problem *problem_ptr, const std::string &fe_name, FEMethod &method, int lower_rank, int upper_rank, MoFEMTypes bh, int verb=DEFAULT_VERBOSITY)
Definition: DeprecatedCoreInterface.cpp:579
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:441
DEMKOWICZ_JACOBI_BASE
@ DEMKOWICZ_JACOBI_BASE
Definition: definitions.h:158
SkeletonFE
Definition: continuity_check_on_skeleton_with_simple_2d_for_h1.cpp:46
help
static char help[]
Definition: continuity_check_on_skeleton_3d.cpp:28
MoFEM::BitRefManager
Managing BitRefLevels.
Definition: BitRefManager.hpp:34
AINSWORTH_LEGENDRE_BASE
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:152
HCURL
@ HCURL
field with continuous tangents
Definition: definitions.h:178
MoFEM::DataForcesAndSourcesCore::EntData::getN
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
Definition: DataStructures.hpp:1288
MoFEM::CoreInterface::build_fields
virtual MoFEMErrorCode build_fields(int verb=DEFAULT_VERBOSITY)=0
FieldApproximationBase
FieldApproximationBase
approximation base
Definition: definitions.h:150
MOFEM_ATOM_TEST_INVALID
@ MOFEM_ATOM_TEST_INVALID
Definition: definitions.h:132
MoFEM::DataForcesAndSourcesCore::EntData
Data on single entity (This is passed as argument to DataOperator::doWork)
Definition: DataStructures.hpp:60
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:509
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.
i
FTensor::Index< 'i', 3 > i
Definition: matrix_function.cpp:18
MoFEM::ProblemsManager::partitionGhostDofs
MoFEMErrorCode partitionGhostDofs(const std::string name, int verb=VERBOSE)
determine ghost nodes
Definition: ProblemsManager.cpp:2599
SkeletonFE::doWork
MoFEMErrorCode doWork(int side, EntityType type, DataForcesAndSourcesCore::EntData &data)
Definition: continuity_check_on_skeleton_with_simple_2d_for_h1.cpp:109
MoFEM::CoreInterface::add_problem
virtual MoFEMErrorCode add_problem(const std::string &name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add problem.
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:485
HDIV
@ HDIV
field with continuous normal traction
Definition: definitions.h:179
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1129
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(const MOFEMuuid &uuid, IFACE *&iface) const
Get interface by uuid and return reference to pointer of interface.
Definition: UnknownInterface.hpp:130
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:415
MoFEM::ProblemsManager::partitionProblem
MoFEMErrorCode partitionProblem(const std::string name, int verb=VERBOSE)
partition problem dofs (collective)
Definition: ProblemsManager.cpp:1954