PeriDEM 0.2.0
PeriDEM -- Peridynamics-based high-fidelity model for granular media
Loading...
Searching...
No Matches
rw::reader::VtkReader Class Reference

A class to read VTK (.vtu) mesh files. More...

#include <vtkReader.h>

Collaboration diagram for rw::reader::VtkReader:

Public Member Functions

 VtkReader (const std::string &filename)
 Constructor.
 
void readMesh (size_t dim, std::vector< util::Point > *nodes, size_t &element_type, size_t &num_elems, 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.
 
void readNodes (std::vector< util::Point > *nodes)
 Reads nodal position.
 
void readCells (size_t dim, size_t &element_type, size_t &num_elems, std::vector< size_t > *enc, std::vector< std::vector< size_t > > *nec)
 Reads cell data, i.e. element-node connectivity and node-element connectivity.
 
bool readPointData (const std::string &name, std::vector< uint8_t > *data)
 reads point data from .vtu file
 
bool readPointData (const std::string &name, std::vector< size_t > *data)
 reads point data from .vtu file
 
bool readPointData (const std::string &name, std::vector< int > *data)
 reads point data from .vtu file
 
bool readPointData (const std::string &name, std::vector< float > *data)
 reads point data from .vtu file
 
bool readPointData (const std::string &name, std::vector< double > *data)
 reads point data from .vtu file
 
bool readPointData (const std::string &name, std::vector< util::Point > *data)
 reads point data from .vtu file
 
bool readPointData (const std::string &name, std::vector< util::SymMatrix3 > *data)
 reads point data from .vtu file
 
bool readPointData (const std::string &name, std::vector< util::Matrix3 > *data)
 reads point data from .vtu file
 
bool readCellData (const std::string &name, std::vector< float > *data)
 reads cell data from .vtu file
 
bool readCellData (const std::string &name, std::vector< double > *data)
 reads cell data from .vtu file
 
bool readCellData (const std::string &name, std::vector< util::Point > *data)
 reads cell data from .vtu file
 
bool readCellData (const std::string &name, std::vector< util::SymMatrix3 > *data)
 reads cell data from .vtu file
 
bool readCellData (const std::string &name, std::vector< util::Matrix3 > *data)
 reads cell data from .vtu file
 
void close ()
 Close the file.
 

Private Attributes

vtkSmartPointer< vtkXMLUnstructuredGridReader > d_reader_p
 XML unstructured grid writer.
 
vtkSmartPointer< vtkUnstructuredGrid > d_grid_p
 Unstructured grid.
 

Static Private Attributes

static size_t d_count = 0
 Counter.
 

Detailed Description

A class to read VTK (.vtu) mesh files.

Note
Depends on VTK library.

Definition at line 32 of file vtkReader.h.

Constructor & Destructor Documentation

◆ VtkReader()

rw::reader::VtkReader::VtkReader ( const std::string &  filename)
explicit

Constructor.

Parameters
filenameName of mesh file (with/without .vtu extension)

Definition at line 30 of file vtkReader.cpp.

30 {
31 d_count++;
32
33 auto f = util::io::checkAndCreateNewFilename(filename, "vtu");
34
35 // Append the extension vtu to file_name
36 d_reader_p = vtkSmartPointer<vtkXMLUnstructuredGridReader>::New();
37 d_reader_p->SetFileName(const_cast<char *>(f.c_str()));
38 d_reader_p->Update();
39}
static size_t d_count
Counter.
Definition vtkReader.h:133
vtkSmartPointer< vtkXMLUnstructuredGridReader > d_reader_p
XML unstructured grid writer.
Definition vtkReader.h:136
std::string checkAndCreateNewFilename(std::string const &filename, std::string filename_ext)
Check for extension and if possible create new filename from a given filename and a given extension.
Definition io.h:443

References util::io::checkAndCreateNewFilename(), d_count, and d_reader_p.

Here is the call graph for this function:

Member Function Documentation

◆ close()

void rw::reader::VtkReader::close ( )

Close the file.

Definition at line 616 of file vtkReader.cpp.

616 {
617 // delete d_reader_p;
618 // delete d_grid_p;
619}

◆ readCellData() [1/5]

bool rw::reader::VtkReader::readCellData ( const std::string &  name,
std::vector< double > *  data 
)

reads cell data from .vtu file

Parameters
nameName of data
dataPointer to the vector of data
Returns
status True if data is found otherwise false

Definition at line 494 of file vtkReader.cpp.

495 {
496 // read point field data
497 d_grid_p = d_reader_p->GetOutput();
498 vtkCellData *c_field = d_grid_p->GetCellData();
499
500 // handle for displacement, fixity and node element connectivity
501 if (c_field->HasArray(name.c_str()) == 0) return false;
502
503 vtkDataArray *array = c_field->GetArray(name.c_str());
504
505 // Below is not efficient. Later this can be improved.
506 // declare another array data to hold the ux,uy,uz
507 auto data_a = vtkSmartPointer<vtkDoubleArray>::New();
508 data_a->SetNumberOfComponents(1);
509 data_a->Allocate(1, 1); // allocate memory
510
511 (*data).resize(array->GetNumberOfTuples());
512 for (size_t i = 0; i < array->GetNumberOfTuples(); i++) {
513 array->GetTuples(i, i, data_a);
514 (*data)[i] = data_a->GetValue(0);
515 }
516
517 return true;
518}
vtkSmartPointer< vtkUnstructuredGrid > d_grid_p
Unstructured grid.
Definition vtkReader.h:139

◆ readCellData() [2/5]

bool rw::reader::VtkReader::readCellData ( const std::string &  name,
std::vector< float > *  data 
)

reads cell data from .vtu file

Parameters
nameName of data
dataPointer to the vector of data
Returns
status True if data is found otherwise false

Definition at line 468 of file vtkReader.cpp.

469 {
470 // read point field data
471 d_grid_p = d_reader_p->GetOutput();
472 vtkCellData *c_field = d_grid_p->GetCellData();
473
474 // handle for displacement, fixity and node element connectivity
475 if (c_field->HasArray(name.c_str()) == 0) return false;
476
477 vtkDataArray *array = c_field->GetArray(name.c_str());
478
479 // Below is not efficient. Later this can be improved.
480 // declare another array data to hold the ux,uy,uz
481 auto data_a = vtkSmartPointer<vtkDoubleArray>::New();
482 data_a->SetNumberOfComponents(1);
483 data_a->Allocate(1, 1); // allocate memory
484
485 (*data).resize(array->GetNumberOfTuples());
486 for (size_t i = 0; i < array->GetNumberOfTuples(); i++) {
487 array->GetTuples(i, i, data_a);
488 (*data)[i] = data_a->GetValue(0);
489 }
490
491 return true;
492}

◆ readCellData() [3/5]

bool rw::reader::VtkReader::readCellData ( const std::string &  name,
std::vector< util::Matrix3 > *  data 
)

reads cell data from .vtu file

Parameters
nameName of data
dataPointer to the vector of data
Returns
status True if data is found otherwise false

Definition at line 575 of file vtkReader.cpp.

576 {
577 // read point field data
578 d_grid_p = d_reader_p->GetOutput();
579 vtkCellData *c_field = d_grid_p->GetCellData();
580
581 // handle for displacement, fixity and node element connectivity
582 if (c_field->HasArray(name.c_str()) == 0) return false;
583
584 vtkDataArray *array = c_field->GetArray(name.c_str());
585
586 // Below is not efficient. Later this can be improved.
587 // declare another array data to hold the ux,uy,uz
588 auto data_a = vtkSmartPointer<vtkDoubleArray>::New();
589 data_a->SetNumberOfComponents(6);
590 data_a->Allocate(6, 1); // allocate memory
591
592 (*data).resize(array->GetNumberOfTuples());
593 for (size_t i = 0; i < array->GetNumberOfTuples(); i++) {
594 array->GetTuples(i, i, data_a);
595
597 m(0,0) = data_a->GetValue(0);
598 m(1,1) = data_a->GetValue(1);
599 m(2,2) = data_a->GetValue(2);
600 m(1,2) = data_a->GetValue(3);
601 m(0,2) = data_a->GetValue(4);
602 m(0,1) = data_a->GetValue(5);
603
604 // symmetrize
605 m(2,1) = m(1,2);
606 m(2,0) = m(0,2);
607 m(1,0) = m(0,1);
608
609 // add to the data
610 (*data)[i] = m;
611 }
612
613 return true;
614}
A structure to represent 3d matrices.
Definition matrix.h:19

◆ readCellData() [4/5]

bool rw::reader::VtkReader::readCellData ( const std::string &  name,
std::vector< util::Point > *  data 
)

reads cell data from .vtu file

Parameters
nameName of data
dataPointer to the vector of data
Returns
status True if data is found otherwise false

Definition at line 520 of file vtkReader.cpp.

521 {
522 // read point field data
523 d_grid_p = d_reader_p->GetOutput();
524 vtkCellData *c_field = d_grid_p->GetCellData();
525
526 // handle for displacement, fixity and node element connectivity
527 if (c_field->HasArray(name.c_str()) == 0) return false;
528
529 vtkDataArray *array = c_field->GetArray(name.c_str());
530
531 // Below is not efficient. Later this can be improved.
532 // declare another array data to hold the ux,uy,uz
533 auto data_a = vtkSmartPointer<vtkDoubleArray>::New();
534 data_a->SetNumberOfComponents(3);
535 data_a->Allocate(3, 1); // allocate memory
536
537 (*data).resize(array->GetNumberOfTuples());
538 for (size_t i = 0; i < array->GetNumberOfTuples(); i++) {
539 array->GetTuples(i, i, data_a);
540 (*data)[i] = util::Point(data_a->GetValue(0), data_a->GetValue(1),
541 data_a->GetValue(2));
542 }
543
544 return true;
545}
A structure to represent 3d vectors.
Definition point.h:30

◆ readCellData() [5/5]

bool rw::reader::VtkReader::readCellData ( const std::string &  name,
std::vector< util::SymMatrix3 > *  data 
)

reads cell data from .vtu file

Parameters
nameName of data
dataPointer to the vector of data
Returns
status True if data is found otherwise false

Definition at line 547 of file vtkReader.cpp.

548 {
549 // read point field data
550 d_grid_p = d_reader_p->GetOutput();
551 vtkCellData *c_field = d_grid_p->GetCellData();
552
553 // handle for displacement, fixity and node element connectivity
554 if (c_field->HasArray(name.c_str()) == 0) return false;
555
556 vtkDataArray *array = c_field->GetArray(name.c_str());
557
558 // Below is not efficient. Later this can be improved.
559 // declare another array data to hold the ux,uy,uz
560 auto data_a = vtkSmartPointer<vtkDoubleArray>::New();
561 data_a->SetNumberOfComponents(6);
562 data_a->Allocate(6, 1); // allocate memory
563
564 (*data).resize(array->GetNumberOfTuples());
565 for (size_t i = 0; i < array->GetNumberOfTuples(); i++) {
566 array->GetTuples(i, i, data_a);
567 (*data)[i] = util::SymMatrix3({data_a->GetValue(0), data_a->GetValue(1),
568 data_a->GetValue(2), data_a->GetValue(3),
569 data_a->GetValue(4), data_a->GetValue(5)});
570 }
571
572 return true;
573}
A structure to represent 3d matrices.
Definition matrix.h:258

◆ readCells()

void rw::reader::VtkReader::readCells ( size_t  dim,
size_t &  element_type,
size_t &  num_elems,
std::vector< size_t > *  enc,
std::vector< std::vector< size_t > > *  nec 
)

Reads cell data, i.e. element-node connectivity and node-element connectivity.

Parameters
dimDimension
element_typeType of element
num_elemsNumber of elements
encElement-node connectivity
necNode-element connectivity

Definition at line 183 of file vtkReader.cpp.

186 {
187 d_grid_p = d_reader_p->GetOutput();
188 vtkIdType num_nodes = d_grid_p->GetNumberOfPoints();
189 vtkIdType num_elems = d_grid_p->GetNumberOfCells();
190
191 // set number of elements
192 num_elem = num_elems;
193
194 // resize data (we resize enc data inside loop, see below)
195 nec->clear();
196 nec->resize(num_nodes);
197
198 // get cell types array
199 vtkUnsignedCharArray *elem_types = d_grid_p->GetCellTypesArray();
200
201 // read cell field data
202 vtkCellData *c_field = d_grid_p->GetCellData();
203
204 // to hold integer
205 auto val = vtkSmartPointer<vtkUnsignedIntArray>::New();
206 val->SetNumberOfComponents(1);
207 val->Allocate(1, 1); // allocate memory
208
209 int nds_per_el = 0;
210 for (size_t i = 0; i < num_elems; i++) {
211 vtkIdType id = i;
212 vtkIdType num_ids;
213 vtkIdType const *nodes_ids;
214
215 d_grid_p->GetCellPoints(id, num_ids, nodes_ids);
216
217 // resize element-node connectivity data in the beginning of the loop
218 if (i == 0) {
219 elem_types->GetTuples(i, i, val);
220 element_type = val->GetValue(0);
221
222 nds_per_el = util::vtk_map_element_to_num_nodes[element_type];
223
224 // std::cout << "nodes per element = " << nds_per_el << "\n";
225 // std::cout << "Num elements = " << num_elems << "\n";
226
227 // resize the element-node connectivity
228 enc->clear();
229 enc->resize(nds_per_el * num_elems);
230 }
231
232 for (size_t j = 0; j < num_ids; j++) {
233 // read element-node connectivity data
234 (*enc)[nds_per_el * i + j] = nodes_ids[j];
235
236 // update node-element connectivity data
237 (*nec)[nodes_ids[j]].push_back(i);
238 }
239 }
240}
static int vtk_map_element_to_num_nodes[16]
Map from element type to number of nodes (for vtk)

References util::vtk_map_element_to_num_nodes.

◆ readMesh()

void rw::reader::VtkReader::readMesh ( size_t  dim,
std::vector< util::Point > *  nodes,
size_t &  element_type,
size_t &  num_elems,
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.

Parameters
dimDimension
nodesVector of nodes data
element_typeType of element
num_elemsNumber of elements
encVector holding element-node connectivity
necVector holding node-element connectivity
volumesVector holding volume of the nodes
is_fdFlag indicating if this mesh is for finite_difference simulation

Definition at line 41 of file vtkReader.cpp.

46 {
47 //
48 // For nodes, we read following data
49 // 1. nodes initial configuration
50 // 2. nodes volume (if provided)
51 // 3. nodes element connectivity (if provided)
52 //
53 // For elements, we read following data
54 // 1. element node connectivity
55 // 2. element type
56 //
57 d_grid_p = d_reader_p->GetOutput();
58 vtkIdType num_nodes = d_grid_p->GetNumberOfPoints();
59 vtkIdType num_elems = d_grid_p->GetNumberOfCells();
60
61 // read point field data
62 vtkPointData *p_field = d_grid_p->GetPointData();
63
64 // check if file contains nodal volumes
65 bool has_volume = true;
66 if (p_field->HasArray("Node_Volume") == 0 and
67 p_field->HasArray("Volume") == 0)
68 has_volume = false;
69 vtkDataArray *vol_array;
70 if (has_volume) {
71 if (p_field->HasArray("Node_Volume"))
72 vol_array = p_field->GetArray("Node_Volume");
73 else if (p_field->HasArray("Volume"))
74 vol_array = p_field->GetArray("Volume");
75 }
76
77 // resize data
78 nodes->resize(num_nodes);
79 if (has_volume) volumes->resize(num_nodes);
80
81 // Below is not efficient and need improvement
82 // declare another array data to hold the ux,uy,uz
83 auto u_a = vtkSmartPointer<vtkDoubleArray>::New();
84 u_a->SetNumberOfComponents(3);
85 u_a->Allocate(3, 1); // allocate memory
86
87 auto fix_a = vtkSmartPointer<vtkUnsignedIntArray>::New();
88 fix_a->SetNumberOfComponents(1);
89 fix_a->Allocate(1, 1); // allocate memory
90
91 auto vol_a = vtkSmartPointer<vtkDoubleArray>::New();
92 vol_a->SetNumberOfComponents(1);
93 vol_a->Allocate(1, 1); // allocate memory
94
95 auto con_a = vtkSmartPointer<vtkIntArray>::New();
96 con_a->SetNumberOfComponents(11);
97 con_a->Allocate(11, 1); // allocate memory
98
99 for (size_t i = 0; i < num_nodes; i++) {
100 vtkIdType id = i;
101
102 double x[3];
103 d_grid_p->GetPoint(id, x);
104
105 util::Point i_node = util::Point(x[0], x[1], x[2]);
106 (*nodes)[i] = i_node;
107
108 // read node volume
109 if (has_volume) {
110 vol_array->GetTuples(i, i, vol_a);
111 (*volumes)[i] = vol_a->GetValue(0);
112 }
113 }
114
115 // if mesh is for finite difference simulation and if we have read the
116 // volume data then we do not need to go further
117 if (is_fd and has_volume) return;
118
119 // read elements
120 // to resize element-node connectivity, we need to know the number of
121 // vertex in the given element type
122 num_elem = num_elems;
123
124 // resize node-element connectivity
125 nec->resize(num_nodes);
126
127 // get cell types array
128 vtkUnsignedCharArray *elem_types = d_grid_p->GetCellTypesArray();
129
130 // read cell field data
131 vtkCellData *c_field = d_grid_p->GetCellData();
132
133 int nds_per_el = 0;
134 for (size_t i = 0; i < num_elems; i++) {
135 vtkIdType id = i;
136 vtkIdType num_ids;
137 vtkIdType const *nodes_ids;
138
139 d_grid_p->GetCellPoints(id, num_ids, nodes_ids);
140
141 // resize element-node connectivity data in the beginning of the loop
142 if (i == 0) {
143 elem_types->GetTuples(i, i, fix_a);
144 element_type = fix_a->GetValue(0);
145
146 nds_per_el = util::vtk_map_element_to_num_nodes[element_type];
147
148 // std::cout << "nodes per element = " << nds_per_el << "\n";
149 // std::cout << "Num elements = " << num_elems << "\n";
150
151 // resize the element-node connectivity
152 enc->resize(nds_per_el * num_elems);
153 }
154
155 for (size_t j = 0; j < num_ids; j++) {
156 // read element-node connectivity data
157 (*enc)[nds_per_el * i + j] = nodes_ids[j];
158
159 // update node-element connectivity data
160 (*nec)[nodes_ids[j]].push_back(i);
161 }
162 }
163}

References util::vtk_map_element_to_num_nodes.

◆ readNodes()

void rw::reader::VtkReader::readNodes ( std::vector< util::Point > *  nodes)

Reads nodal position.

Parameters
nodesVector of nodal coordinates

Definition at line 165 of file vtkReader.cpp.

165 {
166 d_grid_p = d_reader_p->GetOutput();
167 vtkIdType num_nodes = d_grid_p->GetNumberOfPoints();
168
169 // resize data
170 nodes->clear();
171 nodes->resize(num_nodes);
172 for (size_t i = 0; i < num_nodes; i++) {
173 vtkIdType id = i;
174
175 double x[3];
176 d_grid_p->GetPoint(id, x);
177
178 util::Point i_node = util::Point(x[0], x[1], x[2]);
179 (*nodes)[i] = i_node;
180 }
181}

◆ readPointData() [1/8]

bool rw::reader::VtkReader::readPointData ( const std::string &  name,
std::vector< double > *  data 
)

reads point data from .vtu file

Parameters
nameName of data
dataPointer to the vector of data
Returns
status True if data is found otherwise false

Definition at line 346 of file vtkReader.cpp.

347 {
348 // read point field data
349 d_grid_p = d_reader_p->GetOutput();
350 vtkPointData *p_field = d_grid_p->GetPointData();
351
352 // handle for displacement, fixity and node element connectivity
353 if (p_field->HasArray(name.c_str()) == 0) return false;
354
355 vtkDataArray *array = p_field->GetArray(name.c_str());
356
357 // Below is not efficient. Later this can be improved.
358 // declare another array data to hold the ux,uy,uz
359 auto data_a = vtkSmartPointer<vtkDoubleArray>::New();
360 data_a->SetNumberOfComponents(1);
361 data_a->Allocate(1, 1); // allocate memory
362
363 (*data).resize(array->GetNumberOfTuples());
364 for (size_t i = 0; i < array->GetNumberOfTuples(); i++) {
365 array->GetTuples(i, i, data_a);
366 (*data)[i] = data_a->GetValue(0);
367 }
368
369 return true;
370}

◆ readPointData() [2/8]

bool rw::reader::VtkReader::readPointData ( const std::string &  name,
std::vector< float > *  data 
)

reads point data from .vtu file

Parameters
nameName of data
dataPointer to the vector of data
Returns
status True if data is found otherwise false

Definition at line 320 of file vtkReader.cpp.

321 {
322 // read point field data
323 d_grid_p = d_reader_p->GetOutput();
324 vtkPointData *p_field = d_grid_p->GetPointData();
325
326 // handle for displacement, fixity and node element connectivity
327 if (p_field->HasArray(name.c_str()) == 0) return false;
328
329 vtkDataArray *array = p_field->GetArray(name.c_str());
330
331 // Below is not efficient. Later this can be improved.
332 // declare another array data to hold the ux,uy,uz
333 auto data_a = vtkSmartPointer<vtkDoubleArray>::New();
334 data_a->SetNumberOfComponents(1);
335 data_a->Allocate(1, 1); // allocate memory
336
337 (*data).resize(array->GetNumberOfTuples());
338 for (size_t i = 0; i < array->GetNumberOfTuples(); i++) {
339 array->GetTuples(i, i, data_a);
340 (*data)[i] = data_a->GetValue(0);
341 }
342
343 return true;
344}

◆ readPointData() [3/8]

bool rw::reader::VtkReader::readPointData ( const std::string &  name,
std::vector< int > *  data 
)

reads point data from .vtu file

Parameters
nameName of data
dataPointer to the vector of data
Returns
status True if data is found otherwise false

Definition at line 294 of file vtkReader.cpp.

295 {
296 // read point field data
297 d_grid_p = d_reader_p->GetOutput();
298 vtkPointData *p_field = d_grid_p->GetPointData();
299
300 // handle for displacement, fixity and node element connectivity
301 if (p_field->HasArray(name.c_str()) == 0) return false;
302
303 vtkDataArray *array = p_field->GetArray(name.c_str());
304
305 // Below is not efficient. Later this can be improved.
306 // declare another array data to hold the ux,uy,uz
307 auto data_a = vtkSmartPointer<vtkDoubleArray>::New();
308 data_a->SetNumberOfComponents(1);
309 data_a->Allocate(1, 1); // allocate memory
310
311 (*data).resize(array->GetNumberOfTuples());
312 for (size_t i = 0; i < array->GetNumberOfTuples(); i++) {
313 array->GetTuples(i, i, data_a);
314 (*data)[i] = data_a->GetValue(0);
315 }
316
317 return true;
318}

◆ readPointData() [4/8]

bool rw::reader::VtkReader::readPointData ( const std::string &  name,
std::vector< size_t > *  data 
)

reads point data from .vtu file

Parameters
nameName of data
dataPointer to the vector of data
Returns
status True if data is found otherwise false

Definition at line 268 of file vtkReader.cpp.

269 {
270 // read point field data
271 d_grid_p = d_reader_p->GetOutput();
272 vtkPointData *p_field = d_grid_p->GetPointData();
273
274 // handle for displacement, fixity and node element connectivity
275 if (p_field->HasArray(name.c_str()) == 0) return false;
276
277 vtkDataArray *array = p_field->GetArray(name.c_str());
278
279 // Below is not efficient. Later this can be improved.
280 // declare another array data to hold the ux,uy,uz
281 auto data_a = vtkSmartPointer<vtkDoubleArray>::New();
282 data_a->SetNumberOfComponents(1);
283 data_a->Allocate(1, 1); // allocate memory
284
285 (*data).resize(array->GetNumberOfTuples());
286 for (size_t i = 0; i < array->GetNumberOfTuples(); i++) {
287 array->GetTuples(i, i, data_a);
288 (*data)[i] = data_a->GetValue(0);
289 }
290
291 return true;
292}

◆ readPointData() [5/8]

bool rw::reader::VtkReader::readPointData ( const std::string &  name,
std::vector< uint8_t > *  data 
)

reads point data from .vtu file

Parameters
nameName of data
dataPointer to the vector of data
Returns
status True if data is found otherwise false

Definition at line 242 of file vtkReader.cpp.

243 {
244 // read point field data
245 d_grid_p = d_reader_p->GetOutput();
246 vtkPointData *p_field = d_grid_p->GetPointData();
247
248 // handle for displacement, fixity and node element connectivity
249 if (p_field->HasArray(name.c_str()) == 0) return false;
250
251 vtkDataArray *array = p_field->GetArray(name.c_str());
252
253 // Below is not efficient. Later this can be improved.
254 // declare another array data to hold the ux,uy,uz
255 auto data_a = vtkSmartPointer<vtkDoubleArray>::New();
256 data_a->SetNumberOfComponents(1);
257 data_a->Allocate(1, 1); // allocate memory
258
259 (*data).resize(array->GetNumberOfTuples());
260 for (size_t i = 0; i < array->GetNumberOfTuples(); i++) {
261 array->GetTuples(i, i, data_a);
262 (*data)[i] = data_a->GetValue(0);
263 }
264
265 return true;
266}

◆ readPointData() [6/8]

bool rw::reader::VtkReader::readPointData ( const std::string &  name,
std::vector< util::Matrix3 > *  data 
)

reads point data from .vtu file

Parameters
nameName of data
dataPointer to the vector of data
Returns
status True if data is found otherwise false

Definition at line 427 of file vtkReader.cpp.

428 {
429 // read point field data
430 d_grid_p = d_reader_p->GetOutput();
431 vtkPointData *p_field = d_grid_p->GetPointData();
432
433 // handle for displacement, fixity and node element connectivity
434 if (p_field->HasArray(name.c_str()) == 0) return false;
435
436 vtkDataArray *array = p_field->GetArray(name.c_str());
437
438 // Below is not efficient. Later this can be improved.
439 // declare another array data to hold the ux,uy,uz
440 auto data_a = vtkSmartPointer<vtkDoubleArray>::New();
441 data_a->SetNumberOfComponents(6);
442 data_a->Allocate(6, 1); // allocate memory
443
444 (*data).resize(array->GetNumberOfTuples());
445 for (size_t i = 0; i < array->GetNumberOfTuples(); i++) {
446 array->GetTuples(i, i, data_a);
447
449 m(0,0) = data_a->GetValue(0);
450 m(1,1) = data_a->GetValue(1);
451 m(2,2) = data_a->GetValue(2);
452 m(1,2) = data_a->GetValue(3);
453 m(0,2) = data_a->GetValue(4);
454 m(0,1) = data_a->GetValue(5);
455
456 // symmetrize
457 m(2,1) = m(1,2);
458 m(2,0) = m(0,2);
459 m(1,0) = m(0,1);
460
461 // add to the data
462 (*data)[i] = m;
463 }
464
465 return true;
466}

◆ readPointData() [7/8]

bool rw::reader::VtkReader::readPointData ( const std::string &  name,
std::vector< util::Point > *  data 
)

reads point data from .vtu file

Parameters
nameName of data
dataPointer to the vector of data
Returns
status True if data is found otherwise false

Definition at line 372 of file vtkReader.cpp.

373 {
374 // read point field data
375 d_grid_p = d_reader_p->GetOutput();
376 vtkPointData *p_field = d_grid_p->GetPointData();
377
378 // handle for displacement, fixity and node element connectivity
379 if (p_field->HasArray(name.c_str()) == 0) return false;
380
381 vtkDataArray *array = p_field->GetArray(name.c_str());
382
383 // Below is not efficient. Later this can be improved.
384 // declare another array data to hold the ux,uy,uz
385 auto data_a = vtkSmartPointer<vtkDoubleArray>::New();
386 data_a->SetNumberOfComponents(3);
387 data_a->Allocate(3, 1); // allocate memory
388
389 (*data).resize(array->GetNumberOfTuples());
390 for (size_t i = 0; i < array->GetNumberOfTuples(); i++) {
391 array->GetTuples(i, i, data_a);
392 (*data)[i] = util::Point(data_a->GetValue(0), data_a->GetValue(1),
393 data_a->GetValue(2));
394 }
395
396 return true;
397}

◆ readPointData() [8/8]

bool rw::reader::VtkReader::readPointData ( const std::string &  name,
std::vector< util::SymMatrix3 > *  data 
)

reads point data from .vtu file

Parameters
nameName of data
dataPointer to the vector of data
Returns
status True if data is found otherwise false

Definition at line 399 of file vtkReader.cpp.

400 {
401 // read point field data
402 d_grid_p = d_reader_p->GetOutput();
403 vtkPointData *p_field = d_grid_p->GetPointData();
404
405 // handle for displacement, fixity and node element connectivity
406 if (p_field->HasArray(name.c_str()) == 0) return false;
407
408 vtkDataArray *array = p_field->GetArray(name.c_str());
409
410 // Below is not efficient. Later this can be improved.
411 // declare another array data to hold the ux,uy,uz
412 auto data_a = vtkSmartPointer<vtkDoubleArray>::New();
413 data_a->SetNumberOfComponents(6);
414 data_a->Allocate(6, 1); // allocate memory
415
416 (*data).resize(array->GetNumberOfTuples());
417 for (size_t i = 0; i < array->GetNumberOfTuples(); i++) {
418 array->GetTuples(i, i, data_a);
419 (*data)[i] = util::SymMatrix3({data_a->GetValue(0), data_a->GetValue(1),
420 data_a->GetValue(2), data_a->GetValue(3),
421 data_a->GetValue(4), data_a->GetValue(5)});
422 }
423
424 return true;
425}

Field Documentation

◆ d_count

size_t rw::reader::VtkReader::d_count = 0
staticprivate

Counter.

Definition at line 133 of file vtkReader.h.

Referenced by VtkReader().

◆ d_grid_p

vtkSmartPointer<vtkUnstructuredGrid> rw::reader::VtkReader::d_grid_p
private

Unstructured grid.

Definition at line 139 of file vtkReader.h.

◆ d_reader_p

vtkSmartPointer<vtkXMLUnstructuredGridReader> rw::reader::VtkReader::d_reader_p
private

XML unstructured grid writer.

Definition at line 136 of file vtkReader.h.

Referenced by VtkReader().


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