v0.13.1
Functions | Variables
check_base_functions_derivatives_on_tet.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"
 
static const double eps = 1e-6
 
static const double eps_diff = 1e-8
 

Function Documentation

◆ main()

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

Definition at line 28 of file check_base_functions_derivatives_on_tet.cpp.

28 {
29
30 MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
31
32 try {
33
34 // Select space
35 enum spaces { H1TET, HDIVTET, HCURLTET, LASTSPACEOP };
36
37 const char *list_spaces[] = {"h1tet", "hdivtet", "hcurltet"};
38
39 PetscBool flg;
40 PetscInt choise_space_value = H1TET;
41 CHKERR PetscOptionsGetEList(PETSC_NULL, NULL, "-space", list_spaces,
42 LASTSPACEOP, &choise_space_value, &flg);
43
44 if (flg != PETSC_TRUE) {
45 SETERRQ(PETSC_COMM_SELF, MOFEM_IMPOSIBLE_CASE, "space not set");
46 }
47
48 FieldSpace space = LASTSPACE;
49 if (choise_space_value == H1TET) {
50 space = H1;
51 } else if (choise_space_value == HDIVTET) {
52 space = HDIV;
53 } else if (choise_space_value == HCURLTET) {
54 space = HCURL;
55 }
56
57 // Select base
58 enum bases { AINSWORTH, DEMKOWICZ, LASBASETOP };
59
60 const char *list_bases[] = {"ainsworth", "demkowicz"};
61
62 PetscInt choice_base_value = AINSWORTH;
63 CHKERR PetscOptionsGetEList(PETSC_NULL, NULL, "-base", list_bases,
64 LASBASETOP, &choice_base_value, &flg);
65
66 if (flg != PETSC_TRUE) {
67 SETERRQ(PETSC_COMM_SELF, MOFEM_IMPOSIBLE_CASE, "base not set");
68 }
69
71 if (choice_base_value == AINSWORTH) {
73 } else if (choice_base_value == DEMKOWICZ) {
75 }
76
77 moab::Core mb_instance;
78 moab::Interface &moab = mb_instance;
79 int rank;
80 MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
81
82 // create one tet
83 double tet_coords[] = {0, 0, 0, 0.5, 0, 0, 0, 0.5, 0, 0, 0, 0.5};
84 EntityHandle nodes[4];
85 for (int nn = 0; nn < 4; nn++) {
86 CHKERR moab.create_vertex(&tet_coords[3 * nn], nodes[nn]);
87 }
88 EntityHandle tet;
89 CHKERR moab.create_element(MBTET, nodes, 4, tet);
90
91 // Create MoFEM database
92 MoFEM::Core core(moab);
93 MoFEM::Interface &m_field = core;
94
95 // set entities bit level
96 BitRefLevel bit_level0;
97 bit_level0.set(0);
98 CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(
99 0, 3, bit_level0);
100
101
102 // Fields
103 CHKERR m_field.add_field("FIELD", space, base, 1);
104
105 // FE TET
106 CHKERR m_field.add_finite_element("TET_FE");
107
108 // Define rows/cols and element data
109 CHKERR m_field.modify_finite_element_add_field_row("TET_FE", "FIELD");
110 CHKERR m_field.modify_finite_element_add_field_col("TET_FE", "FIELD");
111 CHKERR m_field.modify_finite_element_add_field_data("TET_FE", "FIELD");
112
113
114 // Problem
115 CHKERR m_field.add_problem("TEST_PROBLEM");
116
117
118 // set finite elements for problem
119 CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM", "TET_FE");
120
121 // set refinement level for problem
122 CHKERR m_field.modify_problem_ref_level_add_bit("TEST_PROBLEM", bit_level0);
123
124
125 // meshset consisting all entities in mesh
126 EntityHandle root_set = moab.get_root_set();
127 // add entities to field
128 CHKERR m_field.add_ents_to_field_by_type(root_set, MBTET, "FIELD");
129
130 // add entities to finite element
131 CHKERR
132 m_field.add_ents_to_finite_element_by_type(root_set, MBTET, "TET_FE");
133
134
135 // set app. order
136 int order = 5;
137 if (space == H1) {
138 CHKERR m_field.set_field_order(root_set, MBTET, "FIELD", order);
139 CHKERR m_field.set_field_order(root_set, MBTRI, "FIELD", order);
140 CHKERR m_field.set_field_order(root_set, MBEDGE, "FIELD", order);
141 CHKERR m_field.set_field_order(root_set, MBVERTEX, "FIELD", 1);
142
143 }
144 if (space == HCURL) {
145 CHKERR m_field.set_field_order(root_set, MBTET, "FIELD", order);
146 CHKERR m_field.set_field_order(root_set, MBTRI, "FIELD", order);
147 CHKERR m_field.set_field_order(root_set, MBEDGE, "FIELD", order);
148
149 }
150 if (space == HDIV) {
151 CHKERR m_field.set_field_order(root_set, MBTET, "FIELD", order);
152 CHKERR m_field.set_field_order(root_set, MBTRI, "FIELD", order);
153 }
154
155 /****/
156 // build database
157 // build field
158 CHKERR m_field.build_fields();
159
160 // build finite elemnts
162
163 // build adjacencies
164 CHKERR m_field.build_adjacencies(bit_level0);
165
166
167 /****/
168 ProblemsManager *prb_mng_ptr;
169 CHKERR m_field.getInterface(prb_mng_ptr);
170
171 // build problem
172 CHKERR prb_mng_ptr->buildProblem("TEST_PROBLEM", true);
173
174 // mesh partitioning
175 // partition
176 CHKERR prb_mng_ptr->partitionSimpleProblem("TEST_PROBLEM");
177 CHKERR prb_mng_ptr->partitionFiniteElements("TEST_PROBLEM");
178
179 // what are ghost nodes, see Petsc Manual
180 CHKERR prb_mng_ptr->partitionGhostDofs("TEST_PROBLEM");
181
182 typedef tee_device<std::ostream, std::ofstream> TeeDevice;
183 typedef stream<TeeDevice> TeeStream;
184
185 std::ofstream ofs("forces_and_sources_checking_derivatives.txt");
186 TeeDevice my_tee(std::cout, ofs);
187 TeeStream my_split(my_tee);
188
189 struct OpCheckingDerivatives
191
192 TeeStream &mySplit;
193 OpCheckingDerivatives(TeeStream &my_split)
195 "FIELD", UserDataOperator::OPROW),
196 mySplit(my_split) {}
197
198 MoFEMErrorCode doWork(int side, EntityType type,
201
202 if (data.getFieldData().size() == 0)
204
205 mySplit << "Type " << type << " side " << side << endl;
206
207 if (data.getFieldDofs()[0]->getSpace() == H1) {
208
209 // mySplit << std::fixed << data.getN() << std::endl;
210 // mySplit << std::fixed << data.getDiffN() << std::endl;
211
212 const int nb_dofs = data.getN().size2();
213 for (int dd = 0; dd != nb_dofs; dd++) {
214 const double dksi = (data.getN()(1, dd) - data.getN()(0, dd)) / eps;
215 const double deta = (data.getN()(3, dd) - data.getN()(2, dd)) / eps;
216 const double dzeta =
217 (data.getN()(5, dd) - data.getN()(4, dd)) / eps;
218 mySplit << "DKsi " << dksi - data.getDiffN()(6, 3 * dd + 0)
219 << std::endl;
220 mySplit << "DEta " << deta - data.getDiffN()(6, 3 * dd + 1)
221 << std::endl;
222 mySplit << "DZeta " << dzeta - data.getDiffN()(6, 3 * dd + 2)
223 << std::endl;
224 if (fabs(dksi - data.getDiffN()(6, 3 * dd + 0)) > eps_diff) {
225 SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
226 "H1 inconsistent dKsi derivative");
227 }
228 if (fabs(deta - data.getDiffN()(6, 3 * dd + 1)) > eps_diff) {
229 SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
230 "H1 inconsistent dEta derivative");
231 }
232 if (fabs(dzeta - data.getDiffN()(6, 3 * dd + 2)) > eps_diff) {
233 SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
234 "H1 inconsistent dZeta derivative");
235 }
236 }
237 }
238
239 if (data.getFieldDofs()[0]->getSpace() == HDIV ||
240 data.getFieldDofs()[0]->getSpace() == HCURL) {
241
243 &data.getN()(0, 0), &data.getN()(0, 1), &data.getN()(0, 2), 3);
245 &data.getN()(1, 0), &data.getN()(1, 1), &data.getN()(1, 2), 3);
247 &data.getN()(2, 0), &data.getN()(2, 1), &data.getN()(2, 2), 3);
249 &data.getN()(3, 0), &data.getN()(3, 1), &data.getN()(3, 2), 3);
251 &data.getN()(4, 0), &data.getN()(4, 1), &data.getN()(4, 2), 3);
253 &data.getN()(5, 0), &data.getN()(5, 1), &data.getN()(5, 2), 3);
254
256 &data.getDiffN()(6, 0), &data.getDiffN()(6, 3),
257 &data.getDiffN()(6, 6), &data.getDiffN()(6, 1),
258 &data.getDiffN()(6, 4), &data.getDiffN()(6, 7),
259 &data.getDiffN()(6, 2), &data.getDiffN()(6, 5),
260 &data.getDiffN()(6, 8), 9);
261
266
267 const int nb_dofs = data.getN().size2() / 3;
268 for (int dd = 0; dd != nb_dofs; dd++) {
270 dksi(i) = (base_ksi_p(i) - base_ksi_m(i)) / eps;
272 deta(i) = (base_eta_p(i) - base_eta_m(i)) / eps;
274 dzeta(i) = (base_zeta_p(i) - base_zeta_m(i)) / eps;
275
276 dksi(i) -= diff_base(i, N0);
277 deta(i) -= diff_base(i, N1);
278 dzeta(i) -= diff_base(i, N2);
279
280 mySplit << "dKsi " << dksi(0) << " " << dksi(1) << " " << dksi(2)
281 << " " << sqrt(dksi(i) * dksi(i)) << endl;
282 mySplit << "dEta " << deta(0) << " " << deta(1) << " " << deta(2)
283 << " " << sqrt(deta(i) * deta(i)) << endl;
284 mySplit << "dZeta " << dzeta(0) << " " << dzeta(1) << " "
285 << dzeta(2) << " " << sqrt(dzeta(i) * dzeta(i)) << endl;
286
287 if (sqrt(dksi(i) * dksi(i)) > eps_diff) {
288 SETERRQ2(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
289 "%s inconsistent dKsi derivative for type %d",
290 FieldSpaceNames[data.getFieldDofs()[0]->getSpace()],
291 type);
292 }
293 if (sqrt(deta(i) * deta(i)) > eps_diff) {
294 SETERRQ2(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
295 "%s inconsistent dEta derivative for type %d",
296 FieldSpaceNames[data.getFieldDofs()[0]->getSpace()],
297 type);
298 }
299 if (sqrt(dzeta(i) * dzeta(i)) > eps_diff) {
300 SETERRQ2(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
301 "%s inconsistent dZeta derivative for type %d",
302 FieldSpaceNames[data.getFieldDofs()[0]->getSpace()],
303 type);
304 }
305
306 ++base_ksi_m;
307 ++base_ksi_p;
308 ++base_eta_m;
309 ++base_eta_p;
310 ++base_zeta_m;
311 ++base_zeta_p;
312 ++diff_base;
313 }
314 }
315
317 }
318 };
319
320 struct MyFE : public VolumeElementForcesAndSourcesCore {
321
322 MyFE(MoFEM::Interface &m_field)
324 int getRule(int order) { return -1; };
325
328
329 const double ksi = G_TET_X1[0];
330 const double eta = G_TET_Y1[0];
331 const double zeta = G_TET_Z1[0];
332
333 gaussPts.resize(4, 7);
334 gaussPts.clear();
335
336 gaussPts(0, 0) = ksi - eps;
337 gaussPts(1, 0) = eta;
338 gaussPts(2, 0) = zeta;
339 gaussPts(0, 1) = ksi + eps;
340 gaussPts(1, 1) = eta;
341 gaussPts(2, 1) = zeta;
342
343 gaussPts(0, 2) = ksi;
344 gaussPts(1, 2) = eta - eps;
345 gaussPts(2, 2) = zeta;
346 gaussPts(0, 3) = ksi;
347 gaussPts(1, 3) = eta + eps;
348 gaussPts(2, 3) = zeta;
349
350 gaussPts(0, 4) = ksi;
351 gaussPts(1, 4) = eta;
352 gaussPts(2, 4) = zeta - eps;
353 gaussPts(0, 5) = ksi;
354 gaussPts(1, 5) = eta;
355 gaussPts(2, 5) = zeta + eps;
356
357 gaussPts(0, 6) = ksi;
358 gaussPts(1, 6) = eta;
359 gaussPts(2, 6) = zeta;
360
361 for (unsigned int ii = 0; ii != gaussPts.size2(); ii++) {
362 gaussPts(3, ii) = 1;
363 }
364
365 cerr << gaussPts << endl;
366
368 }
369 };
370
371 MyFE tet_fe(m_field);
372
373 tet_fe.getOpPtrVector().push_back(new OpCheckingDerivatives(my_split));
374 CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "TET_FE", tet_fe);
375
376 }
378
380
381}
static Number< 2 > N2
static Number< 1 > N1
static Number< 0 > N0
ForcesAndSourcesCore::UserDataOperator UserDataOperator
static const double eps
static const double eps_diff
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:385
FieldApproximationBase
approximation base
Definition: definitions.h:71
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:73
@ NOBASE
Definition: definitions.h:72
@ DEMKOWICZ_JACOBI_BASE
Definition: definitions.h:79
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:460
FieldSpace
approximation spaces
Definition: definitions.h:95
@ LASTSPACE
FieldSpace in [ 0, LASTSPACE )
Definition: definitions.h:102
@ H1
continuous field
Definition: definitions.h:98
@ HCURL
field with continuous tangents
Definition: definitions.h:99
@ HDIV
field with continuous normal traction
Definition: definitions.h:100
static const char *const FieldSpaceNames[]
Definition: definitions.h:105
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
@ MOFEM_IMPOSIBLE_CASE
Definition: definitions.h:48
@ MOFEM_ATOM_TEST_INVALID
Definition: definitions.h:53
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
#define CHKERR
Inline error check.
Definition: definitions.h:548
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:453
static const double G_TET_X1[]
Definition: fem_tools.h:1122
static const double G_TET_Y1[]
Definition: fem_tools.h:1123
static const double G_TET_Z1[]
Definition: fem_tools.h:1124
tee_device< std::ostream, std::ofstream > TeeDevice
constexpr double eta
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 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 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_col(const std::string &fe_name, const std::string &name_row)=0
set field col 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 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
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
FTensor::Index< 'i', SPACE_DIM > i
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:67
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition: Types.hpp:51
CoreTmp< 0 > Core
Definition: Core.hpp:1096
PetscErrorCode PetscOptionsGetEList(PetscOptions *, const char pre[], const char name[], const char *const *list, PetscInt next, PetscInt *value, PetscBool *set)
DeprecatedCoreInterface Interface
Definition: Interface.hpp:1965
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:92
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
Definition: Core.cpp:85
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:125
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
const VectorDofs & getFieldDofs() const
get dofs data stature FEDofEntity
virtual MoFEMErrorCode setGaussPts(int order_row, int order_col, int order_data)
set user specific integration rule
virtual int getRule(int order_row, int order_col, int order_data)
another variant of getRule
MatrixDouble gaussPts
Matrix of integration points.
Problem manager is used to build and partition problems.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.

Variable Documentation

◆ eps

const double eps = 1e-6
static

◆ eps_diff

const double eps_diff = 1e-8
static

Definition at line 26 of file check_base_functions_derivatives_on_tet.cpp.

◆ help

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

Definition at line 23 of file check_base_functions_derivatives_on_tet.cpp.