PeriDEM 0.2.0
PeriDEM -- Peridynamics-based high-fidelity model for granular media
Loading...
Searching...
No Matches
vtkParticleReader.cpp
Go to the documentation of this file.
1/*
2 * -------------------------------------------
3 * Copyright (c) 2021 - 2024 Prashant K. Jha
4 * -------------------------------------------
5 * PeriDEM https://github.com/prashjha/PeriDEM
6 *
7 * Distributed under the Boost Software License, Version 1.0. (See accompanying
8 * file LICENSE)
9 */
10
11#include "vtkParticleReader.h"
12#include "util/feElementDefs.h"
13#include "util/io.h"
14#include "fe/mesh.h"
15#include "util/methods.h"
16#include "model/modelData.h"
18
19#include <cstdint>
20
21#include <vtkAbstractArray.h>
22#include <vtkCellArray.h>
23#include <vtkCellData.h>
24#include <vtkDoubleArray.h>
25#include <vtkIdList.h>
26#include <vtkIntArray.h>
27#include <vtkPointData.h>
28#include <vtkPoints.h>
29#include <vtkUnsignedCharArray.h>
30#include <vtkUnsignedIntArray.h>
31
33
34 auto f = util::io::checkAndCreateNewFilename(filename, "vtu");
35
36 d_reader_p = vtkSmartPointer<vtkXMLUnstructuredGridReader>::New();
37 d_reader_p->SetFileName(const_cast<char *>(f.c_str()));
38 d_reader_p->Update();
39}
40
43
44 if (model->d_x.size() == 0)
45 return;
46
47 //
48 // For nodes, we read following data
49 // 1. nodes current configuration
50 // 2. nodes displacement (if not then compute from current and reference
51 // config)
52 // 3. nodes velocity (if not then produce error)
53 //
54 d_grid_p = d_reader_p->GetOutput();
55 vtkIdType num_nodes = d_grid_p->GetNumberOfPoints();
56
57 if (num_nodes != model->d_x.size()) {
58 std::cerr << "Error: Number of points data = " << num_nodes <<
59 " in file does not match the number of points = "
60 << model->d_x.size() << ". Aborting reading.\n";
61 exit(1);
62 }
63
64 // read current configuration
65 for (size_t i = 0; i < num_nodes; i++) {
66 vtkIdType id = i;
67
68 double x[3];
69 d_grid_p->GetPoint(id, x);
70
71 util::Point i_node = util::Point(x[0], x[1], x[2]);
72 model->setX(i, i_node);
73 }
74
75 // read point field data
76 vtkPointData *p_field = d_grid_p->GetPointData();
77
78 std::vector<std::string> tags = {"Displacement", "Velocity"};
79 for (auto tag: tags) {
80
81 if (p_field->HasArray(tag.c_str())) {
82 vtkDataArray *array = p_field->GetArray(tag.c_str());
83
84 auto data_a = vtkSmartPointer<vtkDoubleArray>::New();
85 data_a->SetNumberOfComponents(3);
86 data_a->Allocate(3, 1); // allocate memory
87
88 for (size_t i = 0; i < array->GetNumberOfTuples(); i++) {
89 array->GetTuples(i, i, data_a);
90 auto xi = util::Point(data_a->GetValue(0), data_a->GetValue(1),
91 data_a->GetValue(2));
92 if (tag == "Displacement")
93 model->setU(i, xi);
94 else if (tag == "Velocity")
95 model->setV(i, xi);
96 }
97 } else {
98
99 if (tag == "Velocity") {
100 std::cerr << "Error: Restart file does not have <Velocity> data. "
101 "Aborting reading.\n";
102 exit(1);
103 } else if (tag == "Displacement") {
104 // compute from current and reference configuration
105 for (size_t i=0; i<model->d_x.size(); i++)
106 model->setU(i, model->getX(i) - model->getXRef(i));
107 }
108 }
109 }
110}
111
112
114 // delete d_reader_p;
115 // delete d_grid_p;
116}
A class to store model data.
Definition modelData.h:47
VtkParticleReader(const std::string &filename)
Constructor.
void close()
Closes the file and store it to the hard disk.
void readNodes(model::ModelData *model)
Writes the nodes to the file.
vtkSmartPointer< vtkXMLUnstructuredGridReader > d_reader_p
XML unstructured grid writer.
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
A structure to represent 3d vectors.
Definition point.h:30