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