|
| v0.14.0
|
Example of steady thermal analysis.
More...
Go to the source code of this file.
|
int | main (int argc, char *argv[]) |
|
|
static char | help [] = "...\n\n" |
|
Example of steady thermal analysis.
TODO:
- Todo:
- Make it work in distributed meshes with multigird solver. At the moment it is not working efficient as can.
Definition in file thermal_steady.cpp.
◆ main()
int main |
( |
int |
argc, |
|
|
char * |
argv[] |
|
) |
| |
Definition at line 21 of file thermal_steady.cpp.
23 const string default_options =
"-ksp_type fgmres \n"
25 "-pc_factor_mat_solver_type mumps \n"
26 "-mat_mumps_icntl_20 0 \n"
29 string param_file =
"param_file.petsc";
30 if (!
static_cast<bool>(ifstream(param_file))) {
31 std::ofstream file(param_file.c_str(), std::ios::ate);
33 file << default_options;
45 MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
47 PetscBool flg = PETSC_TRUE;
51 if (flg != PETSC_TRUE) {
53 "*** ERROR -my_file (MESH FILE NEEDED)");
68 CHKERR moab.create_meshset(MESHSET_SET, meshset_level0);
93 if (flg != PETSC_TRUE) {
105 "MESH_NODE_POSITIONS");
112 CHKERR thermal_elements.addThermalElements(
"TEMP");
113 CHKERR thermal_elements.addThermalFluxElement(
"TEMP");
114 CHKERR thermal_elements.addThermalConvectionElement(
"TEMP");
121 "THERMAL_CONVECTION_FE");
138 "MESH_NODE_POSITIONS");
156 ->createMPIAIJWithArrays<PetscGlobalIdx_mi_tag>(
"THERMAL_PROBLEM", &
A);
161 true,
false,
false,
false);
163 true,
false,
false,
false);
165 CHKERR thermal_elements.setThermalFiniteElementRhsOperators(
"TEMP",
F);
166 CHKERR thermal_elements.setThermalFiniteElementLhsOperators(
"TEMP",
A);
167 CHKERR thermal_elements.setThermalFluxFiniteElementRhsOperators(
"TEMP",
F);
168 CHKERR thermal_elements.setThermalConvectionFiniteElementLhsOperators(
170 CHKERR thermal_elements.setThermalConvectionFiniteElementRhsOperators(
174 CHKERR VecGhostUpdateBegin(T, INSERT_VALUES, SCATTER_FORWARD);
175 CHKERR VecGhostUpdateEnd(T, INSERT_VALUES, SCATTER_FORWARD);
177 CHKERR VecGhostUpdateBegin(
F, INSERT_VALUES, SCATTER_FORWARD);
178 CHKERR VecGhostUpdateEnd(
F, INSERT_VALUES, SCATTER_FORWARD);
185 "THERMAL_PROBLEM",
ROW, T, INSERT_VALUES, SCATTER_REVERSE);
188 thermal_elements.getLoopFeRhs());
190 thermal_elements.getLoopFeLhs());
192 thermal_elements.getLoopFeFlux());
194 "THERMAL_PROBLEM",
"THERMAL_CONVECTION_FE",
195 thermal_elements.getLoopFeConvectionRhs());
197 "THERMAL_PROBLEM",
"THERMAL_CONVECTION_FE",
198 thermal_elements.getLoopFeConvectionLhs());
204 CHKERR VecGhostUpdateBegin(
F, ADD_VALUES, SCATTER_REVERSE);
205 CHKERR VecGhostUpdateEnd(
F, ADD_VALUES, SCATTER_REVERSE);
208 CHKERR MatAssemblyBegin(
A, MAT_FINAL_ASSEMBLY);
209 CHKERR MatAssemblyEnd(
A, MAT_FINAL_ASSEMBLY);
215 CHKERR KSPCreate(PETSC_COMM_WORLD, &solver);
216 CHKERR KSPSetOperators(solver,
A,
A);
217 CHKERR KSPSetFromOptions(solver);
221 CHKERR VecGhostUpdateBegin(T, INSERT_VALUES, SCATTER_FORWARD);
222 CHKERR VecGhostUpdateEnd(T, INSERT_VALUES, SCATTER_FORWARD);
229 "THERMAL_PROBLEM",
ROW, T, INSERT_VALUES, SCATTER_REVERSE);
232 CHKERR moab.write_file(
"solution.h5m");
238 ent_method_on_10nodeTet.setNodes =
false;
243 CHKERR moab.create_meshset(MESHSET_SET, out_meshset);
245 "THERMAL_PROBLEM",
"THERMAL_FE", out_meshset);
246 CHKERR moab.write_file(
"out.vtk",
"VTK",
"", &out_meshset, 1);
247 CHKERR moab.delete_entities(&out_meshset, 1);
253 CHKERR KSPDestroy(&solver);
◆ help
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
virtual MoFEMErrorCode problem_basic_method_postProcess(const Problem *problem_ptr, BasicMethod &method, int verb=DEFAULT_VERBOSITY)=0
Set data for BasicMethod.
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.
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.
MoFEMErrorCode buildProblem(const std::string name, const bool square_matrix, int verb=VERBOSE)
build problem data structures
Problem manager is used to build and partition problems.
MoFEMErrorCode addHOOpsVol(const std::string field, E &e, bool h1, bool hcurl, bool hdiv, bool l2)
virtual int get_comm_rank() const =0
virtual MoFEMErrorCode get_problem_finite_elements_entities(const std::string name, const std::string &fe_name, const EntityHandle meshset)=0
add finite elements to the meshset
Projection of edge entities with one mid-node on hierarchical basis.
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
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.
Deprecated interface functions.
DeprecatedCoreInterface Interface
virtual MoFEMErrorCode problem_basic_method_preProcess(const Problem *problem_ptr, BasicMethod &method, int verb=DEFAULT_VERBOSITY)=0
Set data for BasicMethod.
#define CHKERR
Inline error check.
virtual MoFEMErrorCode build_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elements.
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 modify_problem_ref_level_add_bit(const std::string &name_problem, const BitRefLevel &bit)=0
add ref level to problem
Vector manager is used to create vectors \mofem_vectors.
Matrix manager is used to build and partition problems.
structure grouping operators and data used for thermal problems
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
#define CATCH_ERRORS
Catch errors.
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
const FTensor::Tensor2< T, Dim, Dim > Vec
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
virtual MoFEMErrorCode build_fields(int verb=DEFAULT_VERBOSITY)=0
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
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
virtual MoFEMErrorCode build_adjacencies(const Range &ents, int verb=DEFAULT_VERBOSITY)=0
build adjacencies
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.
MoFEMErrorCode partitionGhostDofs(const std::string name, int verb=VERBOSE)
determine ghost nodes
virtual MoFEMErrorCode add_problem(const std::string &name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add problem.
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
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.
MoFEMErrorCode partitionProblem(const std::string name, int verb=VERBOSE)
partition problem dofs (collective)