v0.14.0
forces_and_sources_getting_higher_order_skin_normals_atom_tets.cpp
Go to the documentation of this file.
1 /** \file forces_and_sources_getting_higher_order_skin_normals_atom_tets.cpp
2 
3  \brief Atom test to calculate normals on faces approximated with HO geometry
4  representation
5 
6 */
7 
8 
9 
10 #include <MoFEM.hpp>
11 
12 namespace bio = boost::iostreams;
13 using bio::stream;
14 using bio::tee_device;
15 
16 using namespace MoFEM;
17 
18 static char help[] = "...\n\n";
19 
20 int main(int argc, char *argv[]) {
21 
22  MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
23 
24  try {
25 
26  moab::Core mb_instance;
27  moab::Interface &moab = mb_instance;
28  int rank;
29  MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
30 
31  PetscBool flg = PETSC_TRUE;
32  char mesh_file_name[255];
33 #if PETSC_VERSION_GE(3, 6, 4)
34  CHKERR PetscOptionsGetString(PETSC_NULL, "", "-my_file", mesh_file_name,
35  255, &flg);
36 #else
37  CHKERR PetscOptionsGetString(PETSC_NULL, PETSC_NULL, "-my_file",
38  mesh_file_name, 255, &flg);
39 #endif
40  if (flg != PETSC_TRUE) {
41  SETERRQ(PETSC_COMM_SELF, 1, "*** ERROR -my_file (MESH FILE NEEDED)");
42  }
43 
44  // Create MoFEM (Joseph) database
45  MoFEM::Core core(moab);
46  MoFEM::Interface &m_field = core;
47 
48  const char *option;
49  option = "";
50  CHKERR moab.load_file(mesh_file_name, 0, option);
51 
52  // set entitities bit level
53  BitRefLevel bit_level0;
54  bit_level0.set(0);
55  EntityHandle meshset_level0;
56  CHKERR moab.create_meshset(MESHSET_SET, meshset_level0);
57  CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(
58  0, 3, bit_level0);
59 
60  // Fields
61  CHKERR m_field.add_field("FIELD1", H1, AINSWORTH_LEGENDRE_BASE, 3);
62 
63  // FE
64  CHKERR m_field.add_finite_element("TEST_FE");
65 
66  // Define rows/cols and element data
67  CHKERR m_field.modify_finite_element_add_field_row("TEST_FE", "FIELD1");
68  CHKERR m_field.modify_finite_element_add_field_col("TEST_FE", "FIELD1");
69  CHKERR m_field.modify_finite_element_add_field_data("TEST_FE", "FIELD1");
70 
71  // Problem
72  CHKERR m_field.add_problem("TEST_PROBLEM");
73 
74  // set finite elements for problem
75  CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM", "TEST_FE");
76  // set refinement level for problem
77  CHKERR m_field.modify_problem_ref_level_add_bit("TEST_PROBLEM", bit_level0);
78 
79  // meshset consisting all entities in mesh
80  EntityHandle root_set = moab.get_root_set();
81  // add entities to field
82  CHKERR m_field.add_ents_to_field_by_type(root_set, MBTET, "FIELD1");
83  // add entities to finite element
84  Range tets;
85  CHKERR moab.get_entities_by_type(0, MBTET, tets, false);
86  Skinner skin(&m_field.get_moab());
87  Range tets_skin;
88  CHKERR skin.find_skin(0, tets, false, tets_skin);
89  CHKERR m_field.add_ents_to_finite_element_by_type(tets_skin, MBTRI,
90  "TEST_FE");
91 
92  // set app. order
93  // see Hierarchic Finite Element Bases on Unstructured Tetrahedral Meshes
94  // (Mark Ainsworth & Joe Coyle)
95  int order = 3;
96  CHKERR m_field.set_field_order(root_set, MBTET, "FIELD1", order);
97  CHKERR m_field.set_field_order(root_set, MBTRI, "FIELD1", order);
98  CHKERR m_field.set_field_order(root_set, MBEDGE, "FIELD1", order);
99  CHKERR m_field.set_field_order(root_set, MBVERTEX, "FIELD1", 1);
100 
101  /****/
102  // build database
103  // build field
104  CHKERR m_field.build_fields();
105  // set FIELD1 from positions of 10 node tets
106  Projection10NodeCoordsOnField ent_method(m_field, "FIELD1");
107  CHKERR m_field.loop_dofs("FIELD1", ent_method);
108  // build finite elemnts
109  CHKERR m_field.build_finite_elements();
110  // build adjacencies
111  CHKERR m_field.build_adjacencies(bit_level0);
112  // build problem
113  ProblemsManager *prb_mng_ptr;
114  CHKERR m_field.getInterface(prb_mng_ptr);
115  CHKERR prb_mng_ptr->buildProblem("TEST_PROBLEM", true);
116 
117  /****/
118  // mesh partitioning
119  // partition
120  CHKERR prb_mng_ptr->partitionSimpleProblem("TEST_PROBLEM");
121  CHKERR prb_mng_ptr->partitionFiniteElements("TEST_PROBLEM");
122  // what are ghost nodes, see Petsc Manual
123  CHKERR prb_mng_ptr->partitionGhostDofs("TEST_PROBLEM");
124 
125  struct ForcesAndSourcesCore_TestFE
127 
128  typedef tee_device<std::ostream, std::ofstream> TeeDevice;
129  typedef stream<TeeDevice> TeeStream;
130  std::ofstream ofs;
131  TeeDevice my_tee;
132  TeeStream my_split;
133 
134  ForcesAndSourcesCore_TestFE(MoFEM::Interface &m_field)
136  ofs("forces_and_sources_getting_higher_order_skin_normals_atom."
137  "txt"),
138  my_tee(std::cout, ofs), my_split(my_tee){};
139 
140  MoFEMErrorCode preProcess() {
143  }
144 
145  MoFEMErrorCode operator()() {
147 
149 
150  my_split.precision(3);
151  my_split << "coords: " << coordsAtGaussPts << std::endl;
152  my_split << "normals: " << normalsAtGaussPts << std::endl;
153  my_split << "tangent1: " << tangentOneAtGaussPts << std::endl;
154  my_split << "tangent2: " << tangentTwoAtGaussPts << std::endl;
155 
157  }
158 
159  MoFEMErrorCode postProcess() {
162  }
163  };
164 
165  ForcesAndSourcesCore_TestFE fe1(m_field);
166 
167  fe1.getRuleHook = [](int, int, int approx_order) { return -1; };
168  fe1.setRuleHook = [](ForcesAndSourcesCore *fe_ptr, int ro, int co, int ao) {
170  auto &gauss_pts = fe_ptr->gaussPts;
171  gauss_pts.resize(3, 4, false);
172  for (int gg = 0; gg < 4; gg++) {
173  gauss_pts(0, gg) = G_TRI_X4[gg];
174  gauss_pts(1, gg) = G_TRI_Y4[gg];
175  gauss_pts(2, gg) = 0;
176  }
178  };
179 
180  fe1.getOpPtrVector().push_back(new OpCalculateHOCoords("FIELD1"));
181  fe1.getOpPtrVector().push_back(new OpGetHONormalsOnFace("FIELD1"));
182 
183  CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "TEST_FE", fe1);
184  }
185  CATCH_ERRORS;
186 
188 
189  return 0;
190 }
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.
Definition: UnknownInterface.hpp:93
TeeStream
stream< TeeDevice > TeeStream
Definition: forces_and_sources_testing_contact_prism_element.cpp:10
MoFEM::CoreInterface::loop_finite_elements
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.
MoFEM::CoreInterface::loop_dofs
virtual MoFEMErrorCode loop_dofs(const Problem *problem_ptr, const std::string &field_name, RowColData rc, DofMethod &method, int lower_rank, int upper_rank, int verb=DEFAULT_VERBOSITY)=0
Make a loop over dofs.
MoFEM::ProblemsManager::buildProblem
MoFEMErrorCode buildProblem(const std::string name, const bool square_matrix, int verb=VERBOSE)
build problem data structures
Definition: ProblemsManager.cpp:87
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
H1
@ H1
continuous field
Definition: definitions.h:85
EntityHandle
MoFEM::ProblemsManager
Problem manager is used to build and partition problems.
Definition: ProblemsManager.hpp:21
MoFEM::FaceElementForcesAndSourcesCore::operator()
MoFEMErrorCode operator()()
function is run for every finite element
Definition: FaceElementForcesAndSourcesCore.cpp:401
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
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
help
static char help[]
Definition: forces_and_sources_getting_higher_order_skin_normals_atom_tets.cpp:18
MoFEM.hpp
MoFEM::Projection10NodeCoordsOnField
Projection of edge entities with one mid-node on hierarchical basis.
Definition: Projection10NodeCoordsOnField.hpp:24
MoFEM::CoreTmp< 0 >::Finalize
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:112
main
int main(int argc, char *argv[])
Definition: forces_and_sources_getting_higher_order_skin_normals_atom_tets.cpp:20
MoFEM::OpCalculateHOCoords
Calculate HO coordinates at gauss points.
Definition: HODataOperators.hpp:41
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.
order
constexpr int order
Definition: dg_projection.cpp:18
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
TeeDevice
tee_device< std::ostream, std::ofstream > TeeDevice
Definition: forces_and_sources_testing_contact_prism_element.cpp:9
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:1975
approx_order
static constexpr int approx_order
Definition: prism_polynomial_approximation.cpp:14
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
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:535
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::get_moab
virtual moab::Interface & get_moab()=0
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
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::build_finite_elements
virtual MoFEMErrorCode build_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elements.
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::FaceElementForcesAndSourcesCore
Face finite element.
Definition: FaceElementForcesAndSourcesCore.hpp:23
MoFEM::OpGetHONormalsOnFace
Calculate normals at Gauss points of triangle element.
Definition: HODataOperators.hpp:281
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:2167
G_TRI_X4
static const double G_TRI_X4[]
Definition: fem_tools.h:319
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
mesh_file_name
char mesh_file_name[255]
Definition: mesh_smoothing.cpp:23
G_TRI_Y4
static const double G_TRI_Y4[]
Definition: fem_tools.h:322
MoFEM::ForcesAndSourcesCore
structure to get information form mofem into EntitiesFieldData
Definition: ForcesAndSourcesCore.hpp:22
Range
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:72
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:372
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1094
MoFEM::ProblemsManager::partitionSimpleProblem
MoFEMErrorCode partitionSimpleProblem(const std::string name, int verb=VERBOSE)
partition problem dofs
Definition: ProblemsManager.cpp:1537
AINSWORTH_LEGENDRE_BASE
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
MoFEM::ForcesAndSourcesCore::gaussPts
MatrixDouble gaussPts
Matrix of integration points.
Definition: ForcesAndSourcesCore.hpp:109
MoFEM::PetscOptionsGetString
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
Definition: DeprecatedPetsc.hpp:172
MoFEM::CoreInterface::build_fields
virtual MoFEMErrorCode build_fields(int verb=DEFAULT_VERBOSITY)=0
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
MoFEM::BitRefManager
Managing BitRefLevels.
Definition: BitRefManager.hpp:21
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::Types::BitRefLevel
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition: Types.hpp:40
MoFEM::CoreInterface::build_adjacencies
virtual MoFEMErrorCode build_adjacencies(const Range &ents, int verb=DEFAULT_VERBOSITY)=0
build adjacencies
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:440
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.
MoFEM::ProblemsManager::partitionGhostDofs
MoFEMErrorCode partitionGhostDofs(const std::string name, int verb=VERBOSE)
determine ghost nodes
Definition: ProblemsManager.cpp:2339
MoFEM::CoreInterface::add_problem
virtual MoFEMErrorCode add_problem(const std::string &name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add problem.
convert.int
int
Definition: convert.py:64
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346