v0.14.0
Loading...
Searching...
No Matches
thermal_stress_elem.cpp
Go to the documentation of this file.
1
2
4using namespace MoFEM;
5
6namespace bio = boost::iostreams;
7using bio::stream;
8using bio::tee_device;
9
10static char help[] = "...\n\n";
11
12int 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("DISP", H1, AINSWORTH_LEGENDRE_BASE, 3);
49 CHKERR m_field.add_field("TEMP", H1, AINSWORTH_LEGENDRE_BASE, 1);
50
51 // Problem
52 CHKERR m_field.add_problem("PROB");
53
54 // set refinement level for problem
55 CHKERR m_field.modify_problem_ref_level_add_bit("PROB", 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, "DISP");
62
63 // set app. order
64 // see Hierarchic Finite Element Bases on Unstructured Tetrahedral Meshes
65 // (Mark Ainsworth & Joe Coyle)
66 int order_temp = 2;
67 CHKERR m_field.set_field_order(root_set, MBTET, "TEMP", order_temp);
68 CHKERR m_field.set_field_order(root_set, MBTRI, "TEMP", order_temp);
69 CHKERR m_field.set_field_order(root_set, MBEDGE, "TEMP", order_temp);
70 CHKERR m_field.set_field_order(root_set, MBVERTEX, "TEMP", 1);
71
72 int order_disp = 3;
73 CHKERR m_field.set_field_order(root_set, MBTET, "DISP", order_disp);
74 CHKERR m_field.set_field_order(root_set, MBTRI, "DISP", order_disp);
75 CHKERR m_field.set_field_order(root_set, MBEDGE, "DISP", order_disp);
76 CHKERR m_field.set_field_order(root_set, MBVERTEX, "DISP", 1);
77
78 ThermalStressElement thermal_stress_elem(m_field);
79 CHKERR thermal_stress_elem.addThermalStressElement("ELAS", "DISP", "TEMP");
80 CHKERR m_field.modify_problem_add_finite_element("PROB", "ELAS");
81
82 /****/
83 // build database
84 // build field
85 CHKERR m_field.build_fields();
86 // build finite elemnts
88 // build adjacencies
89 CHKERR m_field.build_adjacencies(bit_level0);
90
91 // build problem
92 ProblemsManager *prb_mng_ptr;
93 CHKERR m_field.getInterface(prb_mng_ptr);
94 CHKERR prb_mng_ptr->buildProblem("PROB", true);
95 // mesh partitioning
96 // partition
97 CHKERR prb_mng_ptr->partitionSimpleProblem("PROB");
98 CHKERR prb_mng_ptr->partitionFiniteElements("PROB");
99 // what are ghost nodes, see Petsc Manual
100 CHKERR prb_mng_ptr->partitionGhostDofs("PROB");
101
102 // set temerature at nodes
104 MBVERTEX, dof)) {
105 EntityHandle ent = dof->get()->getEnt();
106 VectorDouble coords(3);
107 CHKERR moab.get_coords(&ent, 1, &coords[0]);
108 dof->get()->getFieldData() = 1;
109 }
110
111 Vec F;
112 CHKERR m_field.getInterface<VecManager>()->vecCreateGhost("PROB", ROW, &F);
113 CHKERR thermal_stress_elem.setThermalStressRhsOperators("DISP", "TEMP", F,
114 1);
115
117 "PROB", "ELAS", thermal_stress_elem.getLoopThermalStressRhs());
118 CHKERR VecAssemblyBegin(F);
119 CHKERR VecAssemblyEnd(F);
120
121 // PetscViewer viewer;
122 // CHKERR
123 // PetscViewerASCIIOpen(PETSC_COMM_WORLD,"forces_and_sources_thermal_stress_elem.txt",&viewer);
124 // CHKERR VecChop(F,1e-4);
125 // CHKERR VecView(F,viewer);
126 // CHKERR PetscViewerDestroy(&viewer);
127
128 double sum = 0;
129 CHKERR VecSum(F, &sum);
130 CHKERR PetscPrintf(PETSC_COMM_WORLD, "sum = %9.8f\n", sum);
131 double fnorm;
132 CHKERR VecNorm(F, NORM_2, &fnorm);
133 CHKERR PetscPrintf(PETSC_COMM_WORLD, "fnorm = %9.8e\n", fnorm);
134 if (fabs(sum) > 1e-7) {
135 SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID, "Failed to pass test");
136 }
137 if (fabs(fnorm - 2.64638118e+00) > 1e-7) {
138 SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID, "Failed to pass test");
139 }
140
141 CHKERR VecZeroEntries(F);
142 CHKERR VecDestroy(&F);
143 }
145
147
148 return 0;
149}
int main()
@ 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.
@ F
#define _IT_GET_DOFS_FIELD_BY_NAME_AND_TYPE_FOR_LOOP_(MFIELD, NAME, TYPE, IT)
loop over all dofs from a moFEM field and particular field
virtual MoFEMErrorCode build_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elements.
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.
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.
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
char mesh_file_name[255]
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition Types.hpp:40
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
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.
Problem manager is used to build and partition problems.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Vector manager is used to create vectors \mofem_vectors.
Implentation of thermal stress element.
MyVolumeFE & getLoopThermalStressRhs()
MoFEMErrorCode addThermalStressElement(const std::string fe_name, const std::string field_name, const std::string thermal_field_name, const std::string mesh_nodals_positions="MESH_NODE_POSITIONS")
MoFEMErrorCode setThermalStressRhsOperators(string field_name, string thermal_field_name, Vec &F, int verb=0)
static char help[]
int order_disp
int order_temp