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 
30 namespace 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 
46 int 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 
63 int 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  else if (cellType == VTK_QUADRATIC_WEDGE)
116  {
117  std::array<vtkIdType, 15> ogs_nodeIds{};
118  for (unsigned i = 0; i < 15; ++i)
119  {
120  ogs_nodeIds[i] = ptIds->GetId(i);
121  }
122  for (unsigned i = 0; i < 3; ++i)
123  {
124  ptIds->SetId(i, ogs_nodeIds[i + 3]);
125  ptIds->SetId(i + 3, ogs_nodeIds[i]);
126  }
127  for (unsigned i = 0; i < 3; ++i)
128  {
129  ptIds->SetId(6 + i, ogs_nodeIds[8 - i]);
130  }
131  for (unsigned i = 0; i < 3; ++i)
132  {
133  ptIds->SetId(9 + i, ogs_nodeIds[14 - i]);
134  }
135  ptIds->SetId(12, ogs_nodeIds[9]);
136  ptIds->SetId(13, ogs_nodeIds[11]);
137  ptIds->SetId(14, ogs_nodeIds[10]);
138  }
139 
140  output->InsertNextCell(cellType, ptIds);
141  }
142 
143  // Arrays
144  for (auto [name, property] : _mesh->getProperties())
145  {
146  if (auto p = dynamic_cast<PropertyVector<double>*>(property))
147  {
148  addProperty(*p);
149  }
150  else if (auto p = dynamic_cast<PropertyVector<float>*>(property))
151  {
152  addProperty(*p);
153  }
154  else if (auto p = dynamic_cast<PropertyVector<int>*>(property))
155  {
156  addProperty(*p);
157  }
158  else if (auto p = dynamic_cast<PropertyVector<unsigned>*>(property))
159  {
160  addProperty(*p);
161  }
162  else if (auto p = dynamic_cast<PropertyVector<long>*>(property))
163  {
164  addProperty(*p);
165  }
166  else if (auto p = dynamic_cast<PropertyVector<long long>*>(property))
167  {
168  addProperty(*p);
169  }
170  else if (auto p =
171  dynamic_cast<PropertyVector<unsigned long>*>(property))
172  {
173  addProperty(*p);
174  }
175  else if (auto p = dynamic_cast<PropertyVector<unsigned long long>*>(
176  property))
177  {
178  addProperty(*p);
179  }
180  else if (auto p = dynamic_cast<PropertyVector<std::size_t>*>(property))
181  {
182  addProperty(*p);
183  }
184  else if (auto p = dynamic_cast<PropertyVector<char>*>(property))
185  {
186  addProperty(*p);
187  }
188  else if (auto p =
189  dynamic_cast<PropertyVector<unsigned char>*>(property))
190  {
191  addProperty(*p);
192  }
193  else
194  {
195  OGS_FATAL(
196  "Mesh property '{:s}' of unhandled data type '{:s}'. Please "
197  "check the data type of the mesh properties. The available "
198  "data types are:"
199  "\n\t double,"
200  "\n\t float,"
201  "\n\t int,"
202  "\n\t unsigned,"
203  "\n\t long,"
204  "\n\t long long,"
205  "\n\t unsigned long long,"
206  "\n\t char,",
207  "\n\t unsigned char.",
208  property->getPropertyName(),
209  typeid(*property).name());
210  }
211  }
212 
213  output->GetPointData()->ShallowCopy(this->PointData.GetPointer());
214  output->GetCellData()->ShallowCopy(this->CellData.GetPointer());
215  output->GetFieldData()->ShallowCopy(this->FieldData.GetPointer());
216  return 1;
217 }
218 
220  vtkInformation* /*request*/,
221  vtkInformationVector** /*inputVector*/,
222  vtkInformationVector* /*outputVector*/)
223 {
224  this->NumberOfDimensions = 3;
226 
227  return 1;
228 }
229 
230 template <typename T>
232  MeshLib::PropertyVector<T> const& property) const
233 {
234  vtkNew<vtkAOSDataArrayTemplate<T>> dataArray;
235  const bool hasArrayOwnership = false;
236  dataArray->SetArray(const_cast<T*>(property.data()),
237  static_cast<vtkIdType>(property.size()),
238  static_cast<int>(!hasArrayOwnership));
239  dataArray->SetNumberOfComponents(property.getNumberOfGlobalComponents());
240  dataArray->SetName(property.getPropertyName().c_str());
241 
243  {
244  this->PointData->AddArray(dataArray.GetPointer());
245  }
246  else if (property.getMeshItemType() == MeshLib::MeshItemType::Cell)
247  {
248  this->CellData->AddArray(dataArray.GetPointer());
249  }
250  else if (property.getMeshItemType() ==
252  {
253  this->FieldData->AddArray(dataArray.GetPointer());
254  }
255 }
256 } // 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)
Definition: VtkOGSEnum.cpp:17
virtual unsigned getNumberOfNodes() const =0
std::vector< Node * > const & getNodes() const
Get the nodes-vector for the mesh.
Definition: Mesh.h:95
std::vector< Element * > const & getElements() const
Get the element-vector for the mesh.
Definition: Mesh.h:98
const std::string getName() const
Get name of the mesh.
Definition: Mesh.h:92
Properties & getProperties()
Definition: Mesh.h:123
std::size_t getNumberOfNodes() const
Get the number of nodes.
Definition: Mesh.h:89
std::string const & getPropertyName() const
MeshItemType getMeshItemType() const
int getNumberOfGlobalComponents() 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
int ProcessRequest(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
static VtkMappedMeshSource * New()
vtkNew< vtkFieldData > FieldData
int RequestInformation(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
vtkNew< vtkCellData > CellData
void PrintSelf(std::ostream &os, vtkIndent indent) override
vtkStandardNewMacro(VtkMappedMeshSource) void VtkMappedMeshSource