v0.14.0
thermal_with_analytical_bc.cpp
Go to the documentation of this file.
1 
2 
3 #include <boost/iostreams/tee.hpp>
4 #include <boost/iostreams/stream.hpp>
5 #include <fstream>
6 #include <iostream>
7 
8 namespace bio = boost::iostreams;
9 using bio::stream;
10 using bio::tee_device;
11 
12 #include <MoFEM.hpp>
13 using namespace MoFEM;
14 
15 #include <MethodForForceScaling.hpp>
16 #include <DirichletBC.hpp>
17 #include <PostProcOnRefMesh.hpp>
18 #include <ThermalElement.hpp>
19 
21 
22 #include <boost/shared_ptr.hpp>
23 #include <boost/numeric/ublas/vector_proxy.hpp>
24 #include <AnalyticalDirichlet.hpp>
25 
26 #include <boost/iostreams/tee.hpp>
27 #include <boost/iostreams/stream.hpp>
28 #include <fstream>
29 #include <iostream>
30 
31 static int debug = 1;
32 static char help[] = "...\n\n";
33 
35 
36  std::vector<VectorDouble> val;
37 
38  std::vector<VectorDouble> &operator()(double x, double y, double z) {
39  val.resize(1);
40  val[0].resize(1);
41  (val[0])[0] = pow(x, 1);
42  return val;
43  }
44 };
45 
46 int main(int argc, char *argv[]) {
47 
48  MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
49 
50  try {
51 
52  moab::Core mb_instance;
53  moab::Interface &moab = mb_instance;
54  int rank;
55  MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
56 
57  PetscBool flg = PETSC_TRUE;
58  char mesh_file_name[255];
59  CHKERR PetscOptionsGetString(PETSC_NULL, PETSC_NULL, "-my_file",
60  mesh_file_name, 255, &flg);
61  if (flg != PETSC_TRUE) {
62  SETERRQ(PETSC_COMM_SELF, 1, "*** ERROR -my_file (MESH FILE NEEDED)");
63  }
64 
65  const char *option;
66  option = "";
67  CHKERR moab.load_file(mesh_file_name, 0, option);
68 
69  // Create MoFEM (Joseph) database
70  MoFEM::Core core(moab);
71  MoFEM::Interface &m_field = core;
72 
73  // set entitities bit level
74  BitRefLevel bit_level0;
75  bit_level0.set(0);
76  EntityHandle meshset_level0;
77  CHKERR moab.create_meshset(MESHSET_SET, meshset_level0);
78  CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(
79  0, 3, bit_level0);
80 
81  // Fields
82  CHKERR m_field.add_field("TEMP", H1, AINSWORTH_LEGENDRE_BASE, 1);
83 
84  // Problem
85  CHKERR m_field.add_problem("TEST_PROBLEM");
86  CHKERR m_field.add_problem("BC_PROBLEM");
87 
88  // set refinement level for problem
89  CHKERR m_field.modify_problem_ref_level_add_bit("TEST_PROBLEM", bit_level0);
90  CHKERR m_field.modify_problem_ref_level_add_bit("BC_PROBLEM", bit_level0);
91 
92  // meshset consisting all entities in mesh
93  EntityHandle root_set = moab.get_root_set();
94  CHKERR m_field.add_field("MESH_NODE_POSITIONS", H1, AINSWORTH_LEGENDRE_BASE,
95  3);
96  CHKERR m_field.add_ents_to_field_by_type(root_set, MBTET,
97  "MESH_NODE_POSITIONS");
98  CHKERR m_field.set_field_order(0, MBTET, "MESH_NODE_POSITIONS", 2);
99  CHKERR m_field.set_field_order(0, MBTRI, "MESH_NODE_POSITIONS", 2);
100  CHKERR m_field.set_field_order(0, MBEDGE, "MESH_NODE_POSITIONS", 2);
101  CHKERR m_field.set_field_order(0, MBVERTEX, "MESH_NODE_POSITIONS", 1);
102 
103  // add entities to field
104  CHKERR m_field.add_ents_to_field_by_type(root_set, MBTET, "TEMP");
105 
106  // set app. order
107  // see Hierarchic Finite Element Bases on Unstructured Tetrahedral Meshes
108  // (Mark Ainsworth & Joe Coyle)
109  PetscInt order;
110  CHKERR PetscOptionsGetInt(PETSC_NULL, PETSC_NULL, "-my_order", &order,
111  &flg);
112  if (flg != PETSC_TRUE) {
113  order = 2;
114  }
115  CHKERR m_field.set_field_order(root_set, MBTET, "TEMP", order);
116  CHKERR m_field.set_field_order(root_set, MBTRI, "TEMP", order);
117  CHKERR m_field.set_field_order(root_set, MBEDGE, "TEMP", order);
118  CHKERR m_field.set_field_order(root_set, MBVERTEX, "TEMP", 1);
119 
120  ThermalElement thermal_elements(m_field);
121  CHKERR thermal_elements.addThermalElements("TEMP");
122  CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM",
123  "THERMAL_FE");
124 
125  Range bc_tris;
126  for (_IT_CUBITMESHSETS_BY_NAME_FOR_LOOP_(m_field, "ANALYTICAL_BC", it)) {
127  CHKERR moab.get_entities_by_type(it->getMeshset(), MBTRI, bc_tris, true);
128  }
129 
130  AnalyticalDirichletBC analytical_bc(m_field);
131  CHKERR analytical_bc.setFiniteElement(m_field, "BC_FE", "TEMP", bc_tris);
132  CHKERR m_field.modify_problem_add_finite_element("BC_PROBLEM", "BC_FE");
133 
134  /****/
135  // build database
136  // build field
137  CHKERR m_field.build_fields();
138  // build finite elemnts
139  CHKERR m_field.build_finite_elements();
140  // build adjacencies
141  CHKERR m_field.build_adjacencies(bit_level0);
142  // build problem
143  ProblemsManager *prb_mng_ptr;
144  CHKERR m_field.getInterface(prb_mng_ptr);
145  CHKERR prb_mng_ptr->buildProblem("TEST_PROBLEM", true);
146  CHKERR prb_mng_ptr->buildProblem("BC_PROBLEM", true);
147 
148  Projection10NodeCoordsOnField ent_method_material(m_field,
149  "MESH_NODE_POSITIONS");
150  CHKERR m_field.loop_dofs("MESH_NODE_POSITIONS", ent_method_material);
151 
152  /****/
153  // mesh partitioning
154  // partition
155  CHKERR prb_mng_ptr->partitionSimpleProblem("TEST_PROBLEM");
156  CHKERR prb_mng_ptr->partitionFiniteElements("TEST_PROBLEM");
157  // what are ghost nodes, see Petsc Manual
158  CHKERR prb_mng_ptr->partitionGhostDofs("TEST_PROBLEM");
159 
160  CHKERR prb_mng_ptr->partitionSimpleProblem("BC_PROBLEM");
161  CHKERR prb_mng_ptr->partitionFiniteElements("BC_PROBLEM");
162  // what are ghost nodes, see Petsc Manual
163  CHKERR prb_mng_ptr->partitionGhostDofs("BC_PROBLEM");
164 
165  Vec F;
166  CHKERR m_field.getInterface<VecManager>()->vecCreateGhost("TEST_PROBLEM",
167  ROW, &F);
168  Vec T;
169  CHKERR VecDuplicate(F, &T);
170  Mat A;
172  ->createMPIAIJWithArrays<PetscGlobalIdx_mi_tag>("TEST_PROBLEM", &A);
173 
174  CHKERR addHOOpsVol("MESH_NODE_POSITIONS", thermal_elements.getLoopFeRhs(),
175  true, false, false, false);
176  CHKERR addHOOpsVol("MESH_NODE_POSITIONS", thermal_elements.getLoopFeLhs(),
177  true, false, false, false);
178 
179  CHKERR thermal_elements.setThermalFiniteElementRhsOperators("TEMP", F);
180  CHKERR thermal_elements.setThermalFiniteElementLhsOperators("TEMP", A);
181  CHKERR thermal_elements.setThermalFluxFiniteElementRhsOperators("TEMP", F);
182 
183  CHKERR VecZeroEntries(T);
184  CHKERR VecZeroEntries(F);
185  CHKERR MatZeroEntries(A);
186 
187  // analytical Dirichlet bc
188  AnalyticalDirichletBC::DirichletBC analytical_dirichlet_bc(m_field, "TEMP",
189  A, T, F);
190 
191  // solve for dirichlet bc dofs
192  CHKERR analytical_bc.setUpProblem(m_field, "BC_PROBLEM");
193 
194  boost::shared_ptr<AnalyticalFunction> testing_function =
195  boost::shared_ptr<AnalyticalFunction>(new AnalyticalFunction);
196 
197  CHKERR analytical_bc.setApproxOps(m_field, "TEMP", testing_function, 0);
198  CHKERR analytical_bc.solveProblem(m_field, "BC_PROBLEM", "BC_FE",
199  analytical_dirichlet_bc);
200 
201  CHKERR analytical_bc.destroyProblem();
202 
203  // preproc
204  CHKERR m_field.problem_basic_method_preProcess("TEST_PROBLEM",
205  analytical_dirichlet_bc);
206 
207  CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "THERMAL_FE",
208  thermal_elements.getLoopFeRhs());
209  CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "THERMAL_FE",
210  thermal_elements.getLoopFeLhs());
211  if (m_field.check_finite_element("THERMAL_FLUX_FE"))
212  CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "THERMAL_FLUX_FE",
213  thermal_elements.getLoopFeFlux());
214 
215  // postproc
216  CHKERR m_field.problem_basic_method_postProcess("TEST_PROBLEM",
217  analytical_dirichlet_bc);
218 
219  CHKERR VecGhostUpdateBegin(F, ADD_VALUES, SCATTER_REVERSE);
220  CHKERR VecGhostUpdateEnd(F, ADD_VALUES, SCATTER_REVERSE);
221  CHKERR VecAssemblyBegin(F);
222  CHKERR VecAssemblyEnd(F);
223  CHKERR MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);
224  CHKERR MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);
225 
226  CHKERR VecScale(F, -1);
227  // std::string wait;
228  // std::cout << "\n matrix is coming = \n" << std::endl;
229  // CHKERR MatView(A,PETSC_VIEWER_DRAW_WORLD);
230  // std::cin >> wait;
231 
232  // Solver
233  KSP solver;
234  CHKERR KSPCreate(PETSC_COMM_WORLD, &solver);
235  CHKERR KSPSetOperators(solver, A, A);
236  CHKERR KSPSetFromOptions(solver);
237  CHKERR KSPSetUp(solver);
238 
239  CHKERR KSPSolve(solver, F, T);
240  CHKERR VecGhostUpdateBegin(T, INSERT_VALUES, SCATTER_FORWARD);
241  CHKERR VecGhostUpdateEnd(T, INSERT_VALUES, SCATTER_FORWARD);
242 
243  // Save data on mesh
244  // CHKERR
245  // m_field.problem_basic_method_preProcess("TEST_PROBLEM",analytical_dirichlet_bc);
246  CHKERR m_field.getInterface<VecManager>()->setGlobalGhostVector(
247  "TEST_PROBLEM", ROW, T, ADD_VALUES, SCATTER_REVERSE);
248  CHKERR m_field.getInterface<VecManager>()->setLocalGhostVector(
249  "TEST_PROBLEM", ROW, T, INSERT_VALUES, SCATTER_FORWARD);
250 
251  // CHKERR VecView(F,PETSC_VIEWER_STDOUT_WORLD);
252 
253  // CHKERR VecView(T,PETSC_VIEWER_STDOUT_WORLD);
254 
255  PetscReal pointwisenorm;
256  CHKERR VecMax(T, NULL, &pointwisenorm);
257  std::cout << "\n The Global Pointwise Norm of error for this problem is : "
258  << pointwisenorm << std::endl;
259 
260  // PetscViewer viewer;
261  // PetscViewerASCIIOpen(PETSC_COMM_WORLD,"thermal_with_analytical_bc.txt",&viewer);
262  // CHKERR VecChop(T,1e-4);
263  // CHKERR VecView(T,viewer);
264  // CHKERR PetscViewerDestroy(&viewer);
265 
266  double sum = 0;
267  CHKERR VecSum(T, &sum);
268  CHKERR PetscPrintf(PETSC_COMM_WORLD, "sum = %9.8e\n", sum);
269  double fnorm;
270  CHKERR VecNorm(T, NORM_2, &fnorm);
271  CHKERR PetscPrintf(PETSC_COMM_WORLD, "fnorm = %9.8e\n", fnorm);
272  if (fabs(sum + 6.46079983e-01) > 1e-7) {
273  SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID, "Failed to pass test");
274  }
275  if (fabs(fnorm - 4.26080052e+00) > 1e-6) {
276  SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID, "Failed to pass test");
277  }
278 
279  if (debug) {
280 
281  PostProcVolumeOnRefinedMesh post_proc(m_field);
283  CHKERR addHOOpsVol("MESH_NODE_POSITIONS", post_proc, true, false, false,
284  false);
285  CHKERR post_proc.addFieldValuesPostProc("TEMP");
286  CHKERR post_proc.addFieldValuesPostProc("MESH_NODE_POSITIONS");
287  CHKERR post_proc.addFieldValuesGradientPostProc("TEMP");
288  CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "THERMAL_FE",
289  post_proc);
290  CHKERR post_proc.writeFile("out.h5m");
291  }
292 
293  CHKERR MatDestroy(&A);
294  CHKERR VecDestroy(&F);
295  CHKERR VecDestroy(&T);
296  CHKERR KSPDestroy(&solver);
297  }
298  CATCH_ERRORS;
299 
301 
302  return 0;
303 }
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Definition: UnknownInterface.hpp:93
AnalyticalDirichletBC::setFiniteElement
MoFEMErrorCode setFiniteElement(MoFEM::Interface &m_field, string fe, string field, Range &tris, string nodals_positions="MESH_NODE_POSITIONS")
set finite element
Definition: AnalyticalDirichlet.cpp:192
MoFEM::CoreInterface::problem_basic_method_postProcess
virtual MoFEMErrorCode problem_basic_method_postProcess(const Problem *problem_ptr, BasicMethod &method, int verb=DEFAULT_VERBOSITY)=0
Set data for BasicMethod.
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
EntityHandle
MoFEM::ProblemsManager
Problem manager is used to build and partition problems.
Definition: ProblemsManager.hpp:21
MoFEM::addHOOpsVol
MoFEMErrorCode addHOOpsVol(const std::string field, E &e, bool h1, bool hcurl, bool hdiv, bool l2)
Definition: HODataOperators.hpp:674
AnalyticalDirichletBC::destroyProblem
MoFEMErrorCode destroyProblem()
Destroy problem.
Definition: AnalyticalDirichlet.cpp:268
ThermalElement::addThermalElements
MoFEMErrorCode addThermalElements(const std::string field_name, const std::string mesh_nodals_positions="MESH_NODE_POSITIONS")
add thermal element on tets
Definition: ThermalElement.cpp:527
help
static char help[]
Definition: thermal_with_analytical_bc.cpp:32
MoFEM.hpp
A
constexpr AssemblyType A
Definition: operators_tests.cpp:30
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
AnalyticalDirichlet.hpp
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
ThermalElement::getLoopFeFlux
MyTriFE & getLoopFeFlux()
Definition: ThermalElement.hpp:67
ROW
@ ROW
Definition: definitions.h:136
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:2010
AnalyticalDirichletBC::setUpProblem
MoFEMErrorCode setUpProblem(MoFEM::Interface &m_field, string problem)
set problem solver and create matrices and vectors
Definition: AnalyticalDirichlet.cpp:208
MoFEM::CoreInterface::problem_basic_method_preProcess
virtual MoFEMErrorCode problem_basic_method_preProcess(const Problem *problem_ptr, BasicMethod &method, int verb=DEFAULT_VERBOSITY)=0
Set data for BasicMethod.
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
ThermalElement::setThermalFluxFiniteElementRhsOperators
MoFEMErrorCode setThermalFluxFiniteElementRhsOperators(string field_name, Vec &F, const std::string mesh_nodals_positions="MESH_NODE_POSITIONS")
this function is used in case of stationary problem for heat flux terms
Definition: ThermalElement.cpp:724
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
MoFEM::CoreInterface::build_finite_elements
virtual MoFEMErrorCode build_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elements.
AnalyticalDirichletBC::DirichletBC
Structure used to enforce analytical boundary conditions.
Definition: AnalyticalDirichlet.hpp:124
ThermalElement.hpp
Operators and data structures for thermal analysis.
AnalyticalDirichletBC::setApproxOps
MoFEMErrorCode setApproxOps(MoFEM::Interface &m_field, const std::string field_name, boost::shared_ptr< FUNEVAL > function_evaluator, const int field_number=0, const string nodals_positions="MESH_NODE_POSITIONS")
Set operators used to calculate the rhs vector and the lhs matrix.
Definition: AnalyticalDirichlet.hpp:158
AnalyticalDirichletBC::solveProblem
MoFEMErrorCode solveProblem(MoFEM::Interface &m_field, string problem, string fe, DirichletBC &bc, Range &tris)
solve boundary problem
Definition: AnalyticalDirichlet.cpp:222
PostProcTemplateOnRefineMesh::writeFile
MoFEMErrorCode writeFile(const std::string file_name, const char *file_type="MOAB", const char *file_options="PARALLEL=WRITE_PART")
wrote results in (MOAB) format, use "file_name.h5m"
Definition: PostProcOnRefMesh.hpp:253
debug
static int debug
Definition: thermal_with_analytical_bc.cpp:31
Projection10NodeCoordsOnField.hpp
ThermalElement::setThermalFiniteElementRhsOperators
MoFEMErrorCode setThermalFiniteElementRhsOperators(string field_name, Vec &F)
this function is used in case of stationary problem to set elements for rhs
Definition: ThermalElement.cpp:699
PostProcTemplateVolumeOnRefinedMesh::generateReferenceElementMesh
MoFEMErrorCode generateReferenceElementMesh()
Generate reference mesh on single element.
Definition: PostProcOnRefMesh.hpp:301
PostProcTemplateOnRefineMesh::addFieldValuesGradientPostProc
MoFEMErrorCode addFieldValuesGradientPostProc(const std::string field_name, Vec v=PETSC_NULL)
Add operator to post-process L2 or H1 field gradient.
Definition: PostProcOnRefMesh.hpp:195
PostProcOnRefMesh.hpp
Post-process fields on refined mesh.
AnalyticalFunction::val
std::vector< VectorDouble > val
Definition: thermal_with_analytical_bc.cpp:36
main
int main(int argc, char *argv[])
Definition: thermal_with_analytical_bc.cpp:46
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::VecManager
Vector manager is used to create vectors \mofem_vectors.
Definition: VecManager.hpp:23
mesh_file_name
char mesh_file_name[255]
Definition: mesh_smoothing.cpp:22
MoFEM::MatrixManager
Matrix manager is used to build and partition problems.
Definition: MatrixManager.hpp:21
Range
ThermalElement
structure grouping operators and data used for thermal problems
Definition: ThermalElement.hpp:27
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:385
PostProcVolumeOnRefinedMesh
Post processing.
Definition: PostProcOnRefMesh.hpp:955
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
AINSWORTH_LEGENDRE_BASE
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
EigenMatrix::Vec
const FTensor::Tensor2< T, Dim, Dim > Vec
Definition: MatrixFunction.hpp:66
AnalyticalFunction
Definition: thermal_with_analytical_bc.cpp:34
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
_IT_CUBITMESHSETS_BY_NAME_FOR_LOOP_
#define _IT_CUBITMESHSETS_BY_NAME_FOR_LOOP_(MESHSET_MANAGER, NAME, IT)
Iterator that loops over Cubit BlockSet having a particular name.
Definition: MeshsetsManager.hpp:94
AnalyticalDirichletBC
Analytical Dirichlet boundary conditions.
Definition: AnalyticalDirichlet.hpp:18
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_ATOM_TEST_INVALID
@ MOFEM_ATOM_TEST_INVALID
Definition: definitions.h:40
MoFEM::CoreInterface::build_adjacencies
virtual MoFEMErrorCode build_adjacencies(const Range &ents, int verb=DEFAULT_VERBOSITY)=0
build adjacencies
MoFEM::CoreInterface::check_finite_element
virtual bool check_finite_element(const std::string &name) const =0
Check if finite element is in database.
ThermalElement::setThermalFiniteElementLhsOperators
MoFEMErrorCode setThermalFiniteElementLhsOperators(string field_name, Mat A)
this function is used in case of stationary heat conductivity problem for lhs
Definition: ThermalElement.cpp:713
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
AnalyticalFunction::operator()
std::vector< VectorDouble > & operator()(double x, double y, double z)
Definition: thermal_with_analytical_bc.cpp:38
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::PetscOptionsGetInt
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
Definition: DeprecatedPetsc.hpp:142
DirichletBC.hpp
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.
F
@ F
Definition: free_surface.cpp:394
PostProcTemplateOnRefineMesh::addFieldValuesPostProc
MoFEMErrorCode addFieldValuesPostProc(const std::string field_name, Vec v=PETSC_NULL)
Add operator to post-process L2, H1, Hdiv, Hcurl field value.
Definition: PostProcOnRefMesh.hpp:153