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

A vtk writer for simple point data and complex fem mesh data. More...

#include <vtkParticleWriter.h>

Collaboration diagram for rw::writer::VtkParticleWriter:

Public Member Functions

 VtkParticleWriter (const std::string &filename, const std::string &compress_type="")
 Constructor.
 
void addTimeStep (const double &timestep)
 Writes the time step to the file.
 
void close ()
 Closes the file and store it to the hard disk.
 
Mesh data
void appendNodes (const model::ModelData *model, const std::vector< std::string > &tags)
 Writes the nodes to the file.
 
void appendMesh (const model::ModelData *model, const std::vector< std::string > &tags)
 Writes the nodes to the file.
 
void appendContactData (const model::ModelData *model, const std::vector< size_t > *processed_nodes, const std::vector< std::pair< size_t, size_t > > *processed_elems)
 Prepares contact data that is set of nodes in contact and line-element connecting two contacting nodes.
 
void appendStrainStress (const model::ModelData *model)
 Writes strain/stress.
 

Private Attributes

vtkSmartPointer< vtkXMLUnstructuredGridWriter > d_writer_p
 XML unstructured grid writer.
 
vtkSmartPointer< vtkUnstructuredGrid > d_grid_p
 Unstructured grid.
 
std::string d_compressType
 compression_type Specify the compressor (if any)
 

Detailed Description

A vtk writer for simple point data and complex fem mesh data.

Definition at line 28 of file vtkParticleWriter.h.

Constructor & Destructor Documentation

◆ VtkParticleWriter()

rw::writer::VtkParticleWriter::VtkParticleWriter ( const std::string &  filename,
const std::string &  compress_type = "" 
)
explicit

Constructor.

Creates and opens .vtu file of name given by filename. The file remains open till the close() function is invoked.

Parameters
filenameName of file which will be created
compress_typeCompression method (optional)

Definition at line 29 of file vtkParticleWriter.cpp.

31 : d_compressType(compress_type) {
32
33 std::string f = filename + ".vtu";
34
35 d_writer_p = vtkSmartPointer<vtkXMLUnstructuredGridWriter>::New();
36 d_writer_p->SetFileName(const_cast<char *>(f.c_str()));
37}
vtkSmartPointer< vtkXMLUnstructuredGridWriter > d_writer_p
XML unstructured grid writer.
std::string d_compressType
compression_type Specify the compressor (if any)

References d_writer_p.

Member Function Documentation

◆ addTimeStep()

void rw::writer::VtkParticleWriter::addTimeStep ( const double &  timestep)

Writes the time step to the file.

Parameters
timestepCurrent time step of the simulation

Definition at line 341 of file vtkParticleWriter.cpp.

341 {
342
343 auto t = vtkDoubleArray::New();
344 t->SetName("TIME");
345 t->SetNumberOfTuples(1);
346 t->SetTuple1(0, timestep);
347 d_grid_p->GetFieldData()->AddArray(t);
348}
vtkSmartPointer< vtkUnstructuredGrid > d_grid_p
Unstructured grid.

◆ appendContactData()

void rw::writer::VtkParticleWriter::appendContactData ( const model::ModelData model,
const std::vector< size_t > *  processed_nodes,
const std::vector< std::pair< size_t, size_t > > *  processed_elems 
)

Prepares contact data that is set of nodes in contact and line-element connecting two contacting nodes.

Parameters
modelModelData class object
processed_nodesNodes in the list for which we are writing contact data
processed_elemsNumber of line elements (two nodes)

Definition at line 361 of file vtkParticleWriter.cpp.

365 {
366
367 if (processed_nodes->size() == 0)
368 return;
369
370 // write point data
371 auto points = vtkSmartPointer<vtkPoints>::New();
372
373
374 const size_t num_nodes = processed_nodes->size();
375 const size_t num_elems = processed_elems->size();
376
377 // get all the nodes first
378 for (const auto &i : *processed_nodes) {
379 const auto &x = model->d_x[i];
380 points->InsertNextPoint(x.d_x, x.d_y, x.d_z);
381 }
382
383 // write point data
384 d_grid_p = vtkSmartPointer<vtkUnstructuredGrid>::New();
385 d_grid_p->SetPoints(points);
386
387 // now wrtie element data
388 const size_t vtk_element_type = 3; // line element
389 const size_t num_vertex = 2;
390 // element node connectivity
391 auto cells = vtkSmartPointer<vtkCellArray>::New();
392 cells->Allocate(num_vertex, num_elems);
393
394 // element type
395 int cell_types[num_elems];
396
397 vtkIdType ids[num_vertex];
398 for (size_t i = 0; i < num_elems; i++) {
399
400 ids[0] = (*processed_elems)[i].first;
401 ids[1] = (*processed_elems)[i].second;
402
403 cells->InsertNextCell(num_vertex, ids);
404 cell_types[i] = vtk_element_type;
405 }
406
407 // element node connectivity
408 d_grid_p->SetCells(cell_types, cells);
409
410 // write cell data (normal direction)
411 {
412 auto array = vtkSmartPointer < vtkDoubleArray > ::New();
413 array->SetNumberOfComponents(3);
414 array->SetName("Normal");
415 array->SetComponentName(0, "x");
416 array->SetComponentName(1, "y");
417 array->SetComponentName(2, "z");
418
419 double value[3];
420 for (size_t i = 0; i < num_elems; i++) {
421
422 ids[0] = (*processed_elems)[i].first;
423 ids[1] = (*processed_elems)[i].second;
424
425 auto glob_id1 = (*processed_nodes)[ids[0]];
426 auto glob_id2 = (*processed_nodes)[ids[1]];
427
428 const auto &x1 = model->d_x[glob_id1];
429 const auto &x2 = model->d_x[glob_id2];
430
431 auto xd = (x1 - x2)/((x2 - x1).length());
432
433 value[0] = xd[0];
434 value[1] = xd[1];
435 value[2] = xd[2];
436 array->InsertNextTuple(value);
437 }
438
439 d_grid_p->GetCellData()->AddArray(array);
440 }
441}

◆ appendMesh()

void rw::writer::VtkParticleWriter::appendMesh ( const model::ModelData model,
const std::vector< std::string > &  tags 
)

Writes the nodes to the file.

Parameters
modelModelData class object
tagsVector of tags (name of data, e.g., 'Displacement', 'Velocity') to append to the file

Definition at line 275 of file vtkParticleWriter.cpp.

277 {
278
279 if (model->d_x.size() == 0)
280 return;
281
282 // write point data
283 appendNodes(model, tags);
284
285 //
286 // process elements data
287 //
288
289 // get total number of elements and maximum number of vertex in any element
290 size_t num_elems = 0;
291 size_t num_vertex = 0;
292
293 // count number of elements in all particles
294 for (const auto &p : model->d_particlesListTypeAll) {
295 //const auto &rp = p->d_rp_p;
296 num_elems += p->getMeshP()->getNumElements();
297 auto n =
298 util::vtk_map_element_to_num_nodes[p->getMeshP()->getElementType()];
299 if (num_vertex < n)
300 num_vertex = n;
301 }
302
303 // element node connectivity
304 auto cells = vtkSmartPointer<vtkCellArray>::New();
305 cells->Allocate(num_vertex, num_elems);
306
307 // element type
308 int cell_types[num_elems];
309
310 // loop over particles
311 size_t global_elem_counter = 0;
312 for (const auto &p : model->d_particlesListTypeAll) {
313 // get mesh of reference particle in this zone
314 const auto &mesh = p->getMeshP();
315
316 // get element type
317 size_t element_type = mesh->getElementType();
318
319 // loop over elements of this particle
320 size_t num_vertex_p = util::vtk_map_element_to_num_nodes[element_type];
321 vtkIdType ids[num_vertex_p];
322 for (size_t e =0; e < mesh->getNumElements(); e++) {
323 auto elem = mesh->getElementConnectivity(e);
324
325 // assign global ids to the nodes
326 for (size_t n=0; n<elem.size(); n++)
327 ids[n] = elem[n] + p->d_globStart;
328
329 cells->InsertNextCell(num_vertex_p, ids);
330 cell_types[global_elem_counter] = element_type;
331
332 // increment global element counter
333 global_elem_counter++;
334 }
335 }
336
337 // element node connectivity
338 d_grid_p->SetCells(cell_types, cells);
339}
void appendNodes(const model::ModelData *model, const std::vector< std::string > &tags)
Writes the nodes to the file.
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.

◆ appendNodes()

void rw::writer::VtkParticleWriter::appendNodes ( const model::ModelData model,
const std::vector< std::string > &  tags 
)

Writes the nodes to the file.

Parameters
modelModelData class object
tagsVector of tags (name of data, e.g., 'Displacement', 'Velocity') to append to the file

Definition at line 39 of file vtkParticleWriter.cpp.

41 {
42
43 if (model->d_x.size() == 0)
44 return;
45
46 // write point data
47 auto points = vtkSmartPointer<vtkPoints>::New();
48
49 // get all the nodes first
50 for (const auto &x : model->d_x)
51 points->InsertNextPoint(x.d_x, x.d_y, x.d_z);
52
53 // write point data
54 d_grid_p = vtkSmartPointer<vtkUnstructuredGrid>::New();
55 d_grid_p->SetPoints(points);
56
57 // now write data associated to nodes in both particle and wall
58 double value[3];
59 value[0] = 0;
60 value[1] = 0;
61 value[2] = 0;
62 double p_tag[1];
63 p_tag[0] = 0;
64
65 // handle displacement
66 if (util::methods::isTagInList("Displacement", tags)) {
67
68 auto array = vtkSmartPointer<vtkDoubleArray>::New();
69 array->SetNumberOfComponents(3);
70 array->SetName("Displacement");
71
72 array->SetComponentName(0, "x");
73 array->SetComponentName(1, "y");
74 array->SetComponentName(2, "z");
75
76 for (const auto &ui : model->d_u) {
77 value[0] = ui.d_x;
78 value[1] = ui.d_y;
79 value[2] = ui.d_z;
80 array->InsertNextTuple(value);
81 }
82
83 // write
84 d_grid_p->GetPointData()->AddArray(array);
85 } // displacement
86
87 // handle velocity
88 if (util::methods::isTagInList("Velocity", tags)) {
89
90 auto array = vtkSmartPointer<vtkDoubleArray>::New();
91 array->SetNumberOfComponents(3);
92 array->SetName("Velocity");
93
94 array->SetComponentName(0, "x");
95 array->SetComponentName(1, "y");
96 array->SetComponentName(2, "z");
97
98 for (const auto &ui : model->d_v) {
99 value[0] = ui.d_x;
100 value[1] = ui.d_y;
101 value[2] = ui.d_z;
102 array->InsertNextTuple(value);
103 }
104
105 // write
106 d_grid_p->GetPointData()->AddArray(array);
107 } // velocity
108
109 // handle force
110 if (util::methods::isTagInList("Force_Density", tags)) {
111
112 auto array = vtkSmartPointer<vtkDoubleArray>::New();
113 array->SetNumberOfComponents(3);
114 array->SetName("Force_Density");
115
116 array->SetComponentName(0, "x");
117 array->SetComponentName(1, "y");
118 array->SetComponentName(2, "z");
119
120 for (const auto &ui : model->d_f) {
121 value[0] = ui.d_x;
122 value[1] = ui.d_y;
123 value[2] = ui.d_z;
124 array->InsertNextTuple(value);
125 }
126
127 // write
128 d_grid_p->GetPointData()->AddArray(array);
129 } // force
130
131 // handle force
132 if (util::methods::isTagInList("Force", tags)) {
133
134 auto array = vtkSmartPointer<vtkDoubleArray>::New();
135 array->SetNumberOfComponents(3);
136 array->SetName("Force");
137
138 array->SetComponentName(0, "x");
139 array->SetComponentName(1, "y");
140 array->SetComponentName(2, "z");
141
142 size_t i_count = 0;
143 for (const auto &ui : model->d_f) {
144 const auto &voli = model->d_vol[i_count];
145 value[0] = ui.d_x * voli;
146 value[1] = ui.d_y * voli;
147 value[2] = ui.d_z * voli;
148 array->InsertNextTuple(value);
149
150 i_count++;
151 }
152
153 // write
154 d_grid_p->GetPointData()->AddArray(array);
155 } // force
156
157 // handle fixity
158 if (util::methods::isTagInList("Fixity", tags)) {
159
160 auto array = vtkSmartPointer<vtkDoubleArray>::New();
161 array->SetNumberOfComponents(1);
162 array->SetName("Fixity");
163
164 for (const auto &n : model->d_fix) {
165 p_tag[0] = double(n);
166 array->InsertNextTuple(p_tag);
167 }
168
169 // write
170 d_grid_p->GetPointData()->AddArray(array);
171 } // fixity
172
173 // handle Particle ID
174 if (util::methods::isTagInList("Particle_ID", tags)) {
175
176 auto array = vtkSmartPointer<vtkDoubleArray>::New();
177 array->SetNumberOfComponents(1);
178 array->SetName("Particle_ID");
179
180 for (size_t i = 0; i<model->d_x.size(); i++) {
181 auto pi = model->getPtId(i);
182 p_tag[0] = double(model->getParticleFromAllList(pi)->getId());
183 array->InsertNextTuple(p_tag);
184 }
185
186 // write
187 d_grid_p->GetPointData()->AddArray(array);
188 } // Particle ID
189
190 // handle Zone ID
191 if (util::methods::isTagInList("Zone_ID", tags)) {
192
193 auto array = vtkSmartPointer<vtkDoubleArray>::New();
194 array->SetNumberOfComponents(1);
195 array->SetName("Zone_ID");
196
197 for (size_t i = 0; i<model->d_x.size(); i++) {
198 auto pi = model->getPtId(i);
199 p_tag[0] = double(model->getParticleFromAllList(pi)->d_zoneId);
200 array->InsertNextTuple(p_tag);
201 }
202
203 // write
204 d_grid_p->GetPointData()->AddArray(array);
205 } // Zone ID
206
207 // handle force fixity
208 if (util::methods::isTagInList("Force_Fixity", tags)) {
209
210 auto array = vtkSmartPointer<vtkDoubleArray>::New();
211 array->SetNumberOfComponents(1);
212 array->SetName("Force_Fixity");
213
214 for (const auto &n : model->d_forceFixity) {
215 p_tag[0] = double(n);
216 array->InsertNextTuple(p_tag);
217 }
218
219 // write
220 d_grid_p->GetPointData()->AddArray(array);
221 } // force fixity
222
223 // handle nodal volume
224 if (util::methods::isTagInList("Nodal_Volume", tags)) {
225
226 auto array = vtkSmartPointer<vtkDoubleArray>::New();
227 array->SetNumberOfComponents(1);
228 array->SetName("Nodal_Volume");
229
230 for (const auto &n : model->d_vol) {
231 p_tag[0] = double(n);
232 array->InsertNextTuple(p_tag);
233 }
234
235 // write
236 d_grid_p->GetPointData()->AddArray(array);
237 } // nodal volume
238
239 // handle damage_Z
240 if (util::methods::isTagInList("Damage_Z", tags)) {
241
242 auto array = vtkSmartPointer<vtkDoubleArray>::New();
243 array->SetNumberOfComponents(1);
244 array->SetName("Damage_Z");
245
246 for (const auto &n : model->d_Z) {
247 p_tag[0] = double(n);
248 array->InsertNextTuple(p_tag);
249 }
250
251 // write
252 d_grid_p->GetPointData()->AddArray(array);
253 } // damage_Z
254
255 // handle theta
256 if (util::methods::isTagInList("Theta", tags)) {
257
258 if (model->getParticleFromAllList(0)->d_material_p->isStateActive()) {
259
260 auto array = vtkSmartPointer<vtkDoubleArray>::New();
261 array->SetNumberOfComponents(1);
262 array->SetName("Theta");
263
264 for (const auto &n : model->d_thetaX) {
265 p_tag[0] = double(n);
266 array->InsertNextTuple(p_tag);
267 }
268
269 // write
270 d_grid_p->GetPointData()->AddArray(array);
271 }
272 } // Theta
273}
bool isTagInList(const std::string &tag, const std::vector< std::string > &tags)
Returns true if tag is found in the list of tags.
Definition methods.h:279

References util::methods::isTagInList().

Here is the call graph for this function:

◆ appendStrainStress()

void rw::writer::VtkParticleWriter::appendStrainStress ( const model::ModelData model)

Writes strain/stress.

Parameters
modelModelData class object

Definition at line 444 of file vtkParticleWriter.cpp.

445 {
446
447 if (model->d_xQuadCur.size() == 0) {
448 std::cout << "VtkParticleWriter::appendStrainStress: Nothing to write.\n";
449 return;
450 }
451
452 // write point data
453 auto points = vtkSmartPointer<vtkPoints>::New();
454
455 // get all the quadrature points first
456 for (const auto &x : model->d_xQuadCur)
457 points->InsertNextPoint(x.d_x, x.d_y, x.d_z);
458
459 // write point data
460 d_grid_p = vtkSmartPointer<vtkUnstructuredGrid>::New();
461 d_grid_p->SetPoints(points);
462
463 // now write data associated to nodes (in this case, quad points are nodes)
464 double value[3] = {0., 0., 0.};
465 double value_s[6] = {0., 0., 0., 0., 0., 0.};
466 double p_tag[1] = {0.};
467
468 auto array_strain = vtkSmartPointer<vtkDoubleArray>::New();
469 array_strain->SetNumberOfComponents(6);
470 array_strain->SetName("Strain");
471
472 auto array_stress = vtkSmartPointer<vtkDoubleArray>::New();
473 array_stress->SetNumberOfComponents(6);
474 array_stress->SetName("Stress");
475
476 std::vector<std::string> coord_strings = {"xx", "yy", "zz", "yz", "xz", "xy"};
477 for (size_t i =0; i<6; i++) {
478 array_strain->SetComponentName(i, coord_strings[i].c_str());
479 array_stress->SetComponentName(i, coord_strings[i].c_str());
480 }
481
482 for (size_t i=0; i<model->d_strain.size(); i++) {
483
484 model->d_strain[i].copy(value_s);
485 array_strain->InsertNextTuple(value_s);
486
487 model->d_stress[i].copy(value_s);
488 array_stress->InsertNextTuple(value_s);
489 }
490
491
492 // write
493 d_grid_p->GetPointData()->AddArray(array_strain);
494 d_grid_p->GetPointData()->AddArray(array_stress);
495}

◆ close()

void rw::writer::VtkParticleWriter::close ( )

Closes the file and store it to the hard disk.

Definition at line 350 of file vtkParticleWriter.cpp.

350 {
351 d_writer_p->SetInputData(d_grid_p);
352 d_writer_p->SetDataModeToAppended();
353 d_writer_p->EncodeAppendedDataOn();
354 if (d_compressType == "zlib")
355 d_writer_p->SetCompressorTypeToZLib();
356 else
357 d_writer_p->SetCompressor(0);
358 d_writer_p->Write();
359}

Field Documentation

◆ d_compressType

std::string rw::writer::VtkParticleWriter::d_compressType
private

compression_type Specify the compressor (if any)

Definition at line 103 of file vtkParticleWriter.h.

◆ d_grid_p

vtkSmartPointer<vtkUnstructuredGrid> rw::writer::VtkParticleWriter::d_grid_p
private

Unstructured grid.

Definition at line 100 of file vtkParticleWriter.h.

◆ d_writer_p

vtkSmartPointer<vtkXMLUnstructuredGridWriter> rw::writer::VtkParticleWriter::d_writer_p
private

XML unstructured grid writer.

Definition at line 97 of file vtkParticleWriter.h.

Referenced by VtkParticleWriter().


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