v0.13.2
Loading...
Searching...
No Matches
forces_and_sources_calculate_jacobian.cpp
Go to the documentation of this file.
1/** \file forces_and_sources_calculate_jacobian.cpp
2
3 \brief Atom test checking with blessed file how Jacobian's are calculated on
4 elements
5
6*/
7
8
9
10#include <MoFEM.hpp>
11
12namespace bio = boost::iostreams;
13using bio::stream;
14using bio::tee_device;
15
16using namespace MoFEM;
17
18static char help[] = "...\n\n";
19
20int 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 CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(
56 0, 3, bit_level0);
57
58 // Fields
59 CHKERR m_field.add_field("FIELD1", H1, AINSWORTH_LEGENDRE_BASE, 1);
60
61 // FE
62 CHKERR m_field.add_finite_element("TEST_FE");
63
64 // Define rows/cols and element data
65 CHKERR m_field.modify_finite_element_add_field_row("TEST_FE", "FIELD1");
66 CHKERR m_field.modify_finite_element_add_field_col("TEST_FE", "FIELD1");
67 CHKERR m_field.modify_finite_element_add_field_data("TEST_FE", "FIELD1");
68
69 // Problem
70 CHKERR m_field.add_problem("TEST_PROBLEM");
71
72 // set finite elements for problem
73 CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM", "TEST_FE");
74 // set refinement level for problem
75 CHKERR m_field.modify_problem_ref_level_add_bit("TEST_PROBLEM", bit_level0);
76
77 // meshset consisting all entities in mesh
78 EntityHandle root_set = moab.get_root_set();
79 // add entities to field
80 CHKERR m_field.add_ents_to_field_by_type(root_set, MBTET, "FIELD1");
81 // add entities to finite element
82 CHKERR m_field.add_ents_to_finite_element_by_type(root_set, MBTET,
83 "TEST_FE");
84
85 // set app. order
86 // see Hierarchic Finite Element Bases on Unstructured Tetrahedral Meshes
87 // (Mark Ainsworth & Joe Coyle)
88 int order = 5;
89 CHKERR m_field.set_field_order(root_set, MBTET, "FIELD1", order);
90 CHKERR m_field.set_field_order(root_set, MBTRI, "FIELD1", order);
91 CHKERR m_field.set_field_order(root_set, MBEDGE, "FIELD1", order);
92 CHKERR m_field.set_field_order(root_set, MBVERTEX, "FIELD1", 1);
93
94 /****/
95 // build database
96 // build field
97 CHKERR m_field.build_fields();
98 // build finite elemnts
100 // build adjacencies
101 CHKERR m_field.build_adjacencies(bit_level0);
102 // build problem
103 ProblemsManager *prb_mng_ptr;
104 CHKERR m_field.getInterface(prb_mng_ptr);
105 CHKERR prb_mng_ptr->buildProblem("TEST_PROBLEM", true);
106
107 /****/
108 // mesh partitioning
109 // partition
110 CHKERR prb_mng_ptr->partitionSimpleProblem("TEST_PROBLEM");
111 CHKERR prb_mng_ptr->partitionFiniteElements("TEST_PROBLEM");
112 // what are ghost nodes, see Petsc Manual
113 CHKERR prb_mng_ptr->partitionGhostDofs("TEST_PROBLEM");
114
115 struct ForcesAndSourcesCore_TestFE : public ForcesAndSourcesCore {
116
117 typedef tee_device<std::ostream, std::ofstream> TeeDevice;
118 typedef stream<TeeDevice> TeeStream;
119
120 std::ofstream ofs;
121 TeeDevice my_tee;
122 TeeStream my_split;
123
124 struct PrintJacobian : public DataOperator {
125
126 TeeStream &my_split;
127 PrintJacobian(TeeStream &_my_split) : my_split(_my_split) {}
128
129 ~PrintJacobian() { my_split.close(); }
130
131 MoFEMErrorCode doWork(int side, EntityType type,
134 const double eps = 1e-6;
135 for (unsigned int ii = 0; ii != data.getDiffN().size1(); ii++) {
136 for (unsigned int jj = 0; jj != data.getDiffN().size2(); jj++) {
137 if (fabs(data.getDiffN()(ii, jj)) < eps) {
138 data.getDiffN()(ii, jj) = 0;
139 }
140 }
141 }
142 my_split << "side: " << side << " type: " << type << std::fixed
143 << std::setprecision(4) << data.getDiffN() << std::endl;
145 }
146 };
147
149 MatrixDouble dataFIELD1;
150 MatrixDouble dataDiffFIELD1;
151 VectorDouble coords;
152 PrintJacobian opPrintJac;
153 OpSetInvJacH1 opSetInvJac;
154 OpGetDataAndGradient<1, 3> opGetData_FIELD1;
155
156 ForcesAndSourcesCore_TestFE(MoFEM::Interface &_m_field)
157 : ForcesAndSourcesCore(_m_field),
158 ofs("forces_and_sources_calculate_jacobian.txt"),
159 my_tee(std::cout, ofs), my_split(my_tee), invJac(3, 3),
160 opPrintJac(my_split), opSetInvJac(invJac),
161 opGetData_FIELD1(dataFIELD1, dataDiffFIELD1), data(MBTET) {}
162
166 }
167
169
172
174
175 CHKERR getEntitySense<MBEDGE>(data);
176 CHKERR getEntitySense<MBTRI>(data);
177 CHKERR getEntityDataOrder<MBEDGE>(data, H1);
178 CHKERR getEntityDataOrder<MBTRI>(data, H1);
179 CHKERR getEntityDataOrder<MBTET>(data, H1);
180 CHKERR getFaceNodes(data);
181
182 const auto bit_number = mField.get_field_bit_number("FIELD1");
183 CHKERR getRowNodesIndices(data, bit_number);
184 CHKERR getEntityRowIndices(data, bit_number, MBEDGE);
185
186 CHKERR getNodesFieldData(data, bit_number);
187 CHKERR getEntityFieldData(data, bit_number, MBEDGE);
188
189 MatrixDouble gauss_pts(4, 4);
190 for (int gg = 0; gg < 4; gg++) {
191 gauss_pts(0, gg) = G_TET_X4[gg];
192 gauss_pts(1, gg) = G_TET_Y4[gg];
193 gauss_pts(2, gg) = G_TET_Z4[gg];
194 gauss_pts(3, gg) = G_TET_W4[gg];
195 }
197 gauss_pts,
198 boost::shared_ptr<BaseFunctionCtx>(
200
201 EntityHandle ent = numeredEntFiniteElementPtr->getEnt();
202 int num_nodes;
203 const EntityHandle *conn;
204 CHKERR mField.get_moab().get_connectivity(ent, conn, num_nodes, true);
205 coords.resize(num_nodes * 3);
206 CHKERR mField.get_moab().get_coords(conn, num_nodes,
207 &*coords.data().begin());
208
210 &*data.dataOnEntities[MBVERTEX][0].getDiffN().data().begin(),
211 &*coords.begin(), &*invJac.data().begin());
212 CHKERR ShapeInvJacVolume(&*invJac.data().begin());
213
214 CHKERR opSetInvJac.opRhs(data);
215 CHKERR opPrintJac.opRhs(data);
216 CHKERR opGetData_FIELD1.opRhs(data);
217
218 my_split << "data FIELD1:" << std::endl;
219 my_split << dataFIELD1 << std::endl;
220 my_split << "data diff FIELD1:" << std::endl;
221 my_split << dataDiffFIELD1 << std::endl;
222
224 }
225
229 }
230 };
231
232 ForcesAndSourcesCore_TestFE fe1(m_field);
233 CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "TEST_FE", fe1);
234 }
236
238
239 return 0;
240}
MatrixDouble invJac
int main()
Definition: adol-c_atom.cpp:46
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
@ 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
PetscErrorCode ShapeJacMBTET(double *diffN, const double *coords, double *jac)
calculate jacobian
Definition: fem_tools.c:288
static const double G_TET_W4[]
Definition: fem_tools.h:1122
PetscErrorCode ShapeInvJacVolume(double *jac)
Definition: fem_tools.c:39
static const double G_TET_Z4[]
Definition: fem_tools.h:1120
static const double G_TET_Y4[]
Definition: fem_tools.h:1118
static const double G_TET_X4[]
Definition: fem_tools.h:1116
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
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
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
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
base operator to do operations at Gauss Pt. level
virtual MoFEMErrorCode opRhs(EntitiesFieldData &data, const bool error_if_no_base=false)
Deprecated interface functions.
Class used to pass element data to calculate base functions on tet,triangle,edge.
Data on single entity (This is passed as argument to DataOperator::doWork)
MatrixDouble & getDiffN(const FieldApproximationBase base)
get derivatives of base functions
data structure for finite element entity
std::array< boost::ptr_vector< EntData >, MBMAXTYPE > dataOnEntities
structure to get information form mofem into EntitiesFieldData
MoFEMErrorCode getEntityRowIndices(EntitiesFieldData &data, const int bit_number, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const
MoFEMErrorCode getRowNodesIndices(EntitiesFieldData &data, const int bit_number) const
get row node indices from FENumeredDofEntity_multiIndex
MoFEMErrorCode getSpacesAndBaseOnEntities(EntitiesFieldData &data) const
Get field approximation space and base on entities.
virtual MoFEMErrorCode operator()()
function is run for every finite element
MoFEMErrorCode getFaceNodes(EntitiesFieldData &data) const
Get nodes on faces.
MoFEMErrorCode getEntityFieldData(EntitiesFieldData &data, const int bit_number, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const
virtual MoFEMErrorCode preProcess()
function is run at the beginning of loop
MoFEMErrorCode getNodesFieldData(EntitiesFieldData &data, const int bit_number) const
Get data on nodes.
virtual MoFEMErrorCode postProcess()
function is run at the end of loop
Get field values and gradients at Gauss points.
Transform local reference derivatives of shape function to global derivatives.
Problem manager is used to build and partition problems.
Calculate base functions on tetrahedral.
MoFEMErrorCode getValue(MatrixDouble &pts, boost::shared_ptr< BaseFunctionCtx > ctx_ptr)
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.