v0.14.0
Loading...
Searching...
No Matches
thermal_elem_unsteady.cpp File Reference

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 12 of file thermal_elem_unsteady.cpp.

12 {
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
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
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 }
231
233
234 return 0;
235}
@ ROW
#define CATCH_ERRORS
Catch errors.
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base nme:nme847.
Definition definitions.h:60
@ H1
continuous field
Definition definitions.h:85
@ MOFEM_ATOM_TEST_INVALID
Definition definitions.h:40
#define CHKERR
Inline error check.
constexpr int order
@ F
virtual MoFEMErrorCode build_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elements.
virtual MoFEMErrorCode modify_finite_element_add_field_data(const std::string &fe_name, const std::string name_field)=0
set finite element field data
virtual MoFEMErrorCode build_fields(int verb=DEFAULT_VERBOSITY)=0
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.
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.
MoFEMErrorCode partitionGhostDofs(const std::string name, int verb=VERBOSE)
determine ghost nodes
MoFEMErrorCode partitionSimpleProblem(const std::string name, int verb=VERBOSE)
partition problem dofs
MoFEMErrorCode buildProblem(const std::string name, const bool square_matrix, int verb=VERBOSE)
build problem data structures
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
virtual MoFEMErrorCode add_problem(const std::string &name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add problem.
virtual MoFEMErrorCode modify_problem_ref_level_add_bit(const std::string &name_problem, const BitRefLevel &bit)=0
add ref level to problem
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
virtual MoFEMErrorCode load_series_data(const std::string &serie_name, const int step_number)
virtual MoFEMErrorCode finalize_series_recorder(const std::string &serie_name)
virtual MoFEMErrorCode initialize_series_recorder(const std::string &serie_name)
#define _IT_SERIES_STEPS_BY_NAME_FOR_LOOP_(RECORDER, NAME, IT)
loop over recorded series step
virtual MoFEMErrorCode add_series_recorder(const std::string &series_name)
MoFEM::TsCtx * ts_ctx
char mesh_file_name[255]
static MoFEMErrorCodeGeneric< moab::ErrorCode > rval
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition Types.hpp:40
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
PetscErrorCode TsMonitorSet(TS ts, PetscInt step, PetscReal t, Vec u, void *ctx)
Set monitor for TS solver.
Definition TsCtx.cpp:259
PetscErrorCode TsSetIFunction(TS ts, PetscReal t, Vec u, Vec u_t, Vec F, void *ctx)
Set IFunction for TS solver.
Definition TsCtx.cpp:56
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
constexpr AssemblyType A
Managing BitRefLevels.
virtual MoFEMErrorCode build_adjacencies(const Range &ents, int verb=DEFAULT_VERBOSITY)=0
build adjacencies
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.
Core (interface) class.
Definition Core.hpp:82
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
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition Core.cpp:112
Deprecated interface functions.
Matrix manager is used to build and partition problems.
Problem manager is used to build and partition problems.
Interface for Time Stepping (TS) solver.
Definition TsCtx.hpp:17
BasicMethodsSequence & getPostProcessIJacobian()
Get the postProcess to do IJacobian object.
Definition TsCtx.hpp:128
BasicMethodsSequence & getPostProcessIFunction()
Get the postProcess to do IFunction object.
Definition TsCtx.hpp:112
BasicMethodsSequence & getPreProcessIFunction()
Get the preProcess to do IFunction object.
Definition TsCtx.hpp:105
BasicMethodsSequence & getPostProcessMonitor()
Get the postProcess to do Monitor object.
Definition TsCtx.hpp:144
BasicMethodsSequence & getPreProcessIJacobian()
Get the preProcess to do IJacobian object.
Definition TsCtx.hpp:121
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Vector manager is used to create vectors \mofem_vectors.
TS monitore it records temperature at time steps.
this calass is to control time stepping
structure grouping operators and data used for thermal problems
static char help[]

Variable Documentation

◆ help

char help[] = "...\n\n"
static

Definition at line 10 of file thermal_elem_unsteady.cpp.