v0.14.0
forces_and_sources_testing_ftensor.cpp
Go to the documentation of this file.
1 
2 
3 #define FTENSOR_DEBUG
4 
5 #include <MoFEM.hpp>
6 using namespace MoFEM;
7 
8 namespace bio = boost::iostreams;
9 using bio::tee_device;
10 using bio::stream;
11 
12 
13 static char help[] = "...\n\n";
14 
15 int main(int argc, char *argv[]) {
16 
17 
18 
19 
20  MoFEM::Core::Initialize(&argc,&argv,(char *)0,help);
21 
22  try {
23 
24  moab::Core mb_instance;
25  moab::Interface& moab = mb_instance;
26  int rank;
27  MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
28 
29  PetscBool flg = PETSC_TRUE;
30  char mesh_file_name[255];
31  #if PETSC_VERSION_GE(3,6,4)
32  ierr = PetscOptionsGetString(PETSC_NULL,"","-my_file",mesh_file_name,255,&flg); CHKERRG(ierr);
33  #else
34  ierr = PetscOptionsGetString(PETSC_NULL,PETSC_NULL,"-my_file",mesh_file_name,255,&flg); CHKERRG(ierr);
35  #endif
36  if(flg != PETSC_TRUE) {
37  SETERRQ(PETSC_COMM_SELF,MOFEM_DATA_INCONSISTENCY,"*** ERROR -my_file (MESH FILE NEEDED)");
38  }
39 
40  //Create MoFEM (Joseph) database
41  MoFEM::Core core(moab);
42  MoFEM::Interface& m_field = core;
43 
44  const char *option;
45  option = "";
46  rval = moab.load_file(mesh_file_name, 0, option); CHKERRG(rval);
47 
48  //set entitities bit level
49  BitRefLevel bit_level0;
50  bit_level0.set(0);
51  EntityHandle meshset_level0;
52  rval = moab.create_meshset(MESHSET_SET,meshset_level0); CHKERRG(rval);
53  ierr = m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(0,3,bit_level0); CHKERRG(ierr);
54 
55  //Fields
56  ierr = m_field.add_field("FIELD1",H1,AINSWORTH_LEGENDRE_BASE,3); CHKERRG(ierr);
57  ierr = m_field.add_field("FIELD2",H1,AINSWORTH_LEGENDRE_BASE,1); CHKERRG(ierr);
58 
59 
60  //FE
61  ierr = m_field.add_finite_element("TEST_FE1"); CHKERRG(ierr);
62 
63  //Define rows/cols and element data
64  ierr = m_field.modify_finite_element_add_field_row("TEST_FE1","FIELD1"); CHKERRG(ierr);
65  ierr = m_field.modify_finite_element_add_field_col("TEST_FE1","FIELD1"); CHKERRG(ierr);
66  ierr = m_field.modify_finite_element_add_field_data("TEST_FE1","FIELD1"); CHKERRG(ierr);
67  ierr = m_field.modify_finite_element_add_field_data("TEST_FE1","FIELD2"); CHKERRG(ierr);
68 
69  //Problem
70  ierr = m_field.add_problem("TEST_PROBLEM"); CHKERRG(ierr);
71 
72  //set finite elements for problem
73  ierr = m_field.modify_problem_add_finite_element("TEST_PROBLEM","TEST_FE1"); CHKERRG(ierr);
74  //set refinement level for problem
75  ierr = m_field.modify_problem_ref_level_add_bit("TEST_PROBLEM",bit_level0); CHKERRG(ierr);
76 
77 
78  //meshset consisting all entities in mesh
79  EntityHandle root_set = moab.get_root_set();
80  //add entities to field
81  ierr = m_field.add_ents_to_field_by_type(root_set,MBTET,"FIELD1"); CHKERRG(ierr);
82  ierr = m_field.add_ents_to_field_by_type(root_set,MBTET,"FIELD2"); CHKERRG(ierr);
83 
84 
85  //add entities to finite element
86  ierr = m_field.add_ents_to_finite_element_by_type(root_set,MBTET,"TEST_FE1"); CHKERRG(ierr);
87 
88  //set app. order
89  //see Hierarchic Finite Element Bases on Unstructured Tetrahedral Meshes (Mark Ainsworth & Joe Coyle)
90  int order = 4;
91  ierr = m_field.set_field_order(root_set,MBTET,"FIELD1",order); CHKERRG(ierr);
92  ierr = m_field.set_field_order(root_set,MBTRI,"FIELD1",order); CHKERRG(ierr);
93  ierr = m_field.set_field_order(root_set,MBEDGE,"FIELD1",order); CHKERRG(ierr);
94  ierr = m_field.set_field_order(root_set,MBVERTEX,"FIELD1",1); CHKERRG(ierr);
95 
96  ierr = m_field.set_field_order(root_set,MBTET,"FIELD2",order); CHKERRG(ierr);
97  ierr = m_field.set_field_order(root_set,MBTRI,"FIELD2",order); CHKERRG(ierr);
98  ierr = m_field.set_field_order(root_set,MBEDGE,"FIELD2",order); CHKERRG(ierr);
99  ierr = m_field.set_field_order(root_set,MBVERTEX,"FIELD2",1); CHKERRG(ierr);
100 
101  /****/
102  //build database
103  //build field
104  ierr = m_field.build_fields(); CHKERRG(ierr);
105  //build finite elemnts
106  ierr = m_field.build_finite_elements(); CHKERRG(ierr);
107  //build adjacencies
108  ierr = m_field.build_adjacencies(bit_level0); CHKERRG(ierr);
109 
110 
111  ProblemsManager *prb_mng_ptr;
112  ierr = m_field.getInterface(prb_mng_ptr); CHKERRG(ierr);
113  //build problem
114  ierr = prb_mng_ptr->buildProblem("TEST_PROBLEM",true); CHKERRG(ierr);
115  ierr = prb_mng_ptr->partitionSimpleProblem("TEST_PROBLEM"); CHKERRG(ierr);
116  ierr = prb_mng_ptr->partitionFiniteElements("TEST_PROBLEM"); CHKERRG(ierr);
117  ierr = prb_mng_ptr->partitionGhostDofs("TEST_PROBLEM"); CHKERRG(ierr);
118 
119  typedef tee_device<std::ostream, std::ofstream> TeeDevice;
120  typedef stream<TeeDevice> TeeStream;
121 
122  std::ofstream ofs("forces_and_sources_testing_ftensor.txt");
123  TeeDevice my_tee(std::cout, ofs);
124  TeeStream my_split(my_tee);
125 
127 
128  boost::shared_ptr<MatrixDouble> field1ValuesDataPtr;
129  boost::shared_ptr<VectorDouble> field2ValuesDataPtr;
130  boost::shared_ptr<MatrixDouble> grad1ValuesDataPtr;
131  boost::shared_ptr<MatrixDouble> grad2ValuesDataPtr;
132  TeeStream &my_split;
133 
134  MyOp1(
135  boost::shared_ptr<MatrixDouble> field1_values_data_ptr,
136  boost::shared_ptr<VectorDouble> field2_values_data_ptr,
137  boost::shared_ptr<MatrixDouble> grad1_values_data_ptr,
138  boost::shared_ptr<MatrixDouble> grad2_values_data_ptr,
139  TeeStream &_my_split,
140  char type
141  ):
143  field1ValuesDataPtr(field1_values_data_ptr),
144  field2ValuesDataPtr(field2_values_data_ptr),
145  grad1ValuesDataPtr(grad1_values_data_ptr),
146  grad2ValuesDataPtr(grad2_values_data_ptr),
147  my_split(_my_split) {}
148 
149 
150  MoFEMErrorCode doWork(
151  int side,
152  EntityType type,
155 
157  const int nb_gauss_pts = data.getN().size1();
158  const int nb_base_functions = data.getN().size2();
159 
162 
163  auto base_function = data.getFTensor0N();
164  auto diff_base = data.getFTensor1DiffN<3>();
165  // FTensor::Tensor1<double*,3> field_values = getFTensor1FromMat<3>(field1ValuesDataPtr);
166  // #ifdef WITH_ADOL_C
167  // adouble val;
168  // FTensor::Tensor0<adouble*> adouble_t0(&val);
169  // // cerr << endl;
170  // // cerr << base_function << endl;
171  // adouble_t0<<=base_function;
172  // adouble_t0>>=base_function;
173  // // cerr << adouble_t0 << " " << base_function << endl;
174  // if(adouble_t0 - base_function != 0) {
175  // SETERRQ(PETSC_COMM_SELF,MOFEM_DATA_INCONSISTENCY,"Data inconsistency");
176  // }
177  // FTensor::Tensor1<adouble,3> adouble_t1;
178  // adouble_t1(I)<<=diff_base(I);
179  // adouble_t1(I)>>=diff_base(I);
180  // for(int II = 0;II!=3;II++) {
181  // // cerr << adouble_t1(II) << " " << diff_base(II) << endl;
182  // if(adouble_t1(II)-diff_base(II)!=0) {
183  // SETERRQ(PETSC_COMM_SELF,MOFEM_DATA_INCONSISTENCY,"Data inconsistency");
184  // }
185  // }
186  // FTensor::Tensor2<adouble,3,3> adouble_t2;
187  // FTensor::Tensor2<double,3,3> double_t2;
188  // double_t2(I,J) = diff_base(I)*diff_base(J);
189  // adouble_t2(I,J)<<=double_t2(I,J);
190  // adouble_t2(I,J)>>=double_t2(I,J);
191  // cerr << endl;
192  // for(int II = 0;II!=3;II++) {
193  // for(int JJ = 0;JJ!=3;JJ++) {
194  // cerr << adouble_t2(II,JJ) << " " << double_t2(II,JJ) << endl;
195  // if(adouble_t2(II,JJ)-double_t2(II,JJ)!=0) {
196  // SETERRQ(PETSC_COMM_SELF,MOFEM_DATA_INCONSISTENCY,"Data inconsistency");
197  // }
198  // }
199  // }
200  // adouble_t2(0,0)<<=1;
201  // adouble_t1(0)<<=1;
202  // double a;
203  // adouble_t2(0,0)>>=a;
204  // adouble_t1(0)>>=a;
205  // #endif
206 
207  for(int gg = 0;gg!=nb_gauss_pts;gg++) {
208 
209  for(int bb = 0;bb!=nb_base_functions;bb++) {
210 
211  double base_val = base_function;
212  // my_split << base_val << endl;
213  if(base_val!=data.getN(gg)[bb]) {
214  SETERRQ(PETSC_COMM_SELF,MOFEM_DATA_INCONSISTENCY,"Data inconsistency");
215  }
216  ++base_function;
217 
218  for(int dd = 0;dd<3;dd++) {
219  if(diff_base(dd)!=data.getDiffN()(gg,3*bb+dd)) {
220  SETERRQ2(
221  PETSC_COMM_SELF,
223  "Data inconsistency gg = %d bb = %d",
224  gg,
225  bb
226  );
227  }
228  }
229 
230  t2(I,J) = diff_base(I)*diff_base(J);
231 
232  MatrixAdaptor mat = MatrixAdaptor(3,3,ublas::shallow_array_adaptor<double>(9,&t2(0,0)));
233 
234  for(int II = 0;II!=3;II++) {
235  for(int JJ = 0;JJ!=3;JJ++) {
236  if(t2(II,JJ)-mat(II,JJ)!=0) {
237  SETERRQ(PETSC_COMM_SELF,MOFEM_DATA_INCONSISTENCY,"Data inconsistency");
238  }
239  }
240  }
241 
242 
243  ++diff_base;
244 
245  }
246 
247  // VectorAdaptor vec = VectorAdaptor(3,ublas::shallow_array_adaptor<double>(3,&field_values(0)));
248  // my_split << vec << endl;
249 
250  // ++field_values;
251 
252  }
254  }
255 
256  MoFEMErrorCode doWork(
257  int row_side,int col_side,
258  EntityType row_type,EntityType col_type,
259  EntitiesFieldData::EntData &row_data,
261  ) {
263  //
264 
265  const int nb_gauss_pts = row_data.getN().size1();
266  const int nb_base_functions_row = row_data.getN().size2();
267  const int nb_base_functions_col = col_data.getN().size2();
268 
272 
276 
279 
280  for(int br = 0;br!=nb_base_functions_row;br++) {
281  auto base_row = row_data.getFTensor0N(0, br);
282  auto diff_base_row = row_data.getFTensor1DiffN<3>(0, br);
283 
284  for(int bc = 0;bc!=nb_base_functions_col;bc++) {
285  auto base_col = col_data.getFTensor0N(0, bc);
286  auto diff_base_col= row_data.getFTensor1DiffN<3>(0, bc);
287 
288  auto field1_values = getFTensor1FromMat<3>(*field1ValuesDataPtr);
289  auto field2_values = getFTensor0FromVec(*field2ValuesDataPtr);
290  auto grad1_values = getFTensor2FromMat<3,3>(*grad1ValuesDataPtr);
291  auto grad2_values = getFTensor1FromMat<3>(*grad2ValuesDataPtr);
292 
293  for(int gg = 0;gg!=nb_gauss_pts;gg++) {
294  // This make no sense (just do some calculations)
295  // FIXME: Some stuff can be calculated only in loop by integration pts for efficiency (I don't care for purpose of this test)
296  t2(I,J) = diff_base_row(I)*diff_base_col(J)*base_row*base_col;
297  t3(I,J,K) = grad1_values(I,J)*field1_values(K)*field2_values;
298  ++field1_values;
299  ++field2_values;
300  ++grad1_values;
301  ++grad2_values;
302  }
303  ++base_col;
304  ++diff_base_col;
305 
306  }
307  ++base_row;
308  ++diff_base_row;
309 
310  }
311 
313  }
314 
315  };
316 
317  boost::shared_ptr<MatrixDouble> values1_at_gauss_pts_ptr =
318  boost::shared_ptr<MatrixDouble>(new MatrixDouble());
319  boost::shared_ptr<VectorDouble> values2_at_gauss_pts_ptr =
320  boost::shared_ptr<VectorDouble>(new VectorDouble());
321  boost::shared_ptr<MatrixDouble> grad1_at_gauss_pts_ptr =
322  boost::shared_ptr<MatrixDouble>(new MatrixDouble());
323  boost::shared_ptr<MatrixDouble> grad2_at_gauss_pts_ptr =
324  boost::shared_ptr<MatrixDouble>(new MatrixDouble());
325 
327 
328  fe1.getOpPtrVector().push_back(
329  new OpCalculateVectorFieldValues<3>("FIELD1", values1_at_gauss_pts_ptr));
330  fe1.getOpPtrVector().push_back(
331  new OpCalculateScalarFieldValues("FIELD2", values2_at_gauss_pts_ptr));
333  "FIELD1", grad1_at_gauss_pts_ptr));
334  fe1.getOpPtrVector().push_back(
335  new OpCalculateScalarFieldGradient<3>("FIELD2", grad2_at_gauss_pts_ptr));
336 
337  fe1.getOpPtrVector().push_back(
338  new MyOp1(
339  values1_at_gauss_pts_ptr,
340  values2_at_gauss_pts_ptr,
341  grad1_at_gauss_pts_ptr,
342  grad2_at_gauss_pts_ptr,
343  my_split,
345  )
346  );
347  fe1.getOpPtrVector().push_back(
348  new MyOp1(
349  values1_at_gauss_pts_ptr,
350  values2_at_gauss_pts_ptr,
351  grad1_at_gauss_pts_ptr,
352  grad2_at_gauss_pts_ptr,
353  my_split,
355  )
356  );
357 
358  ierr = m_field.loop_finite_elements("TEST_PROBLEM","TEST_FE1",fe1); CHKERRG(ierr);
359 
360 
361  } catch (MoFEMException const &e) {
362  SETERRQ(PETSC_COMM_SELF,e.errorCode,e.errorMessage);
363  }
364 
366 
367  return 0;
368 
369 }
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.
Definition: UnknownInterface.hpp:93
MoFEM::K
VectorDouble K
Definition: Projection10NodeCoordsOnField.cpp:125
MoFEM::EntitiesFieldData::EntData
Data on single entity (This is passed as argument to DataOperator::doWork)
Definition: EntitiesFieldData.hpp:127
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
EntityHandle
MoFEM::ProblemsManager
Problem manager is used to build and partition problems.
Definition: ProblemsManager.hpp:21
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::Exceptions::MoFEMException
Exception to catch.
Definition: Exceptions.hpp:20
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
MoFEM::OpCalculateVectorFieldValues
Get values at integration pts for tensor filed rank 1, i.e. vector field.
Definition: UserDataOperators.hpp:466
MoFEM::Types::MatrixDouble
UBlasMatrix< double > MatrixDouble
Definition: Types.hpp:77
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:1316
J
FTensor::Index< 'J', DIM1 > J
Definition: level_set.cpp:30
MoFEM::ForcesAndSourcesCore::UserDataOperator::OPROWCOL
@ OPROWCOL
operator doWork is executed on FE rows &columns
Definition: ForcesAndSourcesCore.hpp:569
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.
FTensor::Tensor2< double, 3, 3 >
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::OpCalculateScalarFieldGradient
Get field gradients at integration pts for scalar filed rank 0, i.e. vector field.
Definition: UserDataOperators.hpp:1294
I
constexpr IntegrationType I
Definition: operators_tests.cpp:31
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:1975
MoFEM::EntitiesFieldData::EntData::getFTensor0N
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0N(const FieldApproximationBase base)
Get base function as Tensor0.
Definition: EntitiesFieldData.hpp:1489
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
MoFEM::Exceptions::rval
static MoFEMErrorCodeGeneric< moab::ErrorCode > rval
Definition: Exceptions.hpp:74
FTensor::Number< 0 >
FTensor::Tensor3
Definition: Tensor3_value.hpp:12
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.
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.
convert.type
type
Definition: convert.py:64
main
int main(int argc, char *argv[])
Definition: forces_and_sources_testing_ftensor.cpp:15
MoFEM::getFTensor0FromVec
static auto getFTensor0FromVec(ublas::vector< T, A > &data)
Get tensor rank 0 (scalar) form data vector.
Definition: Templates.hpp:135
MoFEM::OpCalculateScalarFieldValues
Get value at integration points for scalar field.
Definition: UserDataOperators.hpp:82
MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator
friend class UserDataOperator
Definition: VolumeElementForcesAndSourcesCore.hpp:105
MoFEM::Exceptions::MoFEMException::errorMessage
char errorMessage[1024]
Definition: Exceptions.hpp:22
MoFEM::VolumeElementForcesAndSourcesCore
Volume finite element base.
Definition: VolumeElementForcesAndSourcesCore.hpp:26
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
MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator
Definition: VolumeElementForcesAndSourcesCore.hpp:108
FTensor::Index
Definition: Index.hpp:23
mesh_file_name
char mesh_file_name[255]
Definition: mesh_smoothing.cpp:23
MoFEM::OpCalculateVectorFieldGradient
Get field gradients at integration pts for scalar filed rank 0, i.e. vector field.
Definition: UserDataOperators.hpp:1536
MoFEM::Exceptions::MoFEMException::errorCode
const int errorCode
Definition: Exceptions.hpp:21
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
FTensor::dd
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
MoFEM::EntitiesFieldData::EntData::getFTensor1DiffN
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1DiffN(const FieldApproximationBase base)
Get derivatives of base functions.
Definition: EntitiesFieldData.cpp:526
help
static char help[]
Definition: forces_and_sources_testing_ftensor.cpp:13
MoFEM::EntitiesFieldData::EntData::getN
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
Definition: EntitiesFieldData.hpp:1305
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1094
MoFEM::ProblemsManager::partitionSimpleProblem
MoFEMErrorCode partitionSimpleProblem(const std::string name, int verb=VERBOSE)
partition problem dofs
Definition: ProblemsManager.cpp:1537
MoFEM::Exceptions::ierr
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Exceptions.hpp:76
AINSWORTH_LEGENDRE_BASE
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
MoFEM::ForcesAndSourcesCore::getOpPtrVector
boost::ptr_deque< UserDataOperator > & getOpPtrVector()
Use to push back operator for row operator.
Definition: ForcesAndSourcesCore.hpp:83
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
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::CoreInterface::modify_finite_element_add_field_data
virtual MoFEMErrorCode modify_finite_element_add_field_data(const std::string &fe_name, const std::string name_filed)=0
set finite element field data
MoFEM::Types::VectorDouble
UBlasVector< double > VectorDouble
Definition: Types.hpp:68
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::MatrixAdaptor
MatrixShallowArrayAdaptor< double > MatrixAdaptor
Matrix adaptor.
Definition: Types.hpp:132
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:440
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.
MoFEM::ProblemsManager::partitionGhostDofs
MoFEMErrorCode partitionGhostDofs(const std::string name, int verb=VERBOSE)
determine ghost nodes
Definition: ProblemsManager.cpp:2339
MoFEM::CoreInterface::add_problem
virtual MoFEMErrorCode add_problem(const std::string &name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add problem.
CHKERRG
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:483
MoFEM::ForcesAndSourcesCore::UserDataOperator::OPROW
@ OPROW
operator doWork function is executed on FE rows
Definition: ForcesAndSourcesCore.hpp:567