OGS
MeshLib::VtkMappedMeshSource Class Referencefinal

Detailed Description

Adapter which maps a MeshLib::Mesh to a vtkUnstructuredGridAlgorithm. Allows for zero-copy access of the mesh from the visualization side.

Definition at line 44 of file VtkMappedMeshSource.h.

#include <VtkMappedMeshSource.h>

Inheritance diagram for MeshLib::VtkMappedMeshSource:
[legend]
Collaboration diagram for MeshLib::VtkMappedMeshSource:
[legend]

Public Member Functions

 vtkTypeMacro (VtkMappedMeshSource, vtkUnstructuredGridAlgorithm)
 
void PrintSelf (std::ostream &os, vtkIndent indent) override
 
 VtkMappedMeshSource (const VtkMappedMeshSource &)=delete
 
void operator= (const VtkMappedMeshSource &)=delete
 
void SetMesh (const MeshLib::Mesh *mesh)
 Sets the mesh. Calling is mandatory.
 
const MeshLib::MeshGetMesh () const
 Returns the mesh.
 

Static Public Member Functions

static VtkMappedMeshSourceNew ()
 

Protected Member Functions

 VtkMappedMeshSource ()
 
int ProcessRequest (vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
 
int RequestData (vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
 
int RequestInformation (vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
 

Private Member Functions

template<typename T >
void addProperty (MeshLib::PropertyVector< T > const &property) const
 Adds a zero-copy vtk array wrapper.
 

Private Attributes

const MeshLib::Mesh_mesh {}
 
int NumberOfDimensions {0}
 
int NumberOfNodes {0}
 
vtkNew< vtkPoints > Points
 
vtkNew< vtkPointData > PointData
 
vtkNew< vtkCellData > CellData
 
vtkNew< vtkFieldData > FieldData
 

Constructor & Destructor Documentation

◆ VtkMappedMeshSource() [1/2]

MeshLib::VtkMappedMeshSource::VtkMappedMeshSource ( const VtkMappedMeshSource & )
delete

◆ VtkMappedMeshSource() [2/2]

MeshLib::VtkMappedMeshSource::VtkMappedMeshSource ( )
protected

Definition at line 40 of file VtkMappedMeshSource.cpp.

42{
43 this->SetNumberOfInputPorts(0);
44}

Member Function Documentation

◆ addProperty()

template<typename T >
void MeshLib::VtkMappedMeshSource::addProperty ( MeshLib::PropertyVector< T > const & property) const
private

Adds a zero-copy vtk array wrapper.

Definition at line 221 of file VtkMappedMeshSource.cpp.

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}
MeshItemType getMeshItemType() const
int getNumberOfGlobalComponents() const
std::string const & getPropertyName() const
std::size_t size() const
vtkNew< vtkPointData > PointData
vtkNew< vtkFieldData > FieldData

References MeshLib::Cell, CellData, FieldData, MeshLib::PropertyVectorBase::getMeshItemType(), MeshLib::PropertyVectorBase::getNumberOfGlobalComponents(), MeshLib::PropertyVectorBase::getPropertyName(), MeshLib::IntegrationPoint, MeshLib::Node, PointData, and MeshLib::PropertyVector< PROP_VAL_TYPE >::size().

Referenced by RequestData().

◆ GetMesh()

const MeshLib::Mesh * MeshLib::VtkMappedMeshSource::GetMesh ( ) const
inline

Returns the mesh.

Definition at line 62 of file VtkMappedMeshSource.h.

62{ return _mesh; }

References _mesh.

Referenced by VtkVisPipeline::addPipelineItem(), MeshItem::getMesh(), and MainWindow::showMeshQualitySelectionDialog().

◆ New()

static VtkMappedMeshSource * MeshLib::VtkMappedMeshSource::New ( )
static

Referenced by MeshItem::MeshItem(), and RequestData().

◆ operator=()

void MeshLib::VtkMappedMeshSource::operator= ( const VtkMappedMeshSource & )
delete

◆ PrintSelf()

void MeshLib::VtkMappedMeshSource::PrintSelf ( std::ostream & os,
vtkIndent indent )
override

◆ ProcessRequest()

int MeshLib::VtkMappedMeshSource::ProcessRequest ( vtkInformation * request,
vtkInformationVector ** inputVector,
vtkInformationVector * outputVector )
overrideprotected

Definition at line 46 of file VtkMappedMeshSource.cpp.

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}
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
int RequestInformation(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override

References RequestData(), and RequestInformation().

◆ RequestData()

int MeshLib::VtkMappedMeshSource::RequestData ( vtkInformation * ,
vtkInformationVector ** ,
vtkInformationVector * outputVector )
overrideprotected

Definition at line 63 of file VtkMappedMeshSource.cpp.

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 =
158 dynamic_cast<PropertyVector<unsigned long long>*>(
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}
#define OGS_FATAL(...)
Definition Error.h:26
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
void addProperty(MeshLib::PropertyVector< T > const &property) const
Adds a zero-copy vtk array wrapper.
static VtkMappedMeshSource * New()

References _mesh, addProperty(), CellData, FieldData, MeshLib::Mesh::getElements(), MeshLib::Mesh::getNodes(), MeshLib::Element::getNumberOfNodes(), MeshLib::Mesh::getProperties(), New(), OGS_FATAL, OGSToVtkCellType(), PointData, and Points.

Referenced by ProcessRequest().

◆ RequestInformation()

int MeshLib::VtkMappedMeshSource::RequestInformation ( vtkInformation * ,
vtkInformationVector ** ,
vtkInformationVector *  )
overrideprotected

Definition at line 209 of file VtkMappedMeshSource.cpp.

213{
214 this->NumberOfDimensions = 3;
216
217 return 1;
218}
std::size_t getNumberOfNodes() const
Get the number of nodes.
Definition Mesh.h:100

References _mesh, MeshLib::Mesh::getNumberOfNodes(), NumberOfDimensions, and NumberOfNodes.

Referenced by ProcessRequest().

◆ SetMesh()

void MeshLib::VtkMappedMeshSource::SetMesh ( const MeshLib::Mesh * mesh)
inline

Sets the mesh. Calling is mandatory.

Definition at line 55 of file VtkMappedMeshSource.h.

56 {
57 this->_mesh = mesh;
58 this->Modified();
59 }

References _mesh.

◆ vtkTypeMacro()

MeshLib::VtkMappedMeshSource::vtkTypeMacro ( VtkMappedMeshSource ,
vtkUnstructuredGridAlgorithm  )

Member Data Documentation

◆ _mesh

const MeshLib::Mesh* MeshLib::VtkMappedMeshSource::_mesh {}
private

◆ CellData

vtkNew<vtkCellData> MeshLib::VtkMappedMeshSource::CellData
private

Definition at line 89 of file VtkMappedMeshSource.h.

Referenced by addProperty(), and RequestData().

◆ FieldData

vtkNew<vtkFieldData> MeshLib::VtkMappedMeshSource::FieldData
private

Definition at line 90 of file VtkMappedMeshSource.h.

Referenced by addProperty(), and RequestData().

◆ NumberOfDimensions

int MeshLib::VtkMappedMeshSource::NumberOfDimensions {0}
private

Definition at line 84 of file VtkMappedMeshSource.h.

84{0};

Referenced by RequestInformation().

◆ NumberOfNodes

int MeshLib::VtkMappedMeshSource::NumberOfNodes {0}
private

Definition at line 85 of file VtkMappedMeshSource.h.

85{0};

Referenced by RequestInformation().

◆ PointData

vtkNew<vtkPointData> MeshLib::VtkMappedMeshSource::PointData
private

Definition at line 88 of file VtkMappedMeshSource.h.

Referenced by addProperty(), and RequestData().

◆ Points

vtkNew<vtkPoints> MeshLib::VtkMappedMeshSource::Points
private

Definition at line 87 of file VtkMappedMeshSource.h.

Referenced by RequestData().


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