v0.14.0
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 
12 namespace bio = boost::iostreams;
13 using bio::stream;
14 using bio::tee_device;
15 
16 using namespace MoFEM;
17 
18 static char help[] = "...\n\n";
19 
20 int 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
99  CHKERR m_field.build_finite_elements();
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 
163  MoFEMErrorCode preProcess() {
166  }
167 
168  EntitiesFieldData data;
169 
170  MoFEMErrorCode operator()() {
172 
173  CHKERR getSpacesAndBaseOnEntities(data);
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>(new EntPolynomialBaseCtx(
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 
226  MoFEMErrorCode postProcess() {
229  }
230  };
231 
232  ForcesAndSourcesCore_TestFE fe1(m_field);
233  CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "TEST_FE", fe1);
234  }
235  CATCH_ERRORS;
236 
238 
239  return 0;
240 }
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:460
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Definition: UnknownInterface.hpp:93
MoFEM::EntitiesFieldData::EntData
Data on single entity (This is passed as argument to DataOperator::doWork)
Definition: EntitiesFieldData.hpp:128
TeeStream
stream< TeeDevice > TeeStream
Definition: forces_and_sources_testing_contact_prism_element.cpp:10
MoFEM::CoreInterface::loop_finite_elements
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.
MoFEM::ProblemsManager::buildProblem
MoFEMErrorCode buildProblem(const std::string name, const bool square_matrix, int verb=VERBOSE)
build problem data structures
Definition: ProblemsManager.cpp:87
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
H1
@ H1
continuous field
Definition: definitions.h:85
MoFEM::DataOperator
base operator to do operations at Gauss Pt. level
Definition: DataOperators.hpp:24
MoFEM::DataOperator::opRhs
virtual MoFEMErrorCode opRhs(EntitiesFieldData &data, const bool error_if_no_base=false)
Definition: DataOperators.cpp:140
EntityHandle
MoFEM::ProblemsManager
Problem manager is used to build and partition problems.
Definition: ProblemsManager.hpp:21
MoFEM::EntPolynomialBaseCtx
Class used to pass element data to calculate base functions on tet,triangle,edge.
Definition: EntPolynomialBaseCtx.hpp:22
G_TET_Z4
static const double G_TET_Z4[]
Definition: fem_tools.h:1120
MoFEM::CoreInterface::modify_finite_element_add_field_row
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
MoFEM::OpSetInvJacH1
Transform local reference derivatives of shape function to global derivatives.
Definition: DataOperators.hpp:131
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
help
static char help[]
Definition: forces_and_sources_calculate_jacobian.cpp:18
MoFEM::Types::MatrixDouble
UBlasMatrix< double > MatrixDouble
Definition: Types.hpp:77
G_TET_Y4
static const double G_TET_Y4[]
Definition: fem_tools.h:1118
MoFEM.hpp
MoFEM::CoreTmp< 0 >::Finalize
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:112
MoFEM::EntitiesFieldData::EntData::getDiffN
MatrixDouble & getDiffN(const FieldApproximationBase base)
get derivatives of base functions
Definition: EntitiesFieldData.hpp:1329
MoFEM::CoreInterface::add_ents_to_field_by_type
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.
order
constexpr int order
Definition: dg_projection.cpp:18
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
TeeDevice
tee_device< std::ostream, std::ofstream > TeeDevice
Definition: forces_and_sources_testing_contact_prism_element.cpp:9
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:2010
MoFEM::CoreInterface::add_ents_to_finite_element_by_type
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
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
MoFEM::CoreInterface::add_finite_element
virtual MoFEMErrorCode add_finite_element(const std::string &fe_name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
add finite element
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
MoFEM::CoreInterface::modify_finite_element_add_field_col
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
MoFEM::CoreInterface::build_finite_elements
virtual MoFEMErrorCode build_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elements.
convert.type
type
Definition: convert.py:64
MoFEM::CoreInterface::modify_finite_element_add_field_data
virtual MoFEMErrorCode modify_finite_element_add_field_data(const std::string &fe_name, const std::string name_field)=0
set finite element field data
MoFEM::ProblemsManager::partitionFiniteElements
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
Definition: ProblemsManager.cpp:2167
MoFEM::CoreInterface::modify_problem_ref_level_add_bit
virtual MoFEMErrorCode modify_problem_ref_level_add_bit(const std::string &name_problem, const BitRefLevel &bit)=0
add ref level to problem
mesh_file_name
char mesh_file_name[255]
Definition: mesh_smoothing.cpp:22
MoFEM::ForcesAndSourcesCore
structure to get information form mofem into EntitiesFieldData
Definition: ForcesAndSourcesCore.hpp:22
G_TET_X4
static const double G_TET_X4[]
Definition: fem_tools.h:1116
G_TET_W4
static const double G_TET_W4[]
Definition: fem_tools.h:1122
MoFEM::CoreTmp< 0 >::Initialize
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
ShapeJacMBTET
PetscErrorCode ShapeJacMBTET(double *diffN, const double *coords, double *jac)
calculate jacobian
Definition: fem_tools.c:288
main
int main(int argc, char *argv[])
Definition: forces_and_sources_calculate_jacobian.cpp:20
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:385
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1148
MoFEM::ProblemsManager::partitionSimpleProblem
MoFEMErrorCode partitionSimpleProblem(const std::string name, int verb=VERBOSE)
partition problem dofs
Definition: ProblemsManager.cpp:1537
eps
static const double eps
Definition: check_base_functions_derivatives_on_tet.cpp:11
AINSWORTH_LEGENDRE_BASE
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
MoFEM::PetscOptionsGetString
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
Definition: DeprecatedPetsc.hpp:172
MoFEM::CoreInterface::build_fields
virtual MoFEMErrorCode build_fields(int verb=DEFAULT_VERBOSITY)=0
MoFEM::TetPolynomialBase::getValue
MoFEMErrorCode getValue(MatrixDouble &pts, boost::shared_ptr< BaseFunctionCtx > ctx_ptr)
Definition: TetPolynomialBase.cpp:2021
MoFEM::Types::MatrixDouble3by3
MatrixBoundedArray< double, 9 > MatrixDouble3by3
Definition: Types.hpp:105
MoFEM::Types::VectorDouble
UBlasVector< double > VectorDouble
Definition: Types.hpp:68
MoFEM::EntitiesFieldData::dataOnEntities
std::array< boost::ptr_vector< EntData >, MBMAXTYPE > dataOnEntities
Definition: EntitiesFieldData.hpp:57
MoFEM::BitRefManager
Managing BitRefLevels.
Definition: BitRefManager.hpp:21
MoFEM::CoreInterface::modify_problem_add_finite_element
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
MoFEM::Types::BitRefLevel
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition: Types.hpp:40
MoFEM::CoreInterface::build_adjacencies
virtual MoFEMErrorCode build_adjacencies(const Range &ents, int verb=DEFAULT_VERBOSITY)=0
build adjacencies
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:453
MoFEM::CoreInterface::set_field_order
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.
ShapeInvJacVolume
PetscErrorCode ShapeInvJacVolume(double *jac)
Definition: fem_tools.c:39
MoFEM::EntitiesFieldData
data structure for finite element entity
Definition: EntitiesFieldData.hpp:40
invJac
MatrixDouble invJac
Definition: HookeElement.hpp:683
MoFEM::ProblemsManager::partitionGhostDofs
MoFEMErrorCode partitionGhostDofs(const std::string name, int verb=VERBOSE)
determine ghost nodes
Definition: ProblemsManager.cpp:2339
CONTINUOUS
@ CONTINUOUS
Regular field.
Definition: definitions.h:100
MoFEM::CoreInterface::add_problem
virtual MoFEMErrorCode add_problem(const std::string &name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add problem.
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
MoFEM::OpGetDataAndGradient
Get field values and gradients at Gauss points.
Definition: DataOperators.hpp:240
MoFEM::CoreInterface::add_field
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.
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
MoFEM::TetPolynomialBase
Calculate base functions on tetrahedral.
Definition: TetPolynomialBase.hpp:17