Loading [MathJax]/extensions/tex2jax.js
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 42 of file VtkMappedMeshSource.cpp.

44{
45 this->SetNumberOfInputPorts(0);
46}

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 226 of file VtkMappedMeshSource.cpp.

228{
229 vtkNew<vtkAOSDataArrayTemplate<T>> dataArray;
230 const bool hasArrayOwnership = false;
231 dataArray->SetArray(const_cast<T*>(property.data()),
232 static_cast<vtkIdType>(property.size()),
233 static_cast<int>(!hasArrayOwnership));
234 dataArray->SetNumberOfComponents(property.getNumberOfGlobalComponents());
235 dataArray->SetName(property.getPropertyName().c_str());
236
238 {
239 this->PointData->AddArray(dataArray.GetPointer());
240 }
241 else if (property.getMeshItemType() == MeshLib::MeshItemType::Cell)
242 {
243 this->CellData->AddArray(dataArray.GetPointer());
244 }
245 else if (property.getMeshItemType() ==
247 {
248 this->FieldData->AddArray(dataArray.GetPointer());
249 }
250}
MeshItemType getMeshItemType() const
int getNumberOfGlobalComponents() const
std::string const & getPropertyName() const
constexpr std::size_t size() const
constexpr const PROP_VAL_TYPE * data() const
vtkNew< vtkPointData > PointData
vtkNew< vtkFieldData > FieldData

References MeshLib::Cell, CellData, MeshLib::PropertyVector< T >::data(), FieldData, MeshLib::PropertyVectorBase::getMeshItemType(), MeshLib::PropertyVectorBase::getNumberOfGlobalComponents(), MeshLib::PropertyVectorBase::getPropertyName(), MeshLib::IntegrationPoint, MeshLib::Node, PointData, and MeshLib::PropertyVector< T >::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().

◆ 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 48 of file VtkMappedMeshSource.cpp.

51{
52 if (request->Has(vtkDemandDrivenPipeline::REQUEST_DATA()))
53 {
54 return this->RequestData(request, inputVector, outputVector);
55 }
56
57 if (request->Has(vtkDemandDrivenPipeline::REQUEST_INFORMATION()))
58 {
59 return this->RequestInformation(request, inputVector, outputVector);
60 }
61
62 return this->Superclass::ProcessRequest(request, inputVector, outputVector);
63}
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 65 of file VtkMappedMeshSource.cpp.

68{
69 vtkSmartPointer<vtkInformation> outInfo =
70 outputVector->GetInformationObject(0);
71 vtkSmartPointer<vtkUnstructuredGrid> output =
72 vtkUnstructuredGrid::SafeDownCast(
73 outInfo->Get(vtkDataObject::DATA_OBJECT()));
74
75 if (outInfo->Get(vtkStreamingDemandDrivenPipeline::UPDATE_PIECE_NUMBER()) >
76 0)
77 {
78 return 1;
79 }
80
81 // Points
82 this->Points->Reset();
83
84 MeshNodalCoordinatesBackend backend(_mesh->getNodes());
85
86 vtkNew<vtkImplicitArray<MeshNodalCoordinatesBackend>> nodeCoords;
87 nodeCoords->ConstructBackend(backend);
88 nodeCoords->SetNumberOfComponents(3);
89 nodeCoords->SetNumberOfTuples(_mesh->getNumberOfNodes());
90 this->Points->SetData(nodeCoords);
91 output->SetPoints(this->Points.GetPointer());
92
93 // Cells
94 auto elems = _mesh->getElements();
95 output->Allocate(elems.size());
96 for (auto& cell : elems)
97 {
98 auto cellType = OGSToVtkCellType(cell->getCellType());
99
100 const MeshLib::Element* const elem = cell;
101 const unsigned numNodes(elem->getNumberOfNodes());
102 const MeshLib::Node* const* nodes = cell->getNodes();
103 vtkSmartPointer<vtkIdList> ptIds = vtkSmartPointer<vtkIdList>::New();
104 ptIds->SetNumberOfIds(numNodes);
105
106 for (unsigned i = 0; i < numNodes; ++i)
107 {
108 ptIds->SetId(i, nodes[i]->getID());
109 }
110
111 if (cellType == VTK_WEDGE)
112 {
113 for (unsigned i = 0; i < 3; ++i)
114 {
115 const auto prism_swap_id = ptIds->GetId(i);
116 ptIds->SetId(i, ptIds->GetId(i + 3));
117 ptIds->SetId(i + 3, prism_swap_id);
118 }
119 }
120
121 output->InsertNextCell(cellType, ptIds);
122 }
123
124 // Arrays
125 for (auto [name, property] : _mesh->getProperties())
126 {
127 if (!property->is_for_output)
128 {
129 continue;
130 }
131 if (auto const* p = dynamic_cast<PropertyVector<double>*>(property))
132 {
133 addProperty(*p);
134 }
135 else if (auto const* p = dynamic_cast<PropertyVector<float>*>(property))
136 {
137 addProperty(*p);
138 }
139 else if (auto const* p = dynamic_cast<PropertyVector<int>*>(property))
140 {
141 addProperty(*p);
142 }
143 else if (auto const* p =
144 dynamic_cast<PropertyVector<unsigned>*>(property))
145 {
146 addProperty(*p);
147 }
148 else if (auto const* p = dynamic_cast<PropertyVector<long>*>(property))
149 {
150 addProperty(*p);
151 }
152 else if (auto const* p =
153 dynamic_cast<PropertyVector<long long>*>(property))
154 {
155 addProperty(*p);
156 }
157 else if (auto const* p =
158 dynamic_cast<PropertyVector<unsigned long>*>(property))
159 {
160 addProperty(*p);
161 }
162 else if (auto const* p =
163 dynamic_cast<PropertyVector<unsigned long long>*>(
164 property))
165 {
166 addProperty(*p);
167 }
168 else if (auto const* p =
169 dynamic_cast<PropertyVector<std::size_t>*>(property))
170 {
171 addProperty(*p);
172 }
173 else if (auto const* p = dynamic_cast<PropertyVector<char>*>(property))
174 {
175 addProperty(*p);
176 }
177 else if (auto const* p =
178 dynamic_cast<PropertyVector<unsigned char>*>(property))
179 {
180 addProperty(*p);
181 }
182 else if (auto const* p =
183 dynamic_cast<PropertyVector<uint8_t>*>(property))
184 {
185 addProperty(*p);
186 }
187 else
188 {
189 OGS_FATAL(
190 "Mesh property '{:s}' of unhandled data type '{:s}'. Please "
191 "check the data type of the mesh properties. The available "
192 "data types are:"
193 "\n\t double,"
194 "\n\t float,"
195 "\n\t int,"
196 "\n\t unsigned,"
197 "\n\t long,"
198 "\n\t long long,"
199 "\n\t unsigned long long,"
200 "\n\t char,",
201 "\n\t unsigned char,",
202 "\n\t uint8_t.",
203 property->getPropertyName(),
204 typeid(*property).name());
205 }
206 }
207
208 output->GetPointData()->ShallowCopy(this->PointData.GetPointer());
209 output->GetCellData()->ShallowCopy(this->CellData.GetPointer());
210 output->GetFieldData()->ShallowCopy(this->FieldData.GetPointer());
211 return 1;
212}
#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:108
std::vector< Element * > const & getElements() const
Get the element-vector for the mesh.
Definition Mesh.h:111
std::size_t getNumberOfNodes() const
Get the number of nodes.
Definition Mesh.h:102
void addProperty(MeshLib::PropertyVector< T > const &property) const
Adds a zero-copy vtk array wrapper.

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

Referenced by ProcessRequest().

◆ RequestInformation()

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

◆ 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: