v0.14.0
thermal_elem_unsteady.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  rval = 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  rval = 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  CHKERR m_field.add_field("TEMP_RATE", H1, AINSWORTH_LEGENDRE_BASE, 1);
50 
51  // Problem
52  CHKERR m_field.add_problem("TEST_PROBLEM");
53 
54  // set refinement level for problem
55  CHKERR m_field.modify_problem_ref_level_add_bit("TEST_PROBLEM", bit_level0);
56 
57  // meshset consisting all entities in mesh
58  EntityHandle root_set = moab.get_root_set();
59  // add entities to field
60  CHKERR m_field.add_ents_to_field_by_type(root_set, MBTET, "TEMP");
61  CHKERR m_field.add_ents_to_field_by_type(root_set, MBTET, "TEMP_RATE");
62 
63  // set app. order
64  // see Hierarchic Finite Element Bases on Unstructured Tetrahedral Meshes
65  // (Mark Ainsworth & Joe Coyle)
66  int order = 2;
67  CHKERR m_field.set_field_order(root_set, MBTET, "TEMP", order);
68  CHKERR m_field.set_field_order(root_set, MBTRI, "TEMP", order);
69  CHKERR m_field.set_field_order(root_set, MBEDGE, "TEMP", order);
70  CHKERR m_field.set_field_order(root_set, MBVERTEX, "TEMP", 1);
71 
72  CHKERR m_field.set_field_order(root_set, MBTET, "TEMP_RATE", order);
73  CHKERR m_field.set_field_order(root_set, MBTRI, "TEMP_RATE", order);
74  CHKERR m_field.set_field_order(root_set, MBEDGE, "TEMP_RATE", order);
75  CHKERR m_field.set_field_order(root_set, MBVERTEX, "TEMP_RATE", 1);
76 
77  ThermalElement thermal_elements(m_field);
78  CHKERR thermal_elements.addThermalElements("TEMP");
79  CHKERR thermal_elements.addThermalFluxElement("TEMP");
80  // add rate of temerature to data field of finite element
81  CHKERR m_field.modify_finite_element_add_field_data("THERMAL_FE",
82  "TEMP_RATE");
83 
84  CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM",
85  "THERMAL_FE");
86  CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM",
87  "THERMAL_FLUX_FE");
88 
89  /****/
90  // build database
91  // build field
92  CHKERR m_field.build_fields();
93  // build finite elemnts
94  CHKERR m_field.build_finite_elements();
95  // build adjacencies
96  CHKERR m_field.build_adjacencies(bit_level0);
97  // build problem
98  ProblemsManager *prb_mng_ptr;
99  CHKERR m_field.getInterface(prb_mng_ptr);
100  CHKERR prb_mng_ptr->buildProblem("TEST_PROBLEM", true);
101 
102  /****/
103  // mesh partitioning
104  // partition
105  CHKERR prb_mng_ptr->partitionSimpleProblem("TEST_PROBLEM");
106  CHKERR prb_mng_ptr->partitionFiniteElements("TEST_PROBLEM");
107  // what are ghost nodes, see Petsc Manual
108  CHKERR prb_mng_ptr->partitionGhostDofs("TEST_PROBLEM");
109 
110  Vec F;
111  CHKERR m_field.getInterface<VecManager>()->vecCreateGhost("TEST_PROBLEM",
112  ROW, &F);
113  Vec T;
114  CHKERR VecDuplicate(F, &T);
115  Mat A;
117  ->createMPIAIJWithArrays<PetscGlobalIdx_mi_tag>("TEST_PROBLEM", &A);
118 
119  // TS
120  TsCtx ts_ctx(m_field, "TEST_PROBLEM");
121  TS ts;
122  CHKERR TSCreate(PETSC_COMM_WORLD, &ts);
123  CHKERR TSSetType(ts, TSBEULER);
124 
125  DirichletTemperatureBc my_dirichlet_bc(m_field, "TEMP", A, T, F);
126  ThermalElement::UpdateAndControl update_velocities(m_field, "TEMP",
127  "TEMP_RATE");
128  ThermalElement::TimeSeriesMonitor monitor(m_field, "THEMP_SERIES", "TEMP");
129 
130  // preprocess
131  ts_ctx.getPreProcessIFunction().push_back(&update_velocities);
132  ts_ctx.getPreProcessIFunction().push_back(&my_dirichlet_bc);
133  ts_ctx.getPreProcessIJacobian().push_back(&my_dirichlet_bc);
134 
135  // and temperature element functions
136  CHKERR thermal_elements.setTimeSteppingProblem(ts_ctx, "TEMP", "TEMP_RATE");
137 
138  // postprocess
139  ts_ctx.getPostProcessIFunction().push_back(&my_dirichlet_bc);
140  ts_ctx.getPostProcessIJacobian().push_back(&my_dirichlet_bc);
141  ts_ctx.getPostProcessMonitor().push_back(&monitor);
142 
143  CHKERR TSSetIFunction(ts, F, TsSetIFunction, &ts_ctx);
144  CHKERR TSSetIJacobian(ts, A, A, TsSetIJacobian, &ts_ctx);
145  CHKERR TSMonitorSet(ts, TsMonitorSet, &ts_ctx, PETSC_NULL);
146 
147  double ftime = 1;
148  CHKERR TSSetDuration(ts, PETSC_DEFAULT, ftime);
149  CHKERR TSSetSolution(ts, T);
150  CHKERR TSSetFromOptions(ts);
151 
152  SeriesRecorder *recorder_ptr;
153  CHKERR m_field.getInterface(recorder_ptr);
154  CHKERR recorder_ptr->add_series_recorder("THEMP_SERIES");
155  CHKERR recorder_ptr->initialize_series_recorder("THEMP_SERIES");
156 
157 #if PETSC_VERSION_GE(3, 7, 0)
158  CHKERR TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER);
159 #endif
160  CHKERR TSSolve(ts, T);
161  CHKERR TSGetTime(ts, &ftime);
162 
163  CHKERR recorder_ptr->finalize_series_recorder("THEMP_SERIES");
164 
165  PetscInt steps, snesfails, rejects, nonlinits, linits;
166  CHKERR TSGetTimeStepNumber(ts, &steps);
167  CHKERR TSGetSNESFailures(ts, &snesfails);
168  CHKERR TSGetStepRejections(ts, &rejects);
169  CHKERR TSGetSNESIterations(ts, &nonlinits);
170  CHKERR TSGetKSPIterations(ts, &linits);
171  PetscPrintf(PETSC_COMM_WORLD,
172  "steps %D (%D rejected, %D SNES fails), ftime %g, nonlinits "
173  "%D, linits %D\n",
174  steps, rejects, snesfails, ftime, nonlinits, linits);
175 
176  // PetscViewer viewer;
177  // PetscViewerASCIIOpen(PETSC_COMM_WORLD,"thermal_elem_unsteady.txt",&viewer);
178 
179  double sum = 0;
180  double fnorm = 0;
181 
182  for (_IT_SERIES_STEPS_BY_NAME_FOR_LOOP_(recorder_ptr, "THEMP_SERIES",
183  sit)) {
184 
185  CHKERR recorder_ptr->load_series_data("THEMP_SERIES",
186  sit->get_step_number());
187  CHKERR m_field.getInterface<VecManager>()->setLocalGhostVector(
188  "TEST_PROBLEM", ROW, T, INSERT_VALUES, SCATTER_FORWARD);
189 
190  double sum0;
191  CHKERR VecSum(T, &sum0);
192  CHKERR PetscPrintf(PETSC_COMM_WORLD, "sum0 = %9.8e\n", sum0);
193  double fnorm0;
194  CHKERR VecNorm(T, NORM_2, &fnorm0);
195  CHKERR PetscPrintf(PETSC_COMM_WORLD, "fnorm0 = %9.8e\n", fnorm0);
196 
197  sum += sum0;
198  fnorm += fnorm0;
199 
200  // CHKERR VecChop(T,1e-4);
201  // CHKERR VecView(T,viewer);
202  }
203  CHKERR PetscPrintf(PETSC_COMM_WORLD, "sum = %9.8e\n", sum);
204  CHKERR PetscPrintf(PETSC_COMM_WORLD, "fnorm = %9.8e\n", fnorm);
205  if (fabs(sum + 1.32314077e+01) > 1e-7) {
206  SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID, "Failed to pass test");
207  }
208  if (fabs(fnorm - 4.59664623e+00) > 1e-6) {
209  SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID, "Failed to pass test");
210  }
211 
212  // CHKERR PetscViewerDestroy(&viewer);
213 
214  /*PostProcVertexMethod ent_method(moab,"TEMP");
215  CHKERR m_field.loop_dofs("TEST_PROBLEM","TEMP",ROW,ent_method);
216  if(pcomm->rank()==0) {
217  EntityHandle out_meshset;
218  rval = moab.create_meshset(MESHSET_SET,out_meshset);
219  CHKERR
220  m_field.get_problem_finite_elements_entities("TEST_PROBLEM","THERMAL_FE",out_meshset);
221  rval = moab.write_file("out.vtk","VTK","",&out_meshset,1);
222  rval = moab.delete_entities(&out_meshset,1);
223  }*/
224 
225  CHKERR TSDestroy(&ts);
226  CHKERR MatDestroy(&A);
227  CHKERR VecDestroy(&F);
228  CHKERR VecDestroy(&T);
229  }
230  CATCH_ERRORS;
231 
233 
234  return 0;
235 }
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Definition: UnknownInterface.hpp:93
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::TsSetIFunction
PetscErrorCode TsSetIFunction(TS ts, PetscReal t, Vec u, Vec u_t, Vec F, void *ctx)
Set IFunction for TS solver.
Definition: TsCtx.cpp:56
MoFEM::SeriesRecorder
Definition: SeriesRecorder.hpp:25
MoFEM::TsCtx::getPostProcessMonitor
BasicMethodsSequence & getPostProcessMonitor()
Get the postProcess to do Monitor object.
Definition: TsCtx.hpp:144
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
MoFEM::TsCtx::getPostProcessIJacobian
BasicMethodsSequence & getPostProcessIJacobian()
Get the postProcess to do IJacobian object.
Definition: TsCtx.hpp:128
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
ts_ctx
MoFEM::TsCtx * ts_ctx
Definition: level_set.cpp:1932
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
ROW
@ ROW
Definition: definitions.h:136
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:2010
MoFEM::TsMonitorSet
PetscErrorCode TsMonitorSet(TS ts, PetscInt step, PetscReal t, Vec u, void *ctx)
Set monitor for TS solver.
Definition: TsCtx.cpp:259
MoFEM::SeriesRecorder::add_series_recorder
virtual MoFEMErrorCode add_series_recorder(const std::string &series_name)
Definition: SeriesRecorder.cpp:90
MoFEM::Exceptions::rval
static MoFEMErrorCodeGeneric< moab::ErrorCode > rval
Definition: Exceptions.hpp:74
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
help
static char help[]
Definition: thermal_elem_unsteady.cpp:10
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.
MoFEM::TsCtx::getPreProcessIJacobian
BasicMethodsSequence & getPreProcessIJacobian()
Get the preProcess to do IJacobian object.
Definition: TsCtx.hpp:121
main
int main(int argc, char *argv[])
Definition: thermal_elem_unsteady.cpp:12
MoFEM::SeriesRecorder::finalize_series_recorder
virtual MoFEMErrorCode finalize_series_recorder(const std::string &serie_name)
Definition: SeriesRecorder.cpp:264
ThermalElement::setTimeSteppingProblem
MoFEMErrorCode setTimeSteppingProblem(string field_name, string rate_name, const std::string mesh_nodals_positions="MESH_NODE_POSITIONS")
set up operators for unsteady heat flux; convection; radiation problem
Definition: ThermalElement.cpp:774
MoFEM::TsCtx
Interface for Time Stepping (TS) solver.
Definition: TsCtx.hpp:17
MoFEM::TsCtx::getPostProcessIFunction
BasicMethodsSequence & getPostProcessIFunction()
Get the postProcess to do IFunction object.
Definition: TsCtx.hpp:112
MoFEM::SeriesRecorder::load_series_data
virtual MoFEMErrorCode load_series_data(const std::string &serie_name, const int step_number)
Definition: SeriesRecorder.cpp:309
_IT_SERIES_STEPS_BY_NAME_FOR_LOOP_
#define _IT_SERIES_STEPS_BY_NAME_FOR_LOOP_(RECORDER, NAME, IT)
loop over recorded series step
Definition: SeriesRecorder.hpp:205
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
MoFEM::TsSetIJacobian
PetscErrorCode TsSetIJacobian(TS ts, PetscReal t, Vec u, Vec u_t, PetscReal a, Mat A, Mat B, void *ctx)
Set function evaluating jacobian in TS solver.
Definition: TsCtx.cpp:165
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
ThermalElement::TimeSeriesMonitor
TS monitore it records temperature at time steps.
Definition: ThermalElement.hpp:579
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::TsCtx::getPreProcessIFunction
BasicMethodsSequence & getPreProcessIFunction()
Get the preProcess to do IFunction object.
Definition: TsCtx.hpp:105
MoFEM::MatrixManager
Matrix manager is used to build and partition problems.
Definition: MatrixManager.hpp:21
ThermalElement::UpdateAndControl
this calass is to control time stepping
Definition: ThermalElement.hpp:562
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
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
MoFEM::SeriesRecorder::initialize_series_recorder
virtual MoFEMErrorCode initialize_series_recorder(const std::string &serie_name)
Definition: SeriesRecorder.cpp:246
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