v0.13.2
Loading...
Searching...
No Matches
Functions | Variables
forces_and_sources_hcurl_approximation_functions.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[] 
)

Definition at line 13 of file forces_and_sources_hcurl_approximation_functions.cpp.

13 {
14
15 MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
16
17 try {
18
19 moab::Core mb_instance;
20 moab::Interface &moab = mb_instance;
21 int rank;
22 MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
23
24 // create one tet
25 double tet_coords[] = {0, 0, 0, 2, 0, 0, 0, 2, 0, 0, 0, 2};
26 EntityHandle nodes[4];
27 for (int nn = 0; nn < 4; nn++) {
28 CHKERR moab.create_vertex(&tet_coords[3 * nn], nodes[nn]);
29 }
30 EntityHandle tet;
31 CHKERR moab.create_element(MBTET, nodes, 4, tet);
32
33 // Create MoFEM (Joseph) database
34 MoFEM::Core core(moab);
35 MoFEM::Interface &m_field = core;
36
37 // set entitities bit level
38 BitRefLevel bit_level0;
39 bit_level0.set(0);
40 CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(
41 0, 3, bit_level0);
42
43 // Fields
44 CHKERR m_field.add_field("HCURL", HCURL, AINSWORTH_LEGENDRE_BASE, 1);
45
46 // FE TET
47 CHKERR m_field.add_finite_element("HCURL_TET_FE");
48 // Define rows/cols and element data
49 CHKERR m_field.modify_finite_element_add_field_row("HCURL_TET_FE", "HCURL");
50 CHKERR m_field.modify_finite_element_add_field_col("HCURL_TET_FE", "HCURL");
51 CHKERR m_field.modify_finite_element_add_field_data("HCURL_TET_FE",
52 "HCURL");
53
54 // FE TRI
55 CHKERR m_field.add_finite_element("HCURL_TRI_FE");
56 // Define rows/cols and element data
57 CHKERR m_field.modify_finite_element_add_field_row("HCURL_TRI_FE", "HCURL");
58 CHKERR m_field.modify_finite_element_add_field_col("HCURL_TRI_FE", "HCURL");
59 CHKERR m_field.modify_finite_element_add_field_data("HCURL_TRI_FE",
60 "HCURL");
61
62 // FE EDGE
63 CHKERR m_field.add_finite_element("HCURL_EDGE_FE");
64 // Define rows/cols and element data
65 CHKERR m_field.modify_finite_element_add_field_row("HCURL_EDGE_FE",
66 "HCURL");
67 CHKERR m_field.modify_finite_element_add_field_col("HCURL_EDGE_FE",
68 "HCURL");
69 CHKERR m_field.modify_finite_element_add_field_data("HCURL_EDGE_FE",
70 "HCURL");
71
72 // Problem
73 CHKERR m_field.add_problem("TEST_PROBLEM");
74
75 // set finite elements for problem
76 CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM",
77 "HCURL_TET_FE");
78 CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM",
79 "HCURL_TRI_FE");
80 CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM",
81 "HCURL_EDGE_FE");
82 // set refinement level for problem
83 CHKERR m_field.modify_problem_ref_level_add_bit("TEST_PROBLEM", bit_level0);
84
85 // meshset consisting all entities in mesh
86 EntityHandle root_set = moab.get_root_set();
87 // add entities to field
88 CHKERR m_field.add_ents_to_field_by_type(root_set, MBTET, "HCURL");
89 // add entities to finite element
90 CHKERR m_field.add_ents_to_finite_element_by_type(root_set, MBTET,
91 "HCURL_TET_FE");
92
93 Range tets;
94 CHKERR m_field.getInterface<BitRefManager>()->getEntitiesByTypeAndRefLevel(
95 BitRefLevel().set(0), BitRefLevel().set(), MBTET, tets);
96 Skinner skin(&moab);
97 Range skin_faces; // skin faces from 3d ents
98 CHKERR skin.find_skin(0, tets, false, skin_faces);
99 CHKERR m_field.add_ents_to_finite_element_by_type(skin_faces, MBTRI,
100 "HCURL_TRI_FE");
101 Range skin_edges;
102 CHKERR moab.get_adjacencies(skin_faces, 1, false, skin_edges,
103 moab::Interface::UNION);
104 CHKERR m_field.add_ents_to_finite_element_by_type(skin_edges, MBEDGE,
105 "HCURL_EDGE_FE");
106
107 // set app. order
108 int order = 4;
109 CHKERR m_field.set_field_order(root_set, MBTET, "HCURL", order);
110 CHKERR m_field.set_field_order(root_set, MBTRI, "HCURL", order);
111 CHKERR m_field.set_field_order(root_set, MBEDGE, "HCURL", order);
112
113 // build database
114
115 // build field
116 CHKERR m_field.build_fields();
117 // build finite elemnts
119 // build adjacencies
120 CHKERR m_field.build_adjacencies(bit_level0);
121 // build problem
122 ProblemsManager *prb_mng_ptr;
123 CHKERR m_field.getInterface(prb_mng_ptr);
124 CHKERR prb_mng_ptr->buildProblem("TEST_PROBLEM", true);
125
126 // mesh partitioning
127
128 // partition
129 CHKERR prb_mng_ptr->partitionSimpleProblem("TEST_PROBLEM");
130 CHKERR prb_mng_ptr->partitionFiniteElements("TEST_PROBLEM");
131 // what are ghost nodes, see Petsc Manual
132 CHKERR prb_mng_ptr->partitionGhostDofs("TEST_PROBLEM");
133
134 typedef tee_device<std::ostream, std::ofstream> TeeDevice;
135 typedef stream<TeeDevice> TeeStream;
136
137 std::ofstream ofs("forces_and_sources_hcurl_approximation_functions.txt");
138 TeeDevice my_tee(std::cout, ofs);
139 TeeStream my_split(my_tee);
140
141 struct OpPrintingHdivApproximationFunctions
142 : public VolumeElementForcesAndSourcesCore::UserDataOperator {
143
144 TeeStream &mySplit;
145 OpPrintingHdivApproximationFunctions(TeeStream &my_split)
147 "HCURL", UserDataOperator::OPROW),
148 mySplit(my_split) {}
149
150 MoFEMErrorCode doWork(int side, EntityType type,
153
154 if (data.getFieldData().size() == 0)
156
157 mySplit << std::endl
158 << "type " << type << " side " << side << std::endl;
159
160 const double eps = 1e-4;
161 for (unsigned int dd = 0; dd < data.getN().data().size(); dd++) {
162 if (std::abs(data.getN().data()[dd]) < eps)
163 data.getN().data()[dd] = 0;
164 }
165 for (unsigned int dd = 0; dd < data.getDiffN().data().size(); dd++) {
166 if (std::abs(data.getDiffN().data()[dd]) < eps)
167 data.getDiffN().data()[dd] = 0;
168 }
169
170 mySplit << std::fixed << std::setprecision(5) << data.getN()
171 << std::endl;
172 mySplit << std::fixed << std::setprecision(5) << data.getDiffN()
173 << std::endl;
174
176 }
177 };
178
179 struct MyFE : public VolumeElementForcesAndSourcesCore {
180
181 MyFE(MoFEM::Interface &m_field)
183 int getRule(int order) { return 1; };
184 };
185
186 struct OpFacePrintingHdivApproximationFunctions
187 : public FaceElementForcesAndSourcesCore::UserDataOperator {
188
189 TeeStream &mySplit;
190 OpFacePrintingHdivApproximationFunctions(TeeStream &my_split)
192 "HCURL", UserDataOperator::OPROW),
193 mySplit(my_split) {}
194
195 MoFEMErrorCode doWork(int side, EntityType type,
198
199 if (data.getFieldData().size() == 0)
201
202 mySplit << std::endl
203 << "type " << type << " side " << side << std::endl;
204 mySplit.precision(5);
205
206 const double eps = 1e-4;
207 for (unsigned int dd = 0; dd < data.getN().data().size(); dd++) {
208 if (std::abs(data.getN().data()[dd]) < eps)
209 data.getN().data()[dd] = 0;
210 }
211 for (unsigned int dd = 0; dd < data.getDiffN().data().size(); dd++) {
212 if (std::abs(data.getDiffN().data()[dd]) < eps)
213 data.getDiffN().data()[dd] = 0;
214 }
215
216 mySplit << std::fixed << std::setprecision(5) << data.getN()
217 << std::endl;
218 mySplit << std::fixed << std::setprecision(5) << data.getDiffN()
219 << std::endl;
220
222 }
223 };
224
225 struct MyFaceFE : public FaceElementForcesAndSourcesCore {
226
227 MyFaceFE(MoFEM::Interface &m_field)
229 int getRule(int order) { return 1; };
230 };
231
232 struct OpEdgePrintingHdivApproximationFunctions
233 : public EdgeElementForcesAndSourcesCore::UserDataOperator {
234
235 TeeStream &mySplit;
236 OpEdgePrintingHdivApproximationFunctions(TeeStream &my_split)
238 "HCURL", UserDataOperator::OPROW),
239 mySplit(my_split) {}
240
241 MoFEMErrorCode doWork(int side, EntityType type,
244
245 if (data.getFieldData().size() == 0)
247
248 mySplit << std::endl
249 << "type " << type << " side " << side << std::endl;
250 mySplit.precision(5);
251
252 const double eps = 1e-4;
253 for (unsigned int dd = 0; dd < data.getN().data().size(); dd++) {
254 if (std::abs(data.getN().data()[dd]) < eps)
255 data.getN().data()[dd] = 0;
256 }
257
258 mySplit << std::fixed << std::setprecision(5) << data.getN()
259 << std::endl;
260
262 }
263 };
264
265 struct MyEdgeFE : public EdgeElementForcesAndSourcesCore {
266
267 MyEdgeFE(MoFEM::Interface &m_field)
269 int getRule(int order) { return 1; };
270 };
271
272 MyFE tet_fe(m_field);
273 MyFaceFE tri_fe(m_field);
274 MyEdgeFE edge_fe(m_field);
275
276 tet_fe.getOpPtrVector().push_back(
277 new OpPrintingHdivApproximationFunctions(my_split));
278 tri_fe.getOpPtrVector().push_back(
280 tri_fe.getOpPtrVector().push_back(
281 new OpFacePrintingHdivApproximationFunctions(my_split));
282 edge_fe.getOpPtrVector().push_back(
284 edge_fe.getOpPtrVector().push_back(
285 new OpEdgePrintingHdivApproximationFunctions(my_split));
286
287 CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "HCURL_TET_FE", tet_fe);
288 CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "HCURL_TRI_FE", tri_fe);
289 CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "HCURL_EDGE_FE",
290 edge_fe);
291
292
293 }
295
297}
ForcesAndSourcesCore::UserDataOperator UserDataOperator
static const double eps
#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
@ HCURL
field with continuous tangents
Definition: definitions.h:86
#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
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_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
const Tensor2_symmetric_Expr< const ddTensor0< T, Dim, i, j >, typename promote< T, double >::V, Dim, i, j > dd(const Tensor0< T * > &a, const Index< i, Dim > index1, const Index< j, Dim > index2, const Tensor1< int, Dim > &d_ijk, const Tensor1< double, Dim > &d_xyz)
Definition: ddTensor0.hpp:33
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
Managing BitRefLevels.
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.
Data on single entity (This is passed as argument to DataOperator::doWork)
MatrixDouble & getDiffN(const FieldApproximationBase base)
get derivatives of base functions
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
const VectorDouble & getFieldData() const
get dofs values
virtual int getRule(int order_row, int order_col, int order_data)
another variant of getRule
transform Hcurl base fluxes from reference element to physical edge
transform Hcurl base fluxes from reference element to physical triangle
Problem manager is used to build and partition problems.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.

Variable Documentation

◆ help

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