OGS
VtkMappedMeshSource.cpp
Go to the documentation of this file.
1
14#include "VtkMappedMeshSource.h"
15
16#include <vtkCellType.h>
17#include <vtkDemandDrivenPipeline.h>
18#include <vtkInformation.h>
19#include <vtkInformationVector.h>
20#include <vtkSmartPointer.h>
21#include <vtkStreamingDemandDrivenPipeline.h>
22
23#include <vector>
24
26#include "MeshLib/Mesh.h"
27#include "MeshLib/VtkOGSEnum.h"
29
30namespace MeshLib
31{
33
34 void VtkMappedMeshSource::PrintSelf(std::ostream& os, vtkIndent indent)
35{
36 this->Superclass::PrintSelf(os, indent);
37 os << indent << "Mesh: " << (_mesh ? _mesh->getName() : "(none)") << endl;
38}
39
41
42{
43 this->SetNumberOfInputPorts(0);
44}
45
46int VtkMappedMeshSource::ProcessRequest(vtkInformation* request,
47 vtkInformationVector** inputVector,
48 vtkInformationVector* outputVector)
49{
50 if (request->Has(vtkDemandDrivenPipeline::REQUEST_DATA()))
51 {
52 return this->RequestData(request, inputVector, outputVector);
53 }
54
55 if (request->Has(vtkDemandDrivenPipeline::REQUEST_INFORMATION()))
56 {
57 return this->RequestInformation(request, inputVector, outputVector);
58 }
59
60 return this->Superclass::ProcessRequest(request, inputVector, outputVector);
61}
62
63int VtkMappedMeshSource::RequestData(vtkInformation* /*request*/,
64 vtkInformationVector** /*inputVector*/,
65 vtkInformationVector* outputVector)
66{
67 vtkSmartPointer<vtkInformation> outInfo =
68 outputVector->GetInformationObject(0);
69 vtkSmartPointer<vtkUnstructuredGrid> output =
70 vtkUnstructuredGrid::SafeDownCast(
71 outInfo->Get(vtkDataObject::DATA_OBJECT()));
72
73 if (outInfo->Get(vtkStreamingDemandDrivenPipeline::UPDATE_PIECE_NUMBER()) >
74 0)
75 {
76 return 1;
77 }
78
79 // Points
80 this->Points->Reset();
81
82 vtkSmartPointer<VtkMeshNodalCoordinatesTemplate<double>> nodeCoords =
83 vtkSmartPointer<VtkMeshNodalCoordinatesTemplate<double>>::New();
84 nodeCoords->SetNodes(_mesh->getNodes());
85 this->Points->SetData(nodeCoords);
86 output->SetPoints(this->Points.GetPointer());
87
88 // Cells
89 auto elems = _mesh->getElements();
90 output->Allocate(elems.size());
91 for (auto& cell : elems)
92 {
93 auto cellType = OGSToVtkCellType(cell->getCellType());
94
95 const MeshLib::Element* const elem = cell;
96 const unsigned numNodes(elem->getNumberOfNodes());
97 const MeshLib::Node* const* nodes = cell->getNodes();
98 vtkSmartPointer<vtkIdList> ptIds = vtkSmartPointer<vtkIdList>::New();
99 ptIds->SetNumberOfIds(numNodes);
100
101 for (unsigned i = 0; i < numNodes; ++i)
102 {
103 ptIds->SetId(i, nodes[i]->getID());
104 }
105
106 if (cellType == VTK_WEDGE)
107 {
108 for (unsigned i = 0; i < 3; ++i)
109 {
110 const auto prism_swap_id = ptIds->GetId(i);
111 ptIds->SetId(i, ptIds->GetId(i + 3));
112 ptIds->SetId(i + 3, prism_swap_id);
113 }
114 }
115
116 output->InsertNextCell(cellType, ptIds);
117 }
118
119 // Arrays
120 for (auto [name, property] : _mesh->getProperties())
121 {
122 if (!property->is_for_output)
123 {
124 continue;
125 }
126 if (auto const* p = dynamic_cast<PropertyVector<double>*>(property))
127 {
128 addProperty(*p);
129 }
130 else if (auto const* p = dynamic_cast<PropertyVector<float>*>(property))
131 {
132 addProperty(*p);
133 }
134 else if (auto const* p = dynamic_cast<PropertyVector<int>*>(property))
135 {
136 addProperty(*p);
137 }
138 else if (auto const* p =
139 dynamic_cast<PropertyVector<unsigned>*>(property))
140 {
141 addProperty(*p);
142 }
143 else if (auto const* p = dynamic_cast<PropertyVector<long>*>(property))
144 {
145 addProperty(*p);
146 }
147 else if (auto const* p =
148 dynamic_cast<PropertyVector<long long>*>(property))
149 {
150 addProperty(*p);
151 }
152 else if (auto const* p =
153 dynamic_cast<PropertyVector<unsigned long>*>(property))
154 {
155 addProperty(*p);
156 }
157 else if (auto const* p =
159 property))
160 {
161 addProperty(*p);
162 }
163 else if (auto const* p =
164 dynamic_cast<PropertyVector<std::size_t>*>(property))
165 {
166 addProperty(*p);
167 }
168 else if (auto const* p = dynamic_cast<PropertyVector<char>*>(property))
169 {
170 addProperty(*p);
171 }
172 else if (auto const* p =
173 dynamic_cast<PropertyVector<unsigned char>*>(property))
174 {
175 addProperty(*p);
176 }
177 else if (auto const* p =
178 dynamic_cast<PropertyVector<uint8_t>*>(property))
179 {
180 addProperty(*p);
181 }
182 else
183 {
184 OGS_FATAL(
185 "Mesh property '{:s}' of unhandled data type '{:s}'. Please "
186 "check the data type of the mesh properties. The available "
187 "data types are:"
188 "\n\t double,"
189 "\n\t float,"
190 "\n\t int,"
191 "\n\t unsigned,"
192 "\n\t long,"
193 "\n\t long long,"
194 "\n\t unsigned long long,"
195 "\n\t char,",
196 "\n\t unsigned char,",
197 "\n\t uint8_t.",
198 property->getPropertyName(),
199 typeid(*property).name());
200 }
201 }
202
203 output->GetPointData()->ShallowCopy(this->PointData.GetPointer());
204 output->GetCellData()->ShallowCopy(this->CellData.GetPointer());
205 output->GetFieldData()->ShallowCopy(this->FieldData.GetPointer());
206 return 1;
207}
208
210 vtkInformation* /*request*/,
211 vtkInformationVector** /*inputVector*/,
212 vtkInformationVector* /*outputVector*/)
213{
214 this->NumberOfDimensions = 3;
216
217 return 1;
218}
219
220template <typename T>
222 MeshLib::PropertyVector<T> const& property) const
223{
224 vtkNew<vtkAOSDataArrayTemplate<T>> dataArray;
225 const bool hasArrayOwnership = false;
226 dataArray->SetArray(const_cast<T*>(property.data()),
227 static_cast<vtkIdType>(property.size()),
228 static_cast<int>(!hasArrayOwnership));
229 dataArray->SetNumberOfComponents(property.getNumberOfGlobalComponents());
230 dataArray->SetName(property.getPropertyName().c_str());
231
233 {
234 this->PointData->AddArray(dataArray.GetPointer());
235 }
236 else if (property.getMeshItemType() == MeshLib::MeshItemType::Cell)
237 {
238 this->CellData->AddArray(dataArray.GetPointer());
239 }
240 else if (property.getMeshItemType() ==
242 {
243 this->FieldData->AddArray(dataArray.GetPointer());
244 }
245}
246} // Namespace MeshLib
Definition of the Element class.
#define OGS_FATAL(...)
Definition Error.h:26
Definition of the Mesh class.
VtkMappedMeshSource is a source class to transform OGS meshes into complete vtkUnstructuredGrids....
VtkMeshNodalCoordinatesTemplate is a adapter for node coordinates of OGS meshes to VTK unstructured g...
int OGSToVtkCellType(MeshLib::CellType ogs)
virtual unsigned getNumberOfNodes() const =0
std::vector< Node * > const & getNodes() const
Get the nodes-vector for the mesh.
Definition Mesh.h:106
std::vector< Element * > const & getElements() const
Get the element-vector for the mesh.
Definition Mesh.h:109
Properties & getProperties()
Definition Mesh.h:134
const std::string getName() const
Get name of the mesh.
Definition Mesh.h:103
std::size_t getNumberOfNodes() const
Get the number of nodes.
Definition Mesh.h:100
MeshItemType getMeshItemType() const
int getNumberOfGlobalComponents() const
std::string const & getPropertyName() const
std::size_t size() const
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
void addProperty(MeshLib::PropertyVector< T > const &property) const
Adds a zero-copy vtk array wrapper.
vtkNew< vtkPointData > PointData
static VtkMappedMeshSource * New()
int ProcessRequest(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
vtkNew< vtkFieldData > FieldData
int RequestInformation(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
void PrintSelf(std::ostream &os, vtkIndent indent) override
vtkStandardNewMacro(VtkMappedMeshSource) void VtkMappedMeshSource