v0.13.2
Loading...
Searching...
No Matches
Functions | Variables
convective_matrix.cpp File Reference
#include <BasicFiniteElements.hpp>

Go to the source code of this file.

Functions

int main (int argc, char *argv[])
 

Variables

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

Detailed Description

Atom test for convective mass element

Definition in file convective_matrix.cpp.

Function Documentation

◆ main()

int main ( int  argc,
char *  argv[] 
)

Definition at line 17 of file convective_matrix.cpp.

17 {
18
19 MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
20
21 try {
22
23 moab::Core mb_instance;
24 moab::Interface &moab = mb_instance;
25 int rank;
26 MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
27
28 PetscBool flg = PETSC_TRUE;
29 char mesh_file_name[255];
30 CHKERR PetscOptionsGetString(PETSC_NULL, PETSC_NULL, "-my_file",
31 mesh_file_name, 255, &flg);
32 if (flg != PETSC_TRUE) {
33 SETERRQ(PETSC_COMM_SELF, 1, "*** ERROR -my_file (MESH FILE NEEDED)");
34 }
35
36 const char *option;
37 option = ""; //"PARALLEL=BCAST;";//;DEBUG_IO";
38 CHKERR moab.load_file(mesh_file_name, 0, option);
39
40 MoFEM::Core core(moab);
41 MoFEM::Interface &m_field = core;
42
43 // ref meshset ref level 0
44 BitRefLevel bit_level0;
45 bit_level0.set(0);
46 EntityHandle meshset_level0;
47 CHKERR moab.create_meshset(MESHSET_SET, meshset_level0);
48 CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(
49 0, 3, bit_level0);
50 CHKERR m_field.getInterface<BitRefManager>()->getEntitiesByRefLevel(
51 bit_level0, BitRefLevel().set(), meshset_level0);
52
53 // Fields
54 CHKERR m_field.add_field("SPATIAL_POSITION", H1, AINSWORTH_LEGENDRE_BASE,
55 3);
56
57 // FE
58 CHKERR m_field.add_finite_element("ELASTIC");
59
60 // Define rows/cols and element data
62 "SPATIAL_POSITION");
64 "SPATIAL_POSITION");
66 "SPATIAL_POSITION");
67
68 // define problems
69 CHKERR m_field.add_problem("ELASTIC_MECHANICS");
70
71 // set finite elements for problems
72 CHKERR m_field.modify_problem_add_finite_element("ELASTIC_MECHANICS",
73 "ELASTIC");
74
75 // set refinement level for problem
76 CHKERR m_field.modify_problem_ref_level_add_bit("ELASTIC_MECHANICS",
77 bit_level0);
78
79 // add entitities (by tets) to the field
80 CHKERR m_field.add_ents_to_field_by_type(0, MBTET, "SPATIAL_POSITION");
81
82 // add finite elements entities
84 bit_level0, BitRefLevel().set(), "ELASTIC", MBTET);
85
86 // set app. order
87 PetscInt order;
88 CHKERR PetscOptionsGetInt(PETSC_NULL, PETSC_NULL, "-my_order", &order,
89 &flg);
90 if (flg != PETSC_TRUE) {
91 order = 1;
92 }
93 CHKERR m_field.set_field_order(0, MBTET, "SPATIAL_POSITION", order);
94 CHKERR m_field.set_field_order(0, MBTRI, "SPATIAL_POSITION", order);
95 CHKERR m_field.set_field_order(0, MBEDGE, "SPATIAL_POSITION", order);
96 CHKERR m_field.set_field_order(0, MBVERTEX, "SPATIAL_POSITION", 1);
97
98 CHKERR m_field.add_finite_element("NEUAMNN_FE");
100 "SPATIAL_POSITION");
101 CHKERR m_field.modify_finite_element_add_field_col("NEUAMNN_FE",
102 "SPATIAL_POSITION");
103 CHKERR m_field.modify_finite_element_add_field_data("NEUAMNN_FE",
104 "SPATIAL_POSITION");
105 CHKERR m_field.modify_problem_add_finite_element("ELASTIC_MECHANICS",
106 "NEUAMNN_FE");
108 it)) {
109 Range tris;
110 CHKERR moab.get_entities_by_type(it->meshset, MBTRI, tris, true);
111 CHKERR m_field.add_ents_to_finite_element_by_type(tris, MBTRI,
112 "NEUAMNN_FE");
113 }
115 m_field, SIDESET | PRESSURESET, it)) {
116 Range tris;
117 CHKERR moab.get_entities_by_type(it->meshset, MBTRI, tris, true);
118 CHKERR m_field.add_ents_to_finite_element_by_type(tris, MBTRI,
119 "NEUAMNN_FE");
120 }
121
122 // Velocity
123 CHKERR m_field.add_field("SPATIAL_VELOCITY", H1, AINSWORTH_LEGENDRE_BASE,
124 3);
125 CHKERR m_field.add_ents_to_field_by_type(0, MBTET, "SPATIAL_VELOCITY");
126 int order_velocity = 1;
127 CHKERR m_field.set_field_order(0, MBTET, "SPATIAL_VELOCITY",
128 order_velocity);
129 CHKERR m_field.set_field_order(0, MBTRI, "SPATIAL_VELOCITY",
130 order_velocity);
131 CHKERR m_field.set_field_order(0, MBEDGE, "SPATIAL_VELOCITY",
132 order_velocity);
133 CHKERR m_field.set_field_order(0, MBVERTEX, "SPATIAL_VELOCITY", 1);
134
135 CHKERR m_field.add_field("DOT_SPATIAL_POSITION", H1,
137 CHKERR m_field.add_ents_to_field_by_type(0, MBTET, "DOT_SPATIAL_POSITION");
138 CHKERR m_field.set_field_order(0, MBTET, "DOT_SPATIAL_POSITION", order);
139 CHKERR m_field.set_field_order(0, MBTRI, "DOT_SPATIAL_POSITION", order);
140 CHKERR m_field.set_field_order(0, MBEDGE, "DOT_SPATIAL_POSITION", order);
141 CHKERR m_field.set_field_order(0, MBVERTEX, "DOT_SPATIAL_POSITION", 1);
142 CHKERR m_field.add_field("DOT_SPATIAL_VELOCITY", H1,
144 CHKERR m_field.add_ents_to_field_by_type(0, MBTET, "DOT_SPATIAL_VELOCITY");
145 CHKERR m_field.set_field_order(0, MBTET, "DOT_SPATIAL_VELOCITY",
146 order_velocity);
147 CHKERR m_field.set_field_order(0, MBTRI, "DOT_SPATIAL_VELOCITY",
148 order_velocity);
149 CHKERR m_field.set_field_order(0, MBEDGE, "DOT_SPATIAL_VELOCITY",
150 order_velocity);
151 CHKERR m_field.set_field_order(0, MBVERTEX, "DOT_SPATIAL_VELOCITY", 1);
152
153 ConvectiveMassElement inertia(m_field, 1);
154 CHKERR inertia.setBlocks();
155 CHKERR inertia.addConvectiveMassElement("MASS_ELEMENT", "SPATIAL_VELOCITY",
156 "SPATIAL_POSITION");
157 CHKERR inertia.addVelocityElement("VELOCITY_ELEMENT", "SPATIAL_VELOCITY",
158 "SPATIAL_POSITION");
159 CHKERR m_field.modify_problem_add_finite_element("ELASTIC_MECHANICS",
160 "MASS_ELEMENT");
161 CHKERR m_field.modify_problem_add_finite_element("ELASTIC_MECHANICS",
162 "VELOCITY_ELEMENT");
163
164 // build field
165 CHKERR m_field.build_fields();
166
167 double scale_positions = 2;
168 {
169 EntityHandle node = 0;
170 double coords[3];
171 for (_IT_GET_DOFS_FIELD_BY_NAME_FOR_LOOP_(m_field, "SPATIAL_POSITION",
172 dof_ptr)) {
173 if (dof_ptr->get()->getEntType() != MBVERTEX)
174 continue;
175 EntityHandle ent = dof_ptr->get()->getEnt();
176 int dof_rank = dof_ptr->get()->getDofCoeffIdx();
177 double &fval = dof_ptr->get()->getFieldData();
178 if (node != ent) {
179 CHKERR moab.get_coords(&ent, 1, coords);
180 node = ent;
181 }
182 fval = scale_positions * coords[dof_rank];
183 }
184 }
185
186 double scale_velocities = 4;
187 {
188 EntityHandle node = 0;
189 double coords[3];
190 for (_IT_GET_DOFS_FIELD_BY_NAME_FOR_LOOP_(m_field, "DOT_SPATIAL_POSITION",
191 dof_ptr)) {
192 if (dof_ptr->get()->getEntType() != MBVERTEX)
193 continue;
194 EntityHandle ent = dof_ptr->get()->getEnt();
195 int dof_rank = dof_ptr->get()->getDofCoeffIdx();
196 double &fval = dof_ptr->get()->getFieldData();
197 if (node != ent) {
198 CHKERR moab.get_coords(&ent, 1, coords);
199 node = ent;
200 }
201 fval = scale_velocities * coords[dof_rank];
202 }
203 }
204
205 // build finite elemnts
207
208 // build adjacencies
209 CHKERR m_field.build_adjacencies(bit_level0);
210
211 // build problem
212
213 ProblemsManager *prb_mng_ptr;
214 CHKERR m_field.getInterface(prb_mng_ptr);
215
216 CHKERR prb_mng_ptr->buildProblem("ELASTIC_MECHANICS", true);
217
218 // partition
219 CHKERR prb_mng_ptr->partitionProblem("ELASTIC_MECHANICS");
220 CHKERR prb_mng_ptr->partitionFiniteElements("ELASTIC_MECHANICS");
221 CHKERR prb_mng_ptr->partitionGhostDofs("ELASTIC_MECHANICS");
222
223 // create matrices
224 Vec F;
225 CHKERR m_field.getInterface<VecManager>()->vecCreateGhost(
226 "ELASTIC_MECHANICS", COL, &F);
227 Vec D;
228 CHKERR VecDuplicate(F, &D);
229 Mat Aij;
231 ->createMPIAIJWithArrays<PetscGlobalIdx_mi_tag>("ELASTIC_MECHANICS",
232 &Aij);
233
234 CHKERR inertia.setConvectiveMassOperators("SPATIAL_VELOCITY",
235 "SPATIAL_POSITION");
236 CHKERR inertia.setVelocityOperators("SPATIAL_VELOCITY", "SPATIAL_POSITION");
237
238 inertia.getLoopFeMassRhs().ts_F = F;
239 inertia.getLoopFeMassRhs().ts_a = 1;
240 inertia.getLoopFeMassLhs().ts_B = Aij;
241 inertia.getLoopFeMassLhs().ts_a = 1;
242
243 inertia.getLoopFeVelRhs().ts_F = F;
244 inertia.getLoopFeVelRhs().ts_a = 1;
245 inertia.getLoopFeVelLhs().ts_B = Aij;
246 inertia.getLoopFeVelLhs().ts_a = 1;
247
248 CHKERR m_field.loop_finite_elements("ELASTIC_MECHANICS", "MASS_ELEMENT",
249 inertia.getLoopFeMassRhs());
250 CHKERR m_field.loop_finite_elements("ELASTIC_MECHANICS", "VELOCITY_ELEMENT",
251 inertia.getLoopFeVelRhs());
252 CHKERR VecGhostUpdateBegin(F, ADD_VALUES, SCATTER_REVERSE);
253 CHKERR VecGhostUpdateEnd(F, ADD_VALUES, SCATTER_REVERSE);
254 CHKERR VecAssemblyBegin(F);
255 CHKERR VecAssemblyEnd(F);
256
257 CHKERR m_field.loop_finite_elements("ELASTIC_MECHANICS", "MASS_ELEMENT",
258 inertia.getLoopFeMassLhs());
259 CHKERR m_field.loop_finite_elements("ELASTIC_MECHANICS", "VELOCITY_ELEMENT",
260 inertia.getLoopFeVelLhs());
261 CHKERR MatAssemblyBegin(Aij, MAT_FINAL_ASSEMBLY);
262 CHKERR MatAssemblyEnd(Aij, MAT_FINAL_ASSEMBLY);
263
264 double sum = 0;
265 CHKERR VecSum(F, &sum);
266 CHKERR PetscPrintf(PETSC_COMM_WORLD, "sum = %9.8e\n", sum);
267 double fnorm;
268 CHKERR VecNorm(F, NORM_2, &fnorm);
269 CHKERR PetscPrintf(PETSC_COMM_WORLD, "fnorm = %9.8e\n", fnorm);
270
271 double mnorm;
272 CHKERR MatNorm(Aij, NORM_1, &mnorm);
273 CHKERR PetscPrintf(PETSC_COMM_WORLD, "mnorm = %9.8e\n", mnorm);
274
275 if (fabs(sum - 6.27285463e+00) > 1e-8) {
276 SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID, "Failed to pass test");
277 }
278 if (fabs(fnorm - 1.28223353e+00) > 1e-6) {
279 SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID, "Failed to pass test");
280 }
281 if (fabs(mnorm - 1.31250000e+00) > 1e-6) {
282 SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID, "Failed to pass test");
283 }
284
285 CHKERR VecDestroy(&F);
286 CHKERR VecDestroy(&D);
287 CHKERR MatDestroy(&Aij);
288 }
290
292
293 return 0;
294}
static char help[]
@ COL
Definition: definitions.h:123
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:372
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
@ H1
continuous field
Definition: definitions.h:85
@ PRESSURESET
Definition: definitions.h:152
@ FORCESET
Definition: definitions.h:151
@ NODESET
Definition: definitions.h:146
@ SIDESET
Definition: definitions.h:147
@ MOFEM_ATOM_TEST_INVALID
Definition: definitions.h:40
#define CHKERR
Inline error check.
Definition: definitions.h:535
#define _IT_GET_DOFS_FIELD_BY_NAME_FOR_LOOP_(MFIELD, NAME, IT)
Definition: Interface.hpp:1843
virtual MoFEMErrorCode add_finite_element(const std::string &fe_name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
add finite element
virtual MoFEMErrorCode add_ents_to_finite_element_by_bit_ref(const BitRefLevel &bit, const BitRefLevel &mask, const std::string &name, EntityType type, int verb=DEFAULT_VERBOSITY)=0
add TET entities from given refinement level to finite element database given by name
virtual MoFEMErrorCode build_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elements.
virtual MoFEMErrorCode modify_finite_element_add_field_col(const std::string &fe_name, const std::string name_row)=0
set field col which finite element use
virtual MoFEMErrorCode add_ents_to_finite_element_by_type(const EntityHandle entities, const EntityType type, const std::string &name, const bool recursive=true)=0
add entities to finite element
virtual MoFEMErrorCode modify_finite_element_add_field_data(const std::string &fe_name, const std::string name_filed)=0
set finite element field data
virtual MoFEMErrorCode modify_finite_element_add_field_row(const std::string &fe_name, const std::string name_row)=0
set field row which finite element use
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.
#define _IT_CUBITMESHSETS_BY_BCDATA_TYPE_FOR_LOOP_(MESHSET_MANAGER, CUBITBCTYPE, IT)
Iterator that loops over a specific Cubit MeshSet in a moFEM field.
MoFEMErrorCode partitionGhostDofs(const std::string name, int verb=VERBOSE)
determine ghost nodes
MoFEMErrorCode buildProblem(const std::string name, const bool square_matrix, int verb=VERBOSE)
build problem data structures
MoFEMErrorCode partitionProblem(const std::string name, int verb=VERBOSE)
partition problem dofs (collective)
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
double D
char mesh_file_name[255]
const FTensor::Tensor2< T, Dim, Dim > Vec
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition: Types.hpp:40
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
structure grouping operators and data used for calculation of mass (convective) element \ nonlinear_e...
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.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.
Vector manager is used to create vectors \mofem_vectors.
Definition: VecManager.hpp:23

Variable Documentation

◆ help

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

Definition at line 15 of file convective_matrix.cpp.