OGS
VtkAppendArrayFilter Class Reference

Detailed Description

Definition at line 26 of file VtkAppendArrayFilter.h.

#include <VtkAppendArrayFilter.h>

Inheritance diagram for VtkAppendArrayFilter:
[legend]
Collaboration diagram for VtkAppendArrayFilter:
[legend]

Public Member Functions

 vtkTypeMacro (VtkAppendArrayFilter, vtkUnstructuredGridAlgorithm)
 
void PrintSelf (ostream &os, vtkIndent indent) override
 Prints the mesh data to an output stream.
 
void SetUserProperty (QString name, QVariant value) override
 Sets user properties.
 
void SetArray (const std::string &array_name, const std::vector< double > &new_array)
 
- Public Member Functions inherited from VtkAlgorithmProperties
 VtkAlgorithmProperties (QObject *parent=nullptr)
 Constructor (sets default values)
 
 ~VtkAlgorithmProperties () override
 
vtkProperty * GetProperties () const
 Returns the vtk properties.
 
vtkTexture * GetTexture ()
 Returns a texture (if one has been assigned).
 
void SetTexture (vtkTexture *t)
 Sets a texture for the VtkVisPipelineItem.
 
vtkLookupTable * GetLookupTable (const QString &array_name)
 Returns the colour lookup table (if one has been assigned).
 
void RemoveLookupTable (const QString &array_name)
 Removes the lookup table for the given scalar.
 
void SetLookUpTable (const QString &array_name, vtkLookupTable *lut)
 Sets a colour lookup table for the given scalar array of the VtkVisPipelineItem.
 
void SetLookUpTable (const QString &array_name, const QString &filename)
 Loads a predefined color lookup table from a file for the specified scalar array.
 
bool GetScalarVisibility () const
 Returns the scalar visibility.
 
void SetScalarVisibility (bool on)
 Sets the scalar visibility.
 
QString GetName () const
 Returns the name. This is set to the file path if it is a source algorithm.
 
void SetName (QString name)
 Sets the name.
 
bool IsRemovable () const
 Is this algorithm removable from the pipeline (view).
 
QMap< QString, QVariant > * GetAlgorithmUserProperties () const
 Returns a map of user properties.
 
QMap< QString, QList< QVariant > > * GetAlgorithmUserVectorProperties () const
 Returns a map of vector user properties.
 
QVariant GetUserProperty (QString name) const
 Returns the value of a user property.
 
virtual void SetUserVectorProperty (QString name, QList< QVariant > values)
 Sets a vector user property. This should be implemented by subclasses.
 
QList< QVariant > GetUserVectorProperty (QString name) const
 Returns a list of values of a vector user property.
 
void SetActiveAttribute (QString name)
 Set the active attribute.
 
QString GetActiveAttribute () const
 Returns the desired active attribute.
 

Static Public Member Functions

static VtkAppendArrayFilterNew ()
 Create new objects with New() because of VTKs object reference counting.
 

Protected Member Functions

 VtkAppendArrayFilter ()
 
 ~VtkAppendArrayFilter () override
 
int RequestData (vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
 The filter logic.
 

Private Attributes

std::vector< double > _array
 
std::string _array_name
 

Additional Inherited Members

- Signals inherited from VtkAlgorithmProperties
void ScalarVisibilityChanged (bool on)
 
- Protected Attributes inherited from VtkAlgorithmProperties
vtkProperty * _property
 
vtkTexture * _texture
 
bool _scalarVisibility
 
std::map< QString, vtkLookupTable * > _lut
 
QString _name
 
QString _activeAttributeName
 
bool _removable
 
QMap< QString, QVariant > * _algorithmUserProperties
 
QMap< QString, QList< QVariant > > * _algorithmUserVectorProperties
 

Constructor & Destructor Documentation

◆ VtkAppendArrayFilter()

VtkAppendArrayFilter::VtkAppendArrayFilter ( )
protecteddefault

◆ ~VtkAppendArrayFilter()

VtkAppendArrayFilter::~VtkAppendArrayFilter ( )
overrideprotecteddefault

Member Function Documentation

◆ New()

static VtkAppendArrayFilter * VtkAppendArrayFilter::New ( )
static

Create new objects with New() because of VTKs object reference counting.

◆ PrintSelf()

void VtkAppendArrayFilter::PrintSelf ( ostream & os,
vtkIndent indent )
override

Prints the mesh data to an output stream.

Definition at line 36 of file VtkAppendArrayFilter.cpp.

37{
38 this->Superclass::PrintSelf(os, indent);
39 os << indent << "== VtkAppendArrayFilter ==" << endl;
40}

◆ RequestData()

int VtkAppendArrayFilter::RequestData ( vtkInformation * request,
vtkInformationVector ** inputVector,
vtkInformationVector * outputVector )
overrideprotected

The filter logic.

Definition at line 42 of file VtkAppendArrayFilter.cpp.

45{
46 if (this->_array.empty())
47 {
48 ERR("VtkAppendArrayFilter::RequestData(): Selection array is empty.");
49 return 0;
50 }
51 vtkInformation* inInfo = inputVector[0]->GetInformationObject(0);
52 vtkUnstructuredGrid* input = vtkUnstructuredGrid::SafeDownCast(
53 inInfo->Get(vtkDataObject::DATA_OBJECT()));
54
55 vtkSmartPointer<vtkDoubleArray> colors =
56 vtkSmartPointer<vtkDoubleArray>::New();
57 colors->SetNumberOfComponents(1);
58 std::size_t arrayLength = this->_array.size();
59 colors->SetNumberOfValues(arrayLength);
60 colors->SetName("Selection");
61
62 std::size_t nCells = input->GetNumberOfCells();
63 if (nCells > arrayLength)
64 WARN(
65 "VtkAppendArrayFilter::RequestData(): Number of cells exceeds "
66 "selection array length. Surplus cells won't be examined.");
67
68 for (std::size_t i = 0; i < arrayLength; i++)
69 {
70 colors->SetValue(i, _array[i]);
71 }
72
73 vtkInformation* outInfo = outputVector->GetInformationObject(0);
74 vtkUnstructuredGrid* output = vtkUnstructuredGrid::SafeDownCast(
75 outInfo->Get(vtkDataObject::DATA_OBJECT()));
76 output->CopyStructure(input);
77 output->GetPointData()->PassData(input->GetPointData());
78 output->GetCellData()->PassData(input->GetCellData());
79 output->GetCellData()->AddArray(colors);
80 output->GetCellData()->SetActiveScalars(_array_name.c_str());
81
82 return 1;
83}
void ERR(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:45
void WARN(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:40
std::vector< double > _array

References _array, _array_name, ERR(), and WARN().

◆ SetArray()

void VtkAppendArrayFilter::SetArray ( const std::string & array_name,
const std::vector< double > & new_array )

Definition at line 85 of file VtkAppendArrayFilter.cpp.

87{
88 this->_array_name = array_name;
89 this->_array = new_array;
90}

References _array, and _array_name.

◆ SetUserProperty()

void VtkAppendArrayFilter::SetUserProperty ( QString name,
QVariant value )
inlineoverridevirtual

Sets user properties.

Reimplemented from VtkAlgorithmProperties.

Definition at line 38 of file VtkAppendArrayFilter.h.

39 {
40 Q_UNUSED(name);
41 Q_UNUSED(value);
42 }

◆ vtkTypeMacro()

VtkAppendArrayFilter::vtkTypeMacro ( VtkAppendArrayFilter ,
vtkUnstructuredGridAlgorithm  )

Member Data Documentation

◆ _array

std::vector<double> VtkAppendArrayFilter::_array
private

Definition at line 57 of file VtkAppendArrayFilter.h.

Referenced by RequestData(), and SetArray().

◆ _array_name

std::string VtkAppendArrayFilter::_array_name
private

Definition at line 58 of file VtkAppendArrayFilter.h.

Referenced by RequestData(), and SetArray().


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