v0.13.2
Loading...
Searching...
No Matches
Functions | Variables
forces_and_sources_getting_higher_order_skin_normals_atom_tets.cpp File Reference

Atom test to calculate normals on faces approximated with HO geometry representation. More...

#include <MoFEM.hpp>

Go to the source code of this file.

Functions

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

Variables

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

Detailed Description

Atom test to calculate normals on faces approximated with HO geometry representation.

Definition in file forces_and_sources_getting_higher_order_skin_normals_atom_tets.cpp.

Function Documentation

◆ main()

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

Definition at line 20 of file forces_and_sources_getting_higher_order_skin_normals_atom_tets.cpp.

20 {
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
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
143 }
144
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
162 }
163 };
164
165 ForcesAndSourcesCore_TestFE fe1(m_field);
166 fe1.meshPositionsFieldName = "FIELD1";
167
168 fe1.getRuleHook = [](int, int, int approx_order) { return -1; };
169 fe1.setRuleHook = [](ForcesAndSourcesCore *fe_ptr, int ro, int co, int ao) {
171 auto &gauss_pts = fe_ptr->gaussPts;
172 gauss_pts.resize(3, 4, false);
173 for (int gg = 0; gg < 4; gg++) {
174 gauss_pts(0, gg) = G_TRI_X4[gg];
175 gauss_pts(1, gg) = G_TRI_Y4[gg];
176 gauss_pts(2, gg) = 0;
177 }
179 };
180
181 fe1.getOpPtrVector().push_back(new OpCalculateHOCoords("FIELD1"));
182 fe1.getOpPtrVector().push_back(new OpGetHONormalsOnFace("FIELD1"));
183
184 CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "TEST_FE", fe1);
185 }
187
189
190 return 0;
191}
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:372
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
@ H1
continuous field
Definition: definitions.h:85
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
#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
static const double G_TRI_X4[]
Definition: fem_tools.h:319
static const double G_TRI_Y4[]
Definition: fem_tools.h:322
tee_device< std::ostream, std::ofstream > TeeDevice
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_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.
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
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
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
Definition: enable_if.hpp:6
static constexpr int approx_order
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.
MoFEMErrorCode operator()()
function is run for every finite element
structure to get information form mofem into EntitiesFieldData
MatrixDouble gaussPts
Matrix of integration points.
virtual MoFEMErrorCode preProcess()
function is run at the beginning of loop
virtual MoFEMErrorCode postProcess()
function is run at the end of loop
Calculate HO coordinates at gauss points.
Calculate normals at Gauss points of triangle element.
Problem manager is used to build and partition problems.
Projection of edge entities with one mid-node on hierarchical basis.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.

Variable Documentation

◆ help

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