|
| v0.14.0
|
Go to the documentation of this file.
3 namespace bio = boost::iostreams;
9 static char help[] =
"...\n\n";
11 static const double eps = 1e-6;
14 int main(
int argc,
char *argv[]) {
21 enum spaces { H1TET, HDIVTET, HCURLTET, LASTSPACEOP };
23 const char *list_spaces[] = {
"h1tet",
"hdivtet",
"hcurltet"};
26 PetscInt choise_space_value = H1TET;
28 LASTSPACEOP, &choise_space_value, &flg);
30 if (flg != PETSC_TRUE) {
35 if (choise_space_value == H1TET) {
37 }
else if (choise_space_value == HDIVTET) {
39 }
else if (choise_space_value == HCURLTET) {
44 enum bases { AINSWORTH, DEMKOWICZ, LASBASETOP };
46 const char *list_bases[] = {
"ainsworth",
"demkowicz"};
48 PetscInt choice_base_value = AINSWORTH;
50 LASBASETOP, &choice_base_value, &flg);
52 if (flg != PETSC_TRUE) {
57 if (choice_base_value == AINSWORTH) {
59 }
else if (choice_base_value == DEMKOWICZ) {
66 MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
69 double tet_coords[] = {0, 0, 0, 0.5, 0, 0, 0, 0.5, 0, 0, 0, 0.5};
71 for (
int nn = 0; nn < 4; nn++) {
72 CHKERR moab.create_vertex(&tet_coords[3 * nn], nodes[nn]);
75 CHKERR moab.create_element(MBTET, nodes, 4, tet);
130 if (space ==
HCURL) {
168 typedef tee_device<std::ostream, std::ofstream>
TeeDevice;
171 std::ofstream ofs(
"forces_and_sources_checking_derivatives.txt");
175 struct OpCheckingDerivatives
179 OpCheckingDerivatives(
TeeStream &my_split)
191 mySplit <<
"Type " <<
type <<
" side " << side << endl;
198 const int nb_dofs = data.
getN().size2();
199 for (
int dd = 0;
dd != nb_dofs;
dd++) {
204 mySplit <<
"DKsi " << dksi - data.
getDiffN()(6, 3 *
dd + 0)
206 mySplit <<
"DEta " << deta - data.
getDiffN()(6, 3 *
dd + 1)
208 mySplit <<
"DZeta " << dzeta - data.
getDiffN()(6, 3 *
dd + 2)
212 "H1 inconsistent dKsi derivative");
216 "H1 inconsistent dEta derivative");
220 "H1 inconsistent dZeta derivative");
229 &data.
getN()(0, 0), &data.
getN()(0, 1), &data.
getN()(0, 2), 3);
231 &data.
getN()(1, 0), &data.
getN()(1, 1), &data.
getN()(1, 2), 3);
233 &data.
getN()(2, 0), &data.
getN()(2, 1), &data.
getN()(2, 2), 3);
235 &data.
getN()(3, 0), &data.
getN()(3, 1), &data.
getN()(3, 2), 3);
237 &data.
getN()(4, 0), &data.
getN()(4, 1), &data.
getN()(4, 2), 3);
239 &data.
getN()(5, 0), &data.
getN()(5, 1), &data.
getN()(5, 2), 3);
253 const int nb_dofs = data.
getN().size2() / 3;
254 for (
int dd = 0;
dd != nb_dofs;
dd++) {
256 dksi(
i) = (base_ksi_p(
i) - base_ksi_m(
i)) /
eps;
258 deta(
i) = (base_eta_p(
i) - base_eta_m(
i)) /
eps;
260 dzeta(
i) = (base_zeta_p(
i) - base_zeta_m(
i)) /
eps;
262 dksi(
i) -= diff_base(
i, N0);
263 deta(
i) -= diff_base(
i, N1);
264 dzeta(
i) -= diff_base(
i, N2);
266 mySplit <<
"dKsi " << dksi(0) <<
" " << dksi(1) <<
" " << dksi(2)
267 <<
" " << sqrt(dksi(
i) * dksi(
i)) << endl;
268 mySplit <<
"dEta " << deta(0) <<
" " << deta(1) <<
" " << deta(2)
269 <<
" " << sqrt(deta(
i) * deta(
i)) << endl;
270 mySplit <<
"dZeta " << dzeta(0) <<
" " << dzeta(1) <<
" "
271 << dzeta(2) <<
" " << sqrt(dzeta(
i) * dzeta(
i)) << endl;
275 "%s inconsistent dKsi derivative for type %d",
281 "%s inconsistent dEta derivative for type %d",
287 "%s inconsistent dZeta derivative for type %d",
310 int getRule(
int order) {
return -1; };
319 gaussPts.resize(4, 7);
322 gaussPts(0, 0) = ksi -
eps;
323 gaussPts(1, 0) =
eta;
324 gaussPts(2, 0) =
zeta;
325 gaussPts(0, 1) = ksi +
eps;
326 gaussPts(1, 1) =
eta;
327 gaussPts(2, 1) =
zeta;
329 gaussPts(0, 2) = ksi;
330 gaussPts(1, 2) =
eta -
eps;
331 gaussPts(2, 2) =
zeta;
332 gaussPts(0, 3) = ksi;
333 gaussPts(1, 3) =
eta +
eps;
334 gaussPts(2, 3) =
zeta;
336 gaussPts(0, 4) = ksi;
337 gaussPts(1, 4) =
eta;
339 gaussPts(0, 5) = ksi;
340 gaussPts(1, 5) =
eta;
343 gaussPts(0, 6) = ksi;
344 gaussPts(1, 6) =
eta;
345 gaussPts(2, 6) =
zeta;
347 for (
unsigned int ii = 0; ii != gaussPts.size2(); ii++) {
351 cerr << gaussPts << endl;
357 MyFE tet_fe(m_field);
359 tet_fe.getOpPtrVector().push_back(
new OpCheckingDerivatives(my_split));
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Data on single entity (This is passed as argument to DataOperator::doWork)
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 buildProblem(const std::string name, const bool square_matrix, int verb=VERBOSE)
build problem data structures
Problem manager is used to build and partition problems.
static const double eps_diff
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
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
const VectorDouble & getFieldData() const
get dofs values
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
MatrixDouble & getDiffN(const FieldApproximationBase base)
get derivatives of base functions
int main(int argc, char *argv[])
double zeta
Viscous hardening.
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.
Deprecated interface functions.
DeprecatedCoreInterface Interface
FieldSpace
approximation spaces
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
#define CHKERR
Inline error check.
virtual MoFEMErrorCode add_finite_element(const std::string &fe_name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
add finite element
implementation of Data Operators for Forces and Sources
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_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elements.
friend class UserDataOperator
@ LASTSPACE
FieldSpace in [ 0, LASTSPACE )
const VectorDofs & getFieldDofs() const
get dofs data stature FEDofEntity
virtual MoFEMErrorCode modify_finite_element_add_field_data(const std::string &fe_name, const std::string name_field)=0
set finite element field data
Volume finite element base.
FTensor::Index< 'i', SPACE_DIM > i
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_ref_level_add_bit(const std::string &name_problem, const BitRefLevel &bit)=0
add ref level to problem
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
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)
const static char *const FieldSpaceNames[]
#define CATCH_ERRORS
Catch errors.
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
MoFEMErrorCode partitionSimpleProblem(const std::string name, int verb=VERBOSE)
partition problem dofs
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
@ HCURL
field with continuous tangents
PetscErrorCode PetscOptionsGetEList(PetscOptions *, const char pre[], const char name[], const char *const *list, PetscInt next, PetscInt *value, PetscBool *set)
virtual MoFEMErrorCode build_fields(int verb=DEFAULT_VERBOSITY)=0
FieldApproximationBase
approximation base
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
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
@ MOFEM_ATOM_TEST_INVALID
virtual MoFEMErrorCode build_adjacencies(const Range &ents, int verb=DEFAULT_VERBOSITY)=0
build adjacencies
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
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.
MoFEMErrorCode partitionGhostDofs(const std::string name, int verb=VERBOSE)
determine ghost nodes
virtual MoFEMErrorCode add_problem(const std::string &name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add problem.
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
@ HDIV
field with continuous normal traction
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.
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ OPROW
operator doWork function is executed on FE rows