PeriDEM 0.2.0
PeriDEM -- Peridynamics-based high-fidelity model for granular media
Loading...
Searching...
No Matches
fe::Mesh Class Reference

A class for mesh data. More...

#include <mesh.h>

Collaboration diagram for fe::Mesh:

Public Member Functions

 Mesh (size_t dim=0)
 Constructor.
 
 Mesh (inp::MeshDeck *deck)
 Constructor.
 
std::string printStr (int nt=0, int lvl=0) const
 Returns the string containing printable information about the object.
 
void print (int nt=0, int lvl=0) const
 Prints the information about the object.
 
Accessor methods
size_t getDimension () const
 Get the dimension of the domain.
 
size_t getNumNodes () const
 Get the number of nodes.
 
size_t getNumElements () const
 Get the number of elements.
 
size_t getNumDofs () const
 Get the number of dofs.
 
size_t getElementType () const
 Get the type of element in mesh.
 
double getMeshSize () const
 Get the mesh size.
 
util::Point getNode (const size_t &i) const
 Get coordinates of node i.
 
double getNodalVolume (const size_t &i) const
 Get nodal volume of node i.
 
const std::vector< util::Point > & getNodes () const
 Get the nodes data.
 
std::vector< util::Point > & getNodes ()
 Get the nodes data.
 
const std::vector< util::Point > * getNodesP () const
 Get the pointer to nodes data.
 
std::vector< util::Point > * getNodesP ()
 Get the pointer to nodes data.
 
const std::vector< uint8_t > * getFixityP () const
 Get the pointer to fixity data.
 
std::vector< uint8_t > * getFixityP ()
 Get the pointer to fixity data.
 
const std::vector< uint8_t > & getFixity () const
 Get the reference to fixity data.
 
std::vector< uint8_t > & getFixity ()
 Get the reference to fixity data.
 
const std::vector< double > & getNodalVolumes () const
 Get the nodal volume data.
 
std::vector< double > & getNodalVolumes ()
 Get the nodal volume data.
 
const std::vector< double > * getNodalVolumesP () const
 Get the pointer to nodal volume data.
 
std::vector< double > * getNodalVolumesP ()
 Get the pointer to nodal volume data.
 
bool isNodeFree (const size_t &i, const unsigned int &dof) const
 Return true if node is free.
 
std::vector< size_tgetElementConnectivity (const size_t &i) const
 Get the connectivity of element.
 
std::vector< util::PointgetElementConnectivityNodes (const size_t &i) const
 Get the vertices of element.
 
const std::vector< size_t > & getElementConnectivities () const
 Get the reference to element-node connectivity data.
 
std::vector< size_t > & getElementConnectivities ()
 Get the reference to element-node connectivity data.
 
const std::vector< size_t > * getElementConnectivitiesP () const
 Get the pointer to element-node connectivity data.
 
std::vector< size_t > * getElementConnectivitiesP ()
 Get the pointer to element-node connectivity data.
 
const std::pair< std::vector< double >, std::vector< double > > & getBoundingBox () const
 Get the bounding box of the mesh.
 
std::pair< std::vector< double >, std::vector< double > > & getBoundingBox ()
 Get the bounding box of the mesh.
 
Setter methods
void setFixity (const size_t &i, const unsigned int &dof, const bool &flag)
 Set the fixity to free (0) or fixed (1)
 
void clearElementData ()
 Clear element-node connectivity data.
 
Utility methods
void createData (const std::string &filename, bool ref_config=false)
 Reads mesh data from the file and populates other data.
 
bool readElementData (const std::string &filename)
 Reads element-node connectivity data from file. This function is meant for cases when mesh was created without element-node connectivity data but later during output, strain/stress were required which needs element-node connectivity data.
 
void computeVol ()
 Compute the nodal volume.
 
void computeBBox ()
 Compute the bounding box

 
void computeMeshSize ()
 Compute the mesh size.
 

Data Fields

size_t d_dim
 Dimension of the mesh.
 
std::string d_spatialDiscretization
 Tag for spatial discretization type.
 
std::string d_filename
 Filename to read mesh data.
 
bool d_needEncData
 Flag that indicates whether we need enc data (set by input mesh deck in constructor)
 
bool d_encDataPopulated
 Flag that indicates whether element-node connectivity data is read from file.
 
size_t d_numDofs
 Number of dofs = (dimension) times (number of nodes)
 
std::vector< size_td_gMap
 Map from global reduced id to default global id.
 
std::vector< intd_gInvMap
 Map from global id to reduced global id.
 
std::pair< std::vector< double >, std::vector< double > > d_bbox
 Bounding box.
 
double d_h
 Mesh size.
 
Mesh data
size_t d_numNodes
 Number of nodes.
 
size_t d_numElems
 Number of elements.
 
size_t d_eType
 Element type.
 
size_t d_eNumVertex
 Number of vertex per element.
 
std::vector< util::Pointd_nodes
 Vector of initial (reference) coordinates of nodes.
 
std::vector< size_td_enc
 Element-node connectivity data.
 
std::vector< std::vector< size_t > > d_nec
 Node-element connectivity data.
 
std::vector< uint8_td_fix
 Vector of fixity mask of each node.
 
std::vector< doubled_vol
 Vector of volume of each node.
 
Mesh data specific to parallel implementation
size_t d_nPart
 Number of partitions.
 
std::string d_partitionMethod
 Partitioning method. It could be either empty string or "metis_recursive" or "metis_kway".
 
std::vector< size_td_nodePartition
 Node partition information. For each node i, d_nodePartition[i] specifies the partition number, i.e., the processor that owns the node in MPI application.
 

Detailed Description

A class for mesh data.

In this class the mesh data such as nodes, element-node connectivity, node-element connectivity are stored. The class also stores fixity mask of nodes which indicates if x-, y-, or z-dof of the node is fixed or free.

We currently only support mesh with only one type of elements, i.e. mesh can not have mix of two types of elements. For example, we can not have mesh with triangle and quadrangle elements together.

This class is used in both finite difference implementation and finite element implementation. For finite difference, we only require nodal volume. If the mesh file contains nodal volume, we skip reading element-node and node-element connectivity, however if mesh file does not have nodal volume data, we read connectivity data and compute the nodal volume.

Definition at line 51 of file mesh.h.

Constructor & Destructor Documentation

◆ Mesh() [1/2]

fe::Mesh::Mesh ( size_t  dim = 0)
explicit

Constructor.

Parameters
dimDimension of the domain

Definition at line 26 of file mesh.cpp.

28 d_h(0.), d_dim(dim), d_encDataPopulated(false), d_needEncData(false),
29 d_nPart(0){}
double d_h
Mesh size.
Definition mesh.h:516
size_t d_nPart
Number of partitions.
Definition mesh.h:449
size_t d_dim
Dimension of the mesh.
Definition mesh.h:468
size_t d_numElems
Number of elements.
Definition mesh.h:378
size_t d_eType
Element type.
Definition mesh.h:389
bool d_encDataPopulated
Flag that indicates whether element-node connectivity data is read from file.
Definition mesh.h:487
size_t d_numNodes
Number of nodes.
Definition mesh.h:375
size_t d_eNumVertex
Number of vertex per element.
Definition mesh.h:404
bool d_needEncData
Flag that indicates whether we need enc data (set by input mesh deck in constructor)
Definition mesh.h:484
size_t d_numDofs
Number of dofs = (dimension) times (number of nodes)
Definition mesh.h:490

◆ Mesh() [2/2]

fe::Mesh::Mesh ( inp::MeshDeck deck)
explicit

Constructor.

The constructor initializes the data using input deck, performs checks on input data, and reads mesh file and populates the mesh related data. The mesh file of **.csv**, **.vtu (VTK)** and **.msh (Gmsh)** are supported.

Parameters
deckInput deck which contains user-specified information

Definition at line 31 of file mesh.cpp.

33 d_h(deck->d_h), d_dim(deck->d_dim),
34 d_spatialDiscretization(deck->d_spatialDiscretization),
35 d_filename(deck->d_filename), d_encDataPopulated(false),
36 d_needEncData(deck->d_populateElementNodeConnectivity),
37 d_nPart(0) {
38
39 // perform check on input data
40 if (d_spatialDiscretization != "finite_difference" and
41 d_spatialDiscretization != "weak_finite_element" and
42 d_spatialDiscretization != "nodal_finite_element" and
43 d_spatialDiscretization != "truss_finite_element") {
44 std::cerr << "Error: Spatial discretization type not known. Check input "
45 "data.\n";
46 std::cerr << deck->printStr() << "\n";
47 exit(1);
48 }
49
51 std::cerr << "Error: Check Dimension in input data.\n";
52 exit(1);
53 }
54
55 if (d_filename.empty()) {
56 std::cerr << "Error: Filename for mesh data not specified.\n";
57 exit(1);
58 }
59
60 // read mesh data from file
62}
std::string d_spatialDiscretization
Tag for spatial discretization type.
Definition mesh.h:478
std::string d_filename
Filename to read mesh data.
Definition mesh.h:481
void createData(const std::string &filename, bool ref_config=false)
Reads mesh data from the file and populates other data.
Definition mesh.cpp:67

References createData(), d_dim, d_filename, d_spatialDiscretization, and inp::MeshDeck::printStr().

Here is the call graph for this function:

Member Function Documentation

◆ clearElementData()

void fe::Mesh::clearElementData ( )

Clear element-node connectivity data.

Definition at line 396 of file mesh.cpp.

396 {
397 if (!d_enc.empty())
398 d_enc.shrink_to_fit();
399 d_numElems = 0;
400 if (!d_nec.empty())
401 d_nec.shrink_to_fit();
402}
std::vector< std::vector< size_t > > d_nec
Node-element connectivity data.
Definition mesh.h:421
std::vector< size_t > d_enc
Element-node connectivity data.
Definition mesh.h:415

◆ computeBBox()

void fe::Mesh::computeBBox ( )

Compute the bounding box

Definition at line 332 of file mesh.cpp.

332 {
333 std::vector<double> p1(3,0.);
334 std::vector<double> p2(3,0.);
335 for (const auto& x : d_nodes) {
336 if (util::isLess(x.d_x, p1[0]))
337 p1[0] = x.d_x;
338 if (util::isLess(x.d_y, p1[1]))
339 p1[1] = x.d_y;
340 if (util::isLess(x.d_z, p1[2]))
341 p1[2] = x.d_z;
342 if (util::isLess(p2[0], x.d_x))
343 p2[0] = x.d_x;
344 if (util::isLess(p2[1], x.d_y))
345 p2[1] = x.d_y;
346 if (util::isLess(p2[2], x.d_z))
347 p2[2] = x.d_z;
348 }
349
350 d_bbox = std::make_pair(p1, p2);
351}
std::vector< util::Point > d_nodes
Vector of initial (reference) coordinates of nodes.
Definition mesh.h:407
std::pair< std::vector< double >, std::vector< double > > d_bbox
Bounding box.
Definition mesh.h:513
bool isLess(const double &a, const double &b)
Returns true if a < b.
Definition function.cpp:20

References util::isLess().

Here is the call graph for this function:

◆ computeMeshSize()

void fe::Mesh::computeMeshSize ( )

Compute the mesh size.

This method searches for minimum distance between any two mesh nodes and stores it as a mesh size.

Definition at line 353 of file mesh.cpp.

353 {
354
355 double guess = 0.;
356 if (d_nodes.size() < 2) {
357 d_h = 0.;
358 return;
359 }
360
361 guess = (d_nodes[0] - d_nodes[1]).length();
362 for (size_t i = 0; i < d_nodes.size(); i++)
363 for (size_t j = 0; j < d_nodes.size(); j++)
364 if (i != j) {
365 double val = d_nodes[i].dist(d_nodes[j]);
366
367 if (util::isLess(val, 1.0E-12)) {
368
369 std::cout << "Check nodes are too close = "
370 << util::io::printStr<util::Point>({d_nodes[i],
371 d_nodes[j]})
372 << "\n";
373 std::cout << "Distance = " << val << ", guess = " << guess << "\n";
374 }
375 if (util::isLess(val, guess))
376 guess = val;
377 }
378
379 d_h = guess;
380}

References util::isLess().

Here is the call graph for this function:

◆ computeVol()

void fe::Mesh::computeVol ( )

Compute the nodal volume.

This method requires element-node connectivity data to compute the nodal volumes. Formula for volume of a node \( i\) is given by

\[ V_i = \sum_{e \in \mathbf{N}_i} \int_{T_e} N_i(x) dx, \]

where \(\mathbf{N}_i\) is a list of elements which have node \( i\) as its vertex, \( T_e\) is the element domain, \( N_i\) is the shape function of the node \( i\) in element e.

Definition at line 250 of file mesh.cpp.

250 {
251
252 // initialize quadrature data
255 quads = new fe::TriElem(2);
256 else if (d_eType == util::vtk_type_quad)
257 quads = new fe::QuadElem(2);
258 else if (d_eType == util::vtk_type_tetra)
259 quads = new fe::TetElem(2);
260
261 // check if we have valid element-node connectivity data for nodal volume
262 // calculations
263 if (d_nec.size() != d_numNodes || d_enc.empty()) {
264 std::cerr << "Error: Can not compute nodal volume for given finite "
265 "element mesh as the element-node connectivity data is "
266 "invalid."
267 << std::endl;
269 }
270
271 if (false) {
272 print(0, 0);
273 std::cout << "\n-------- Node data ----------\n";
274 std::cout << util::io::printStr(d_nodes, 0) << "\n";
275 std::cout << "\n-------- Element data ----------\n";
276 std::cout << util::io::printStr(d_enc, 0) << "\n";
277 }
278
279 //
280 // compute nodal volume
281 //
282 d_vol.resize(d_numNodes);
283
285 tf::Taskflow taskflow;
286
287 taskflow.for_each_index(
288 (std::size_t) 0, this->d_numNodes, (std::size_t) 1, [this, quads](std::size_t i) {
289 double v = 0.0;
290
291 for (auto e : this->d_nec[i]) {
292
293 std::vector<size_t> e_ns = this->getElementConnectivity(e);
294
295 // locate global node i in local list of element el
296 int loc_i = -1;
297 for (size_t l = 0; l < e_ns.size(); l++)
298 if (e_ns[l] == i)
299 loc_i = l;
300
301 if (loc_i == -1) {
302 std::cerr << "Error: Check node element connectivity.\n";
303 exit(1);
304 }
305
306 // get quad data
307 std::vector<util::Point> e_nodes;
308 for (auto k : e_ns)
310
311 // get volume of element
312 double vol = quads->elemSize(e_nodes);
313 double factor = 1.;
314 if (vol < 0.)
315 factor = -1.;
316
317 std::vector<fe::QuadData> qds = quads->getQuadDatas(e_nodes);
318
319 // compute V_e and add it to volume
320 for (auto qd : qds)
321 v += qd.d_shapes[loc_i] * factor * qd.d_w;
322 } // loop over elements
323
324 // update
325 this->d_vol[i] = v;
326 }
327 ); // for_each
328
329 executor.run(taskflow).get();
330}
A base class which provides methods to map points to/from reference element and to compute quadrature...
Definition baseElem.h:84
std::vector< double > d_vol
Vector of volume of each node.
Definition mesh.h:439
std::vector< size_t > getElementConnectivity(const size_t &i) const
Get the connectivity of element.
Definition mesh.h:210
void print(int nt=0, int lvl=0) const
Prints the information about the object.
Definition mesh.h:306
A class for mapping and quadrature related operations for bi-linear quadrangle element.
Definition quadElem.h:64
A class for mapping and quadrature related operations for linear tetrahedron element.
Definition tetElem.h:141
A class for mapping and quadrature related operations for linear triangle element.
Definition triElem.h:91
static const int vtk_type_triangle
Integer flag for triangle element.
static const int vtk_type_quad
Integer flag for quad element.
static const int vtk_type_tetra
Integer flag for tetrahedron element.
std::string printStr(const T &msg, int nt=print_default_tab)
Returns formatted string for output.
Definition io.h:54
unsigned int getNThreads()
Get number of threads to be used by taskflow.

References fe::BaseElem::elemSize(), util::parallel::getNThreads(), fe::BaseElem::getQuadDatas(), util::io::printStr(), util::vtk_type_quad, util::vtk_type_tetra, and util::vtk_type_triangle.

Here is the call graph for this function:

◆ createData()

void fe::Mesh::createData ( const std::string &  filename,
bool  ref_config = false 
)

Reads mesh data from the file and populates other data.

This function calls reader methods in namespace rw::reader to read the mesh file. For finite difference implementation, we support **.csv** mesh file which has nodal coordinates and nodal volumes data.

However, for finite element implementation, we require either **.vtu** or **.msh** file with element-node connectivity data.

Parameters
filenameName of the mesh file
ref_configFlag which specifies if we need to subtract the displacement from nodes obtained from vtu file to get reference position of nodes

Definition at line 67 of file mesh.cpp.

67 {
68
69 util::io::log("Mesh: Reading element data.\n");
70
71 int file_type = -1;
72 // find the extension of file and call correct reader
74 file_type = 0;
76 file_type = 1;
78 file_type = 2;
79 else {
80 std::cerr << "Error: Currently only '.csv', '.msg', and '.vtu' "
81 "files are supported for reading mesh.\n";
83 }
84
85 if (d_spatialDiscretization != "finite_difference" and file_type == 0) {
86
87 std::cerr << "Error: For discretization = " << d_spatialDiscretization
88 << " .vtu or .msh mesh file is required.\n";
89 exit(1);
90 }
91
92 //
93 bool is_fd = false;
94 if (d_spatialDiscretization == "finite_difference")
95 is_fd = true;
96
97 // read node and elements
98 if (file_type == 0)
100 else if (file_type == 1) {
102 &d_enc, &d_nec, &d_vol, false);
103 d_encDataPopulated = true;
104 }
105 else if (file_type == 2) {
106 //
107 // old reading of mesh
108 //
109 // rw::reader::readVtuFile(filename, d_dim, &d_nodes, d_eType, d_numElems,
110 // &d_enc, &d_nec, &d_vol, false);
111
112 //
113 // new reading of mesh
114 // We read the data from file one by one depending on what data we need
115 //
116
117 // read node
119
120 // read volume if required
121 bool found_volume_data = false;
122 if (is_fd) {
125
126 // try another tag for nodal volume
130 }
131
132 // read element data (only if this is fe simulation or if we need
133 // element-node connectivity data for nodal volume calculation)
134 if (!is_fd || !found_volume_data) {
136 &d_nec);
137 d_encDataPopulated = true;
138 }
139
140 // check if file has fixity data
142 }
143
144 // compute data from mesh data
145 d_numNodes = d_nodes.size();
148
149 //
150 // assign default values to fixity
151 //
152 if (d_fix.size() != d_numNodes)
153 d_fix = std::vector<uint8_t>(d_nodes.size(), uint8_t(0));
154
155 //
156 // compute nodal volume if required
157 //
158 bool compute_vol = false;
159 if (is_fd and d_vol.empty()) compute_vol = true;
160
161 // if this is weak finite element simulation then check from policy if
162 // volume is to be computed
163 if (d_spatialDiscretization == "weak_finite_element")
164 compute_vol = false;
165
166 if (compute_vol) {
167 util::io::log("Mesh: Computing nodal volume.\n");
168 computeVol();
169 }
170
171 //
172 // compute bounding box
173 //
174 computeBBox();
175
176 // check if we need to compute mesh size
177 //if (deck->d_computeMeshSize)
179
180 // check nodal volume
181 size_t counter = 0;
182 for (const auto &v : d_vol) {
183
184 if (v < 0.01 * std::pow(d_h, d_dim)) {
185
186 std::cerr << "Error: Check nodal volume " << v
187 << " is less than " << 0.01 * std::pow(d_h, d_dim)
188 << ", Node = " << counter
189 << " at position = " << d_nodes[counter].printStr() << "\n"
190 << "mesh filename = " << filename << "\n"
191 << printStr() << "\n";
192
193 exit(1);
194 }
195
196 counter++;
197 }
198
201}
bool readElementData(const std::string &filename)
Reads element-node connectivity data from file. This function is meant for cases when mesh was create...
Definition mesh.cpp:203
std::vector< uint8_t > d_fix
Vector of fixity mask of each node.
Definition mesh.h:431
void computeBBox()
Compute the bounding box
Definition mesh.cpp:332
void computeMeshSize()
Compute the mesh size.
Definition mesh.cpp:353
void computeVol()
Compute the nodal volume.
Definition mesh.cpp:250
std::string printStr(int nt=0, int lvl=0) const
Returns the string containing printable information about the object.
Definition mesh.cpp:404
static int vtk_map_element_to_num_nodes[16]
Map from element type to number of nodes (for vtk)
void readVtuFileCells(const std::string &filename, size_t dim, size_t &element_type, size_t &num_elem, std::vector< size_t > *enc, std::vector< std::vector< size_t > > *nec)
Reads cell data, i.e. element-node connectivity and node-element connectivity.
Definition reader.cpp:168
bool readVtuFilePointData(const std::string &filename, const std::string &tag, std::vector< uint8_t > *data)
Reads data of specified tag from the vtu file.
Definition reader.cpp:208
void readVtuFileNodes(const std::string &filename, size_t dim, std::vector< util::Point > *nodes, bool ref_config=false)
Reads nodal coordinates.
Definition reader.cpp:137
void readCsvFile(const std::string &filename, size_t dim, std::vector< util::Point > *nodes, std::vector< double > *volumes)
Reads mesh data into node file and element file.
Definition reader.cpp:16
void readMshFile(const std::string &filename, size_t dim, std::vector< util::Point > *nodes, size_t &element_type, size_t &num_elem, std::vector< size_t > *enc, std::vector< std::vector< size_t > > *nec, std::vector< double > *volumes, bool is_fd=false)
Reads mesh data into node file and element file.
Definition reader.cpp:351
std::string getExtensionFromFile(std::string const &filename)
Get extension from the filename.
Definition io.h:428
void log(std::ostringstream &oss, bool screen_out=false, int printMpiRank=print_default_mpi_rank)
Global method to log the message.
Definition io.cpp:38

References util::io::getExtensionFromFile(), util::io::log(), rw::reader::readCsvFile(), rw::reader::readMshFile(), rw::reader::readVtuFileCells(), rw::reader::readVtuFileNodes(), rw::reader::readVtuFilePointData(), and util::vtk_map_element_to_num_nodes.

Referenced by Mesh().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ getBoundingBox() [1/2]

std::pair< std::vector< double >, std::vector< double > > & fe::Mesh::getBoundingBox ( )
inline

Get the bounding box of the mesh.

Returns
box Bounding box

Definition at line 265 of file mesh.h.

265 {
266 return d_bbox;
267 };

References d_bbox.

◆ getBoundingBox() [2/2]

const std::pair< std::vector< double >, std::vector< double > > & fe::Mesh::getBoundingBox ( ) const
inline

Get the bounding box of the mesh.

Returns
box Bounding box

Definition at line 259 of file mesh.h.

260 {
261 return d_bbox;
262 };

References d_bbox.

◆ getDimension()

size_t fe::Mesh::getDimension ( ) const
inline

Get the dimension of the domain.

Returns
N Dimension

Definition at line 82 of file mesh.h.

82{ return d_dim; };

References d_dim.

Referenced by fe::getMaxShearStressAndLoc(), and fe::getStrainStress().

Here is the caller graph for this function:

◆ getElementConnectivities() [1/2]

std::vector< size_t > & fe::Mesh::getElementConnectivities ( )
inline

Get the reference to element-node connectivity data.

Returns
reference Reference

Definition at line 238 of file mesh.h.

238 {
239 return d_enc;
240 };

References d_enc.

◆ getElementConnectivities() [2/2]

const std::vector< size_t > & fe::Mesh::getElementConnectivities ( ) const
inline

Get the reference to element-node connectivity data.

Returns
reference Reference

Definition at line 233 of file mesh.h.

233 {
234 return d_enc;
235 };

References d_enc.

◆ getElementConnectivitiesP() [1/2]

std::vector< size_t > * fe::Mesh::getElementConnectivitiesP ( )
inline

Get the pointer to element-node connectivity data.

Returns
pointer Pointer

Definition at line 251 of file mesh.h.

251 {
252 return &d_enc;
253 };

References d_enc.

◆ getElementConnectivitiesP() [2/2]

const std::vector< size_t > * fe::Mesh::getElementConnectivitiesP ( ) const
inline

Get the pointer to element-node connectivity data.

Returns
pointer Pointer

Definition at line 246 of file mesh.h.

246 {
247 return &d_enc;
248 };

References d_enc.

◆ getElementConnectivity()

std::vector< size_t > fe::Mesh::getElementConnectivity ( const size_t i) const
inline

Get the connectivity of element.

Since we store connectivity in a single vector, we use fe::Mesh::d_eNumVertex to get the connectivity of element. Given element e, the connectivity of e begins from location \( i_0 = e*d\_eNumVertex + 0 \) upto \(i_{n-1} = e*d\_eNumVertex + d\_eNumVertex - 1\).

So the connectivity of e is d_enc[ \(i_0\)], d_enc[ \( i_1 \)], ..., d_end[ \(i_{n-1}\)]

Parameters
iId of an element
Returns
vector Vector of nodal ids

Definition at line 210 of file mesh.h.

210 {
211 return std::vector<size_t>(d_enc.begin() + d_eNumVertex * i,
212 d_enc.begin() + d_eNumVertex * i + d_eNumVertex);
213 };

References d_dim, d_enc, and d_eNumVertex.

Referenced by fe::getCurrentQuadPoints(), fe::getMaxShearStressAndLoc(), and fe::getStrainStress().

Here is the caller graph for this function:

◆ getElementConnectivityNodes()

std::vector< util::Point > fe::Mesh::getElementConnectivityNodes ( const size_t i) const
inline

Get the vertices of element.

Parameters
iId of an element
Returns
vector Vector of vertices

Definition at line 221 of file mesh.h.

222 {
223 std::vector<util::Point> nds;
224 for (size_t k = 0; k < d_eNumVertex; k++)
225 nds.emplace_back(d_nodes[d_enc[d_eNumVertex * i + k]]);
226 return nds;
227 };

References d_dim, d_enc, d_eNumVertex, and d_nodes.

◆ getElementType()

size_t fe::Mesh::getElementType ( ) const
inline

Get the type of element in mesh.

Returns
type Element type (using VTK convention)

Definition at line 106 of file mesh.h.

106{ return d_eType; };

References d_eType.

Referenced by fe::getCurrentQuadPoints(), fe::getMaxShearStressAndLoc(), and fe::getStrainStress().

Here is the caller graph for this function:

◆ getFixity() [1/2]

std::vector< uint8_t > & fe::Mesh::getFixity ( )
inline

Get the reference to fixity data.

Returns
reference Reference to fixity data

Definition at line 162 of file mesh.h.

162{ return d_fix; };

References d_fix.

◆ getFixity() [2/2]

const std::vector< uint8_t > & fe::Mesh::getFixity ( ) const
inline

Get the reference to fixity data.

Returns
reference Reference to fixity data

Definition at line 159 of file mesh.h.

159{ return d_fix; };

References d_fix.

◆ getFixityP() [1/2]

std::vector< uint8_t > * fe::Mesh::getFixityP ( )
inline

Get the pointer to fixity data.

Returns
pointer Pointer to fixity data

Definition at line 153 of file mesh.h.

153{ return &d_fix; };

References d_fix.

◆ getFixityP() [2/2]

const std::vector< uint8_t > * fe::Mesh::getFixityP ( ) const
inline

Get the pointer to fixity data.

Returns
pointer Pointer to fixity data

Definition at line 150 of file mesh.h.

150{ return &d_fix; };

References d_fix.

◆ getMeshSize()

double fe::Mesh::getMeshSize ( ) const
inline

Get the mesh size.

Returns
h Mesh size

Definition at line 112 of file mesh.h.

112{ return d_h; };

References d_h.

◆ getNodalVolume()

double fe::Mesh::getNodalVolume ( const size_t i) const
inline

Get nodal volume of node i.

Parameters
iId of the node
Returns
vol Volume

Definition at line 126 of file mesh.h.

126{ return d_vol[i]; };

References d_dim, and d_vol.

◆ getNodalVolumes() [1/2]

std::vector< double > & fe::Mesh::getNodalVolumes ( )
inline

Get the nodal volume data.

Returns
Vector Vector of nodal volume

Definition at line 171 of file mesh.h.

171{ return d_vol; };

References d_vol.

◆ getNodalVolumes() [2/2]

const std::vector< double > & fe::Mesh::getNodalVolumes ( ) const
inline

Get the nodal volume data.

Returns
Vector Vector of nodal volume

Definition at line 168 of file mesh.h.

168{ return d_vol; };

References d_vol.

◆ getNodalVolumesP() [1/2]

std::vector< double > * fe::Mesh::getNodalVolumesP ( )
inline

Get the pointer to nodal volume data.

Returns
pointer Pointer to nodal volume data

Definition at line 180 of file mesh.h.

180{ return &d_vol; };

References d_vol.

◆ getNodalVolumesP() [2/2]

const std::vector< double > * fe::Mesh::getNodalVolumesP ( ) const
inline

Get the pointer to nodal volume data.

Returns
pointer Pointer to nodal volume data

Definition at line 177 of file mesh.h.

177{ return &d_vol; };

References d_vol.

◆ getNode()

util::Point fe::Mesh::getNode ( const size_t i) const
inline

Get coordinates of node i.

Parameters
iId of the node
Returns
coords Coordinates

Definition at line 119 of file mesh.h.

119{ return d_nodes[i]; };

References d_dim, and d_nodes.

◆ getNodes() [1/2]

std::vector< util::Point > & fe::Mesh::getNodes ( )
inline

Get the nodes data.

Returns
nodes Nodes data

Definition at line 135 of file mesh.h.

135{ return d_nodes; };

References d_nodes.

◆ getNodes() [2/2]

const std::vector< util::Point > & fe::Mesh::getNodes ( ) const
inline

Get the nodes data.

Returns
nodes Nodes data

Definition at line 132 of file mesh.h.

132{ return d_nodes; };

References d_nodes.

◆ getNodesP() [1/2]

std::vector< util::Point > * fe::Mesh::getNodesP ( )
inline

Get the pointer to nodes data.

Returns
pointer Pointer to nodes data

Definition at line 144 of file mesh.h.

144{ return &d_nodes; };

References d_nodes.

◆ getNodesP() [2/2]

const std::vector< util::Point > * fe::Mesh::getNodesP ( ) const
inline

Get the pointer to nodes data.

Returns
pointer Pointer to nodes data

Definition at line 141 of file mesh.h.

141{ return &d_nodes; };

References d_nodes.

◆ getNumDofs()

size_t fe::Mesh::getNumDofs ( ) const
inline

Get the number of dofs.

Returns
N Number of dofs

Definition at line 100 of file mesh.h.

100{ return d_numDofs; };

References d_numDofs.

◆ getNumElements()

size_t fe::Mesh::getNumElements ( ) const
inline

Get the number of elements.

Returns
N Number of elements

Definition at line 94 of file mesh.h.

94{ return d_enc.size()/d_eNumVertex; };

References d_enc, and d_eNumVertex.

Referenced by fe::getCurrentQuadPoints(), fe::getMaxShearStressAndLoc(), and fe::getStrainStress().

Here is the caller graph for this function:

◆ getNumNodes()

size_t fe::Mesh::getNumNodes ( ) const
inline

Get the number of nodes.

Returns
N number of nodes

Definition at line 88 of file mesh.h.

88{ return d_numNodes; };

References d_numNodes.

Referenced by fe::getCurrentQuadPoints(), fe::getMaxShearStressAndLoc(), and fe::getStrainStress().

Here is the caller graph for this function:

◆ isNodeFree()

bool fe::Mesh::isNodeFree ( const size_t i,
const unsigned int dof 
) const
inline

Return true if node is free.

Parameters
iId of node
dofDof to check for
Returns
bool True if dof is free else false

Definition at line 188 of file mesh.h.

188 {
189
190 // below checks if d_fix has 1st bit (if dof=0), 2nd bit (if dof=1), 3rd
191 // bit (if dof=2) is set to 1 or 0. If set to 1, then it means it is fixed,
192 // and therefore it returns false
193 return !(d_fix[i] >> dof & 1UL);
194 };

References d_dim, and d_fix.

◆ print()

void fe::Mesh::print ( int  nt = 0,
int  lvl = 0 
) const
inline

Prints the information about the object.

Parameters
ntNumber of tabs to append before printing
lvlInformation level (higher means more information)

Definition at line 306 of file mesh.h.

306{ std::cout << printStr(nt, lvl); };

References d_dim, and printStr().

Here is the call graph for this function:

◆ printStr()

std::string fe::Mesh::printStr ( int  nt = 0,
int  lvl = 0 
) const

Returns the string containing printable information about the object.

Parameters
ntNumber of tabs to append before printing
lvlInformation level (higher means more information)
Returns
string String containing printable information about the object

Definition at line 404 of file mesh.cpp.

404 {
405
406 auto tabS = util::io::getTabS(nt);
407 std::ostringstream oss;
408 oss << tabS << "------- Mesh --------" << std::endl << std::endl;
409 oss << tabS << "Dimension = " << d_dim << std::endl;
410 oss << tabS << "Spatial discretization type = " << d_spatialDiscretization << std::endl;
411 oss << tabS << "Mesh size = " << d_h << std::endl;
412 oss << tabS << "Num nodes = " << d_numNodes << std::endl;
413 oss << tabS << "Num elements = " << d_numElems << std::endl;
414 oss << tabS << "Element type = " << d_eType << std::endl;
415 oss << tabS << "Num nodes per element = " << d_eNumVertex << std::endl;
416 oss << tabS << "Num nodal vol = " << d_vol.size() << std::endl;
417 oss << tabS << "Bounding box: " << std::endl;
419 oss << tabS << std::endl;
420
421 return oss.str();
422}
std::string printBoxStr(const std::pair< util::Point, util::Point > &box, int nt=print_default_tab)
Returns formatted string for output.
Definition io.h:168
std::string getTabS(int nt)
Returns tab spaces of given size.
Definition io.h:40

References util::io::getTabS(), and util::io::printBoxStr().

Referenced by print().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ readElementData()

bool fe::Mesh::readElementData ( const std::string &  filename)

Reads element-node connectivity data from file. This function is meant for cases when mesh was created without element-node connectivity data but later during output, strain/stress were required which needs element-node connectivity data.

File must be either **.vtu** or **.msh** file and must have element-node connectivity data.

Parameters
filenameName of the mesh file

Definition at line 203 of file mesh.cpp.

203 {
204
205 if (d_encDataPopulated and !d_enc.empty()) {
206 util::io::log("Mesh: Element data is populated already.\n");
207 return false;
208 }
209
210 util::io::log("Mesh: Reading element-node connectivity data.\n");
211
212 int file_type = -1;
213 // find the extension of file and call correct reader
215 file_type = 0;
216 else if (util::io::getExtensionFromFile(filename) == "msh")
217 file_type = 1;
218 else if (util::io::getExtensionFromFile(filename) == "vtu")
219 file_type = 2;
220 else {
221 std::cerr << "Error: Currently only '.csv', '.msg', and '.vtu' "
222 "files are supported for reading mesh.\n";
224 }
225
226 if (file_type == 0) {
227 std::cerr << "Error: readElementData() requires file to be either "
228 ".vtu or .msh mesh file.\n";
230 }
231
232 if (file_type == 1) {
234 &d_enc,
235 &d_nec);
236 d_encDataPopulated = true;
237 return true;
238 }
239 else if (file_type == 2) {
241 &d_enc,
242 &d_nec);
243 d_encDataPopulated = true;
244 return true;
245 }
246
247 return false;
248}
void readMshFileCells(const std::string &filename, size_t dim, size_t &element_type, size_t &num_elem, std::vector< size_t > *enc, std::vector< std::vector< size_t > > *nec)
Reads cell data, i.e. element-node connectivity and node-element connectivity.
Definition reader.cpp:401

References util::io::getExtensionFromFile(), util::io::log(), rw::reader::readMshFileCells(), and rw::reader::readVtuFileCells().

Here is the call graph for this function:

◆ setFixity()

void fe::Mesh::setFixity ( const size_t i,
const unsigned int dof,
const bool flag 
)

Set the fixity to free (0) or fixed (1)

Parameters
iId of node
dofDof which is affected
flagSet fixity to fixed if true or free

Definition at line 385 of file mesh.cpp.

386 {
387
388 // to set i^th bit as true of integer a,
389 // a |= 1UL << (i % 8)
390
391 // to set i^th bit as false of integer a,
392 // a &= ~(1UL << (i % 8))
393
394 flag ? (d_fix[i] |= 1UL << dof) : (d_fix[i] &= ~(1UL << dof));
395}

Field Documentation

◆ d_bbox

std::pair<std::vector<double>, std::vector<double> > fe::Mesh::d_bbox

Bounding box.

Definition at line 513 of file mesh.h.

Referenced by fe::createUniformMesh(), getBoundingBox(), and getBoundingBox().

◆ d_dim

◆ d_enc

std::vector<size_t> fe::Mesh::d_enc

Element-node connectivity data.

Structure: First d_eNumVertex data gives the connectivity of first element, and next d_eNumVertex data gives the connectivity of second element and so on and so fourth.

Definition at line 415 of file mesh.h.

Referenced by fe::createUniformMesh(), getElementConnectivities(), getElementConnectivities(), getElementConnectivitiesP(), getElementConnectivitiesP(), getElementConnectivity(), getElementConnectivityNodes(), and getNumElements().

◆ d_encDataPopulated

bool fe::Mesh::d_encDataPopulated

Flag that indicates whether element-node connectivity data is read from file.

Definition at line 487 of file mesh.h.

◆ d_eNumVertex

size_t fe::Mesh::d_eNumVertex

Number of vertex per element.

This information is useful in getting the connectivity for a given element. We assume that the mesh has only one type of elements and based on that assumption we store the element-node connectivity in a single vector.

The value for different elements are

  • Line element: 2,
  • Triangle element: 3,
  • Quadrilateral element: 4,
  • Tetrahedral element: 4

Definition at line 404 of file mesh.h.

Referenced by fe::createUniformMesh(), fe::getCurrentQuadPoints(), getElementConnectivity(), getElementConnectivityNodes(), fe::getMaxShearStressAndLoc(), getNumElements(), and fe::getStrainStress().

◆ d_eType

size_t fe::Mesh::d_eType

Element type.

We follow VTK convention to identify the elements:

  • Line element = 3,
  • Triangle element = 5,
  • Pixel element = 8,
  • Quadrilateral element = 9,
  • Tetrahedral element = 10

Definition at line 389 of file mesh.h.

Referenced by fe::createUniformMesh(), and getElementType().

◆ d_filename

std::string fe::Mesh::d_filename

Filename to read mesh data.

Definition at line 481 of file mesh.h.

Referenced by Mesh().

◆ d_fix

std::vector<uint8_t> fe::Mesh::d_fix

Vector of fixity mask of each node.

First bit represents x-dof, second represents y-dof, and third represents z-dof. 0 represents free dof and 1 represents fixed dof.

We store data in uint8_t type which can hold 8 bit. At present we only use first 3 bits.

Definition at line 431 of file mesh.h.

Referenced by fe::createUniformMesh(), getFixity(), getFixity(), getFixityP(), getFixityP(), and isNodeFree().

◆ d_gInvMap

std::vector<int> fe::Mesh::d_gInvMap

Map from global id to reduced global id.

This is a inverse of d_gMap

Definition at line 510 of file mesh.h.

◆ d_gMap

std::vector<size_t> fe::Mesh::d_gMap

Map from global reduced id to default global id.

We assign number to each free dof where number ranges from 0 to total number of free dofs. This is referred to set of global reduced id. This new set of ids are subset of set of global ids of all dofs, and therefore, "reduced" word is used.

d_gMap provides a map from global reduced id to global id.

Note
Needed only when the discretization is "weak_finite_element" for the assembly of the mass matrix.

Definition at line 504 of file mesh.h.

◆ d_h

double fe::Mesh::d_h

Mesh size.

Definition at line 516 of file mesh.h.

Referenced by fe::createUniformMesh(), and getMeshSize().

◆ d_nec

std::vector<std::vector<size_t> > fe::Mesh::d_nec

Node-element connectivity data.

At present, this data is never populated.

Definition at line 421 of file mesh.h.

◆ d_needEncData

bool fe::Mesh::d_needEncData

Flag that indicates whether we need enc data (set by input mesh deck in constructor)

Definition at line 484 of file mesh.h.

◆ d_nodePartition

std::vector<size_t> fe::Mesh::d_nodePartition

Node partition information. For each node i, d_nodePartition[i] specifies the partition number, i.e., the processor that owns the node in MPI application.

For uniform square mesh, the volume is simply \( h^2 \) in 2-d and \( h^3\) in 3-d, where \( h\) is the mesh size. For general mesh, the volume is computed using the element-node connectivity of the mesh.

Definition at line 463 of file mesh.h.

Referenced by fe::metisGraphPartition().

◆ d_nodes

std::vector<util::Point> fe::Mesh::d_nodes

Vector of initial (reference) coordinates of nodes.

Definition at line 407 of file mesh.h.

Referenced by fe::createUniformMesh(), getElementConnectivityNodes(), getNode(), getNodes(), getNodes(), getNodesP(), and getNodesP().

◆ d_nPart

size_t fe::Mesh::d_nPart

Number of partitions.

Definition at line 449 of file mesh.h.

Referenced by fe::metisGraphPartition().

◆ d_numDofs

size_t fe::Mesh::d_numDofs

Number of dofs = (dimension) times (number of nodes)

Definition at line 490 of file mesh.h.

Referenced by fe::createUniformMesh(), and getNumDofs().

◆ d_numElems

size_t fe::Mesh::d_numElems

Number of elements.

Definition at line 378 of file mesh.h.

Referenced by fe::createUniformMesh().

◆ d_numNodes

size_t fe::Mesh::d_numNodes

Number of nodes.

Definition at line 375 of file mesh.h.

Referenced by fe::createUniformMesh(), and getNumNodes().

◆ d_partitionMethod

std::string fe::Mesh::d_partitionMethod

Partitioning method. It could be either empty string or "metis_recursive" or "metis_kway".

Definition at line 454 of file mesh.h.

Referenced by fe::metisGraphPartition().

◆ d_spatialDiscretization

std::string fe::Mesh::d_spatialDiscretization

Tag for spatial discretization type.

List of valid values are:

  • finite_difference
  • weak_finite_element
  • nodal_finite_element
  • truss_finite_element

Definition at line 478 of file mesh.h.

Referenced by Mesh().

◆ d_vol

std::vector<double> fe::Mesh::d_vol

Vector of volume of each node.

For uniform square mesh, the volume is simply \( h^2 \) in 2-d and \( h^3\) in 3-d, where \( h\) is the mesh size. For general mesh, the volume is computed using the element-node connectivity of the mesh.

Definition at line 439 of file mesh.h.

Referenced by fe::createUniformMesh(), getNodalVolume(), getNodalVolumes(), getNodalVolumes(), getNodalVolumesP(), and getNodalVolumesP().


The documentation for this class was generated from the following files: