v0.14.0
Functions | Variables
forces_and_sources_testing_triangle_element.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"
 

Function Documentation

◆ main()

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

Definition at line 13 of file forces_and_sources_testing_triangle_element.cpp.

13  {
14 
15  MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
16 
17  try {
18 
19  moab::Core mb_instance;
20  moab::Interface &moab = mb_instance;
21  int rank;
22  MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
23 
24  PetscBool flg = PETSC_TRUE;
25  char mesh_file_name[255];
26 #if PETSC_VERSION_GE(3, 6, 4)
27  CHKERR PetscOptionsGetString(PETSC_NULL, "", "-my_file", mesh_file_name,
28  255, &flg);
29 #else
30  CHKERR PetscOptionsGetString(PETSC_NULL, PETSC_NULL, "-my_file",
31  mesh_file_name, 255, &flg);
32 #endif
33  if (flg != PETSC_TRUE) {
34  SETERRQ(PETSC_COMM_SELF, 1, "*** ERROR -my_file (MESH FILE NEEDED)");
35  }
36 
37  // Create MoFEM (Joseph) database
38  MoFEM::Core core(moab);
39  MoFEM::Interface &m_field = core;
40 
41  const char *option;
42  option = "";
43  CHKERR moab.load_file(mesh_file_name, 0, option);
44 
45  // set entitities bit level
46  BitRefLevel bit_level0;
47  bit_level0.set(0);
48  EntityHandle meshset_level0;
49  CHKERR moab.create_meshset(MESHSET_SET, meshset_level0);
50  CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(
51  0, 3, bit_level0);
52 
53  // Fields
54  CHKERR m_field.add_field("FIELD1", H1, AINSWORTH_LEGENDRE_BASE, 3);
55  CHKERR m_field.add_field("FIELD2", NOFIELD, NOBASE, 3);
56 
57  {
58  // Creating and adding no field entities.
59  const double coords[] = {0, 0, 0};
60  EntityHandle no_field_vertex;
61  CHKERR m_field.get_moab().create_vertex(coords, no_field_vertex);
62  Range range_no_field_vertex;
63  range_no_field_vertex.insert(no_field_vertex);
64  CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevel(
65  range_no_field_vertex, BitRefLevel().set());
66  EntityHandle meshset = m_field.get_field_meshset("FIELD2");
67  CHKERR m_field.get_moab().add_entities(meshset, range_no_field_vertex);
68  }
69 
70  // FE
71  CHKERR m_field.add_finite_element("TEST_FE1");
72  CHKERR m_field.add_finite_element("TEST_FE2");
73 
74  // Define rows/cols and element data
75  CHKERR m_field.modify_finite_element_add_field_row("TEST_FE1", "FIELD1");
76  CHKERR m_field.modify_finite_element_add_field_col("TEST_FE1", "FIELD1");
77  CHKERR m_field.modify_finite_element_add_field_data("TEST_FE1", "FIELD1");
78 
79  CHKERR m_field.modify_finite_element_add_field_row("TEST_FE2", "FIELD2");
80  CHKERR m_field.modify_finite_element_add_field_col("TEST_FE2", "FIELD1");
81  CHKERR m_field.modify_finite_element_add_field_data("TEST_FE2", "FIELD1");
82  CHKERR m_field.modify_finite_element_add_field_data("TEST_FE2", "FIELD2");
83 
84  // Problem
85  CHKERR m_field.add_problem("TEST_PROBLEM");
86 
87  // set finite elements for problem
88  CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM",
89  "TEST_FE1");
90  CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM",
91  "TEST_FE2");
92  // set refinement level for problem
93  CHKERR m_field.modify_problem_ref_level_add_bit("TEST_PROBLEM", bit_level0);
94 
95  // meshset consisting all entities in mesh
96  EntityHandle root_set = moab.get_root_set();
97  // add entities to field
98  CHKERR m_field.add_ents_to_field_by_type(root_set, MBTET, "FIELD1");
99  // add entities to finite element
100  Range tets;
101  CHKERR moab.get_entities_by_type(0, MBTET, tets, false);
102  Skinner skin(&m_field.get_moab());
103  Range tets_skin;
104  CHKERR skin.find_skin(0, tets, false, tets_skin);
105  CHKERR m_field.add_ents_to_finite_element_by_type(tets_skin, MBTRI,
106  "TEST_FE1");
107  CHKERR m_field.add_ents_to_finite_element_by_type(tets_skin, MBTRI,
108  "TEST_FE2");
109 
110  // set app. order
111  // see Hierarchic Finite Element Bases on Unstructured Tetrahedral Meshes
112  // (Mark Ainsworth & Joe Coyle)
113  int order = 3;
114  CHKERR m_field.set_field_order(root_set, MBTET, "FIELD1", order);
115  CHKERR m_field.set_field_order(root_set, MBTRI, "FIELD1", order);
116  CHKERR m_field.set_field_order(root_set, MBEDGE, "FIELD1", order);
117  CHKERR m_field.set_field_order(root_set, MBVERTEX, "FIELD1", 1);
118 
119  /****/
120  // build database
121  // build field
122  CHKERR m_field.build_fields();
123  // set FIELD1 from positions of 10 node tets
124  Projection10NodeCoordsOnField ent_method(m_field, "FIELD1");
125  CHKERR m_field.loop_dofs("FIELD1", ent_method);
126  // build finite elemnts
127  CHKERR m_field.build_finite_elements();
128  // build adjacencies
129  CHKERR m_field.build_adjacencies(bit_level0);
130  // build problem
131  ProblemsManager *prb_mng_ptr;
132  CHKERR m_field.getInterface(prb_mng_ptr);
133  CHKERR prb_mng_ptr->buildProblem("TEST_PROBLEM", true);
134 
135  /****/
136  // mesh partitioning
137  // partition
138  CHKERR prb_mng_ptr->partitionSimpleProblem("TEST_PROBLEM");
139  CHKERR prb_mng_ptr->partitionFiniteElements("TEST_PROBLEM");
140  // what are ghost nodes, see Petsc Manual
141  CHKERR prb_mng_ptr->partitionGhostDofs("TEST_PROBLEM");
142 
143  typedef tee_device<std::ostream, std::ofstream> TeeDevice;
144  typedef stream<TeeDevice> TeeStream;
145 
146  std::ofstream ofs("forces_and_sources_testing_triangle_element.txt");
147  TeeDevice my_tee(std::cout, ofs);
148  TeeStream my_split(my_tee);
149 
151 
152  TeeStream &my_split;
153  MyOp1(TeeStream &_my_split, const char type)
155  "FIELD1", type),
156  my_split(_my_split) {}
157 
158  MoFEMErrorCode doWork(int side, EntityType type,
161 
162  const double eps = 1e-4;
163  for (DoubleAllocator::iterator it = getNormal().data().begin();
164  it != getNormal().data().end(); it++) {
165  *it = fabs(*it) < eps ? 0.0 : *it;
166  }
167 
168  my_split << "NH1" << std::endl;
169  my_split << "side: " << side << " type: " << type << std::endl;
170  my_split << "data: " << data << std::endl;
171  my_split << std::setprecision(3) << getCoords() << std::endl;
172  my_split << std::setprecision(3) << getCoordsAtGaussPts() << std::endl;
173  my_split << std::setprecision(3) << getArea() << std::endl;
174  my_split << std::setprecision(3) << getNormal() << std::endl;
175  my_split << std::setprecision(3) << getNormalsAtGaussPts() << std::endl;
176  my_split << std::setprecision(3) << getTangent1AtGaussPts()
177  << std::endl;
178  my_split << std::setprecision(3) << getTangent2AtGaussPts()
179  << std::endl;
180  my_split << std::endl;
182  }
183 
184  MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type,
185  EntityType col_type,
186  EntitiesFieldData::EntData &row_data,
187  EntitiesFieldData::EntData &col_data) {
188 
190  my_split << "NH1NH1" << std::endl;
191  my_split << "row side: " << row_side << " row_type: " << row_type
192  << std::endl;
193  my_split << row_data << std::endl;
194  my_split << "NH1NH1" << std::endl;
195  my_split << "col side: " << col_side << " col_type: " << col_type
196  << std::endl;
197  my_split << row_data << std::endl;
198 
199  VectorInt row_indices, col_indices;
200  CHKERR getProblemRowIndices("FIELD1", row_type, row_side, row_indices);
201  CHKERR getProblemColIndices("FIELD1", col_type, col_side, col_indices);
202 
203  if (row_indices.size() != row_data.getIndices().size()) {
204  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
205  "row inconsistency");
206  }
207 
208  if (col_indices.size() != col_data.getIndices().size()) {
209  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
210  "col inconsistency");
211  }
212 
213  for (unsigned int rr = 0; rr < row_indices.size(); rr++) {
214  if (row_indices[rr] != row_data.getIndices()[rr]) {
215  std::cerr << row_indices << std::endl;
216  std::cerr << row_data.getIndices() << std::endl;
217  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
218  "row inconsistency");
219  }
220  }
221 
222  for (unsigned int cc = 0; cc < col_indices.size(); cc++) {
223  if (col_indices[cc] != col_data.getIndices()[cc]) {
224  std::cerr << col_indices << std::endl;
225  std::cerr << col_data.getIndices() << std::endl;
226  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
227  "row inconsistency");
228  }
229  }
230 
231  my_split << row_data << std::endl;
232 
234  }
235  };
236 
238 
239  TeeStream &my_split;
240  MyOp2(TeeStream &_my_split, const char type)
242  "FIELD1", type),
243  my_split(_my_split) {
244  sYmm = false;
245  }
246 
247  MoFEMErrorCode doWork(int side, EntityType type,
250 
251  if (type != MBENTITYSET)
253 
254  my_split << "NOFIELD" << std::endl;
255  my_split << "side: " << side << " type: " << type << std::endl;
256  my_split << data << std::endl;
258  }
259 
260  MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type,
261  EntityType col_type,
262  EntitiesFieldData::EntData &row_data,
263  EntitiesFieldData::EntData &col_data) {
265 
266  unSetSymm();
267 
268  if (row_type != MBENTITYSET)
270 
271  my_split << "NOFILEDH1" << std::endl;
272  my_split << "row side: " << row_side << " row_type: " << row_type
273  << std::endl;
274  my_split << row_data << std::endl;
275  my_split << "col side: " << col_side << " col_type: " << col_type
276  << std::endl;
277  my_split << col_data << std::endl;
278 
280  }
281  };
282 
283  FaceElementForcesAndSourcesCore fe1(m_field);
284  fe1.getOpPtrVector().push_back(
285  new MyOp1(my_split, ForcesAndSourcesCore::UserDataOperator::OPROW));
286  fe1.getOpPtrVector().push_back(
287  new MyOp1(my_split, ForcesAndSourcesCore::UserDataOperator::OPROWCOL));
288 
289  FaceElementForcesAndSourcesCore fe2(m_field);
290  fe2.getOpPtrVector().push_back(
291  new MyOp2(my_split, ForcesAndSourcesCore::UserDataOperator::OPROW));
292  fe2.getOpPtrVector().push_back(
293  new MyOp2(my_split, ForcesAndSourcesCore::UserDataOperator::OPROWCOL));
294 
295  CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "TEST_FE1", fe1);
296  CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "TEST_FE2", fe2);
297  }
298  CATCH_ERRORS;
299 
301 
302  return 0;
303 }

Variable Documentation

◆ help

char help[] = "...\n\n"
static
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::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::CoreInterface::loop_dofs
virtual MoFEMErrorCode loop_dofs(const Problem *problem_ptr, const std::string &field_name, RowColData rc, DofMethod &method, int lower_rank, int upper_rank, int verb=DEFAULT_VERBOSITY)=0
Make a loop over dofs.
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::ForcesAndSourcesCore::UserDataOperator::getProblemColIndices
MoFEMErrorCode getProblemColIndices(const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
Get col indices.
Definition: ForcesAndSourcesCore.cpp:1637
EntityHandle
MyOp2::MyOp2
MyOp2(const char type, const char face_type)
Definition: forces_and_sources_testing_contact_prism_element.cpp:133
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
NOBASE
@ NOBASE
Definition: definitions.h:59
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
MoFEM::Projection10NodeCoordsOnField
Projection of edge entities with one mid-node on hierarchical basis.
Definition: Projection10NodeCoordsOnField.hpp:24
MoFEM::CoreTmp< 0 >::Finalize
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:112
MoFEM::CoreInterface::get_field_meshset
virtual EntityHandle get_field_meshset(const std::string name) const =0
get field meshset
help
static char help[]
Definition: forces_and_sources_testing_triangle_element.cpp:11
MoFEM::DataOperator::doWork
virtual MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
Operator for bi-linear form, usually to calculate values on left hand side.
Definition: DataOperators.hpp:40
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:1975
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:535
MoFEM::FaceElementForcesAndSourcesCore::UserDataOperator::getNormalsAtGaussPts
MatrixDouble & getNormalsAtGaussPts()
if higher order geometry return normals at Gauss pts.
Definition: FaceElementForcesAndSourcesCore.hpp:290
MoFEM::FaceElementForcesAndSourcesCore::UserDataOperator::getTangent2AtGaussPts
MatrixDouble & getTangent2AtGaussPts()
if higher order geometry return tangent vector to triangle at Gauss pts.
Definition: FaceElementForcesAndSourcesCore.hpp:310
MoFEM::FaceElementForcesAndSourcesCore::UserDataOperator::getTangent1AtGaussPts
MatrixDouble & getTangent1AtGaussPts()
if higher order geometry return tangent vector to triangle at Gauss pts.
Definition: FaceElementForcesAndSourcesCore.hpp:304
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::CoreInterface::get_moab
virtual moab::Interface & get_moab()=0
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::ForcesAndSourcesCore::UserDataOperator::getProblemRowIndices
MoFEMErrorCode getProblemRowIndices(const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
Get row indices.
Definition: ForcesAndSourcesCore.cpp:1620
MoFEM::CoreInterface::build_finite_elements
virtual MoFEMErrorCode build_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elements.
MoFEM::ForcesAndSourcesCore::UserDataOperator
Definition: ForcesAndSourcesCore.hpp:549
MoFEM::FaceElementForcesAndSourcesCore::UserDataOperator
default operator for TRI element
Definition: FaceElementForcesAndSourcesCore.hpp:94
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
MoFEM::FaceElementForcesAndSourcesCore::UserDataOperator::getNormal
VectorDouble & getNormal()
get triangle normal
Definition: FaceElementForcesAndSourcesCore.hpp:243
MoFEM::EntitiesFieldData::EntData::getIndices
const VectorInt & getIndices() const
Get global indices of dofs on entity.
Definition: EntitiesFieldData.hpp:1201
MoFEM::FaceElementForcesAndSourcesCore
Face finite element.
Definition: FaceElementForcesAndSourcesCore.hpp:23
MoFEM::ForcesAndSourcesCore::UserDataOperator::getCoordsAtGaussPts
MatrixDouble & getCoordsAtGaussPts()
Gauss points and weight, matrix (nb. of points x 3)
Definition: ForcesAndSourcesCore.hpp:1264
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:23
Range
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
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:372
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::Types::VectorInt
UBlasVector< int > VectorInt
Definition: Types.hpp:67
eps
static const double eps
Definition: check_base_functions_derivatives_on_tet.cpp:11
MyOp2
Definition: forces_and_sources_testing_contact_prism_element.cpp:130
AINSWORTH_LEGENDRE_BASE
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
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
MyOp2::doWork
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Definition: forces_and_sources_testing_contact_prism_element.cpp:137
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:440
MoFEM::FaceElementForcesAndSourcesCore::UserDataOperator::getArea
double getArea()
get area of face
Definition: FaceElementForcesAndSourcesCore.hpp:239
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.
MoFEM::DataOperator::sYmm
bool sYmm
If true assume that matrix is symmetric structure.
Definition: DataOperators.hpp:82
MoFEM::DataOperator::unSetSymm
void unSetSymm()
unset if operator is executed for non symmetric problem
Definition: DataOperators.hpp:113
NOFIELD
@ NOFIELD
scalar or vector of scalars describe (no true field)
Definition: definitions.h:84
MoFEM::FaceElementForcesAndSourcesCore::UserDataOperator::getCoords
VectorDouble & getCoords()
get triangle coordinates
Definition: FaceElementForcesAndSourcesCore.hpp:279