v0.14.0
thermal_elem.cpp
Go to the documentation of this file.
1 
2 
4 using namespace MoFEM;
5 
6 namespace bio = boost::iostreams;
7 using bio::stream;
8 using bio::tee_device;
9 
10 static char help[] = "...\n\n";
11 
12 int main(int argc, char *argv[]) {
13 
14  MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
15 
16  try {
17 
18  moab::Core mb_instance;
19  moab::Interface &moab = mb_instance;
20  int rank;
21  MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
22 
23  PetscBool flg = PETSC_TRUE;
24  char mesh_file_name[255];
25  CHKERR PetscOptionsGetString(PETSC_NULL, PETSC_NULL, "-my_file",
26  mesh_file_name, 255, &flg);
27  if (flg != PETSC_TRUE) {
28  SETERRQ(PETSC_COMM_SELF, 1, "*** ERROR -my_file (MESH FILE NEEDED)");
29  }
30 
31  const char *option;
32  option = "";
33  CHKERR moab.load_file(mesh_file_name, 0, option);
34 
35  // Create MoFEM (Joseph) database
36  MoFEM::Core core(moab);
37  MoFEM::Interface &m_field = core;
38 
39  // set entitities bit level
40  BitRefLevel bit_level0;
41  bit_level0.set(0);
42  EntityHandle meshset_level0;
43  CHKERR moab.create_meshset(MESHSET_SET, meshset_level0);
44  CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(
45  0, 3, bit_level0);
46 
47  // Fields
48  CHKERR m_field.add_field("TEMP", H1, AINSWORTH_LEGENDRE_BASE, 1);
49 
50  // Problem
51  CHKERR m_field.add_problem("TEST_PROBLEM");
52 
53  // set refinement level for problem
54  CHKERR m_field.modify_problem_ref_level_add_bit("TEST_PROBLEM", bit_level0);
55 
56  // meshset consisting all entities in mesh
57  EntityHandle root_set = moab.get_root_set();
58  // add entities to field
59  CHKERR m_field.add_ents_to_field_by_type(root_set, MBTET, "TEMP");
60 
61  // set app. order
62  // see Hierarchic Finite Element Bases on Unstructured Tetrahedral Meshes
63  // (Mark Ainsworth & Joe Coyle)
64  int order = 4;
65  CHKERR m_field.set_field_order(root_set, MBTET, "TEMP", order);
66  CHKERR m_field.set_field_order(root_set, MBTRI, "TEMP", order);
67  CHKERR m_field.set_field_order(root_set, MBEDGE, "TEMP", order);
68  CHKERR m_field.set_field_order(root_set, MBVERTEX, "TEMP", 1);
69 
70  ThermalElement thermal_elements(m_field);
71  CHKERR thermal_elements.addThermalElements("TEMP");
72  CHKERR thermal_elements.addThermalFluxElement("TEMP");
73 
74  CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM",
75  "THERMAL_FE");
76  CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM",
77  "THERMAL_FLUX_FE");
78 
79  /****/
80  // build database
81  // build field
82  CHKERR m_field.build_fields();
83  // build finite elemnts
84  CHKERR m_field.build_finite_elements();
85  // build adjacencies
86  CHKERR m_field.build_adjacencies(bit_level0);
87 
88  /****/
89  // mesh partitioning
90  // build problem
91  ProblemsManager *prb_mng_ptr;
92  CHKERR m_field.getInterface(prb_mng_ptr);
93  CHKERR prb_mng_ptr->buildProblem("TEST_PROBLEM", true);
94  // partition
95  CHKERR prb_mng_ptr->partitionSimpleProblem("TEST_PROBLEM");
96  CHKERR prb_mng_ptr->partitionFiniteElements("TEST_PROBLEM");
97  // what are ghost nodes, see Petsc Manual
98  CHKERR prb_mng_ptr->partitionGhostDofs("TEST_PROBLEM");
99 
100  Vec F;
101  CHKERR m_field.getInterface<VecManager>()->vecCreateGhost("TEST_PROBLEM",
102  ROW, &F);
103  Vec T;
104  CHKERR VecDuplicate(F, &T);
105  Mat A;
107  ->createMPIAIJWithArrays<PetscGlobalIdx_mi_tag>("TEST_PROBLEM", &A);
108 
109  DirichletTemperatureBc my_dirichlet_bc(m_field, "TEMP", A, T, F);
110  CHKERR thermal_elements.setThermalFiniteElementRhsOperators("TEMP", F);
111  CHKERR thermal_elements.setThermalFiniteElementLhsOperators("TEMP", A);
112  CHKERR thermal_elements.setThermalFluxFiniteElementRhsOperators("TEMP", F);
113 
114  CHKERR VecZeroEntries(T);
115  CHKERR VecZeroEntries(F);
116  CHKERR MatZeroEntries(A);
117 
118  // preproc
119  CHKERR m_field.problem_basic_method_preProcess("TEST_PROBLEM",
120  my_dirichlet_bc);
121  CHKERR m_field.getInterface<VecManager>()->setGlobalGhostVector(
122  "TEST_PROBLEM", ROW, T, INSERT_VALUES, SCATTER_REVERSE);
123 
124  CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "THERMAL_FE",
125  thermal_elements.getLoopFeRhs());
126  CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "THERMAL_FE",
127  thermal_elements.getLoopFeLhs());
128  CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "THERMAL_FLUX_FE",
129  thermal_elements.getLoopFeFlux());
130 
131  // postproc
132  CHKERR m_field.problem_basic_method_postProcess("TEST_PROBLEM",
133  my_dirichlet_bc);
134 
135  CHKERR VecGhostUpdateBegin(F, ADD_VALUES, SCATTER_REVERSE);
136  CHKERR VecGhostUpdateEnd(F, ADD_VALUES, SCATTER_REVERSE);
137  CHKERR VecAssemblyBegin(F);
138  CHKERR VecAssemblyEnd(F);
139  CHKERR MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);
140  CHKERR MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);
141 
142  CHKERR VecScale(F, -1);
143 
144  // Solver
145  KSP solver;
146  CHKERR KSPCreate(PETSC_COMM_WORLD, &solver);
147  CHKERR KSPSetOperators(solver, A, A);
148  CHKERR KSPSetFromOptions(solver);
149  CHKERR KSPSetUp(solver);
150 
151  CHKERR KSPSolve(solver, F, T);
152  CHKERR VecGhostUpdateBegin(T, INSERT_VALUES, SCATTER_FORWARD);
153  CHKERR VecGhostUpdateEnd(T, INSERT_VALUES, SCATTER_FORWARD);
154 
155  CHKERR m_field.problem_basic_method_preProcess("TEST_PROBLEM",
156  my_dirichlet_bc);
157 
158  // Save data on mesh
159  CHKERR m_field.getInterface<VecManager>()->setGlobalGhostVector(
160  "TEST_PROBLEM", ROW, T, INSERT_VALUES, SCATTER_REVERSE);
161 
162  double sum = 0;
163  CHKERR VecSum(F, &sum);
164  CHKERR PetscPrintf(PETSC_COMM_WORLD, "sum = %9.8f\n", sum);
165  double fnorm;
166  CHKERR VecNorm(F, NORM_2, &fnorm);
167  CHKERR PetscPrintf(PETSC_COMM_WORLD, "fnorm = %9.8e\n", fnorm);
168  if (fabs(sum + 0.59583333) > 1e-7) {
169  SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID, "Failed to pass test");
170  }
171  if (fabs(fnorm - 2.32872499e-01) > 1e-6) {
172  SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID, "Failed to pass test");
173  }
174 
175  CHKERR MatDestroy(&A);
176  CHKERR VecDestroy(&F);
177  CHKERR VecDestroy(&T);
178  CHKERR KSPDestroy(&solver);
179 
180  }
181  CATCH_ERRORS;
182 
184 
185  return 0;
186 }
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Definition: UnknownInterface.hpp:93
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::ProblemsManager::buildProblem
MoFEMErrorCode buildProblem(const std::string name, const bool square_matrix, int verb=VERBOSE)
build problem data structures
Definition: ProblemsManager.cpp:87
help
static char help[]
Definition: thermal_elem.cpp:10
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
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
DirichletTemperatureBc
Definition: DirichletBC.hpp:239
A
constexpr AssemblyType A
Definition: operators_tests.cpp:30
BasicFiniteElements.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::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
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.
ThermalElement::addThermalFluxElement
MoFEMErrorCode addThermalFluxElement(const std::string field_name, const std::string mesh_nodals_positions="MESH_NODE_POSITIONS")
add heat flux element
Definition: ThermalElement.cpp:574
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
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
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
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1148
main
int main(int argc, char *argv[])
Definition: thermal_elem.cpp:12
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
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::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
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
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::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