Skip to content
Snippets Groups Projects
Commit 201ddae2 authored by Nicolas Vuaille's avatar Nicolas Vuaille
Browse files

Introduce vtkFieldDataToDataSetAttribute

* this copy a FieldData array in another DatasetAttribute.
It supports vtkDataArray subclasses and copy only the first component
of the first tuple.
This makes use of ImplicitArrays to save memory.
parent 2668fa35
No related branches found
No related tags found
No related merge requests found
## Introduce vtkFieldDataToDataSetAttribute
With the new `vtkFieldDataToDataSetAttribute` filter VTK provides now
a way to efficiently pass FieldData single-value arrays to other AttributeData.
This is for instance useful with composite data, where FieldData
can be used to store a single scalar, varying at block level only.
This scalar was mostly not usable for computation in other filters.
Moving it, for instance, to PointData, allows to use it in your pipeline.
This is done at low memory cost by using the ImplicitArrays design.
......@@ -49,6 +49,7 @@ set(classes
vtkExtractEdges
vtkFeatureEdges
vtkFieldDataToAttributeDataFilter
vtkFieldDataToDataSetAttribute
vtkFlyingEdges2D
vtkFlyingEdges3D
vtkFlyingEdgesPlaneCutter
......
......@@ -45,6 +45,7 @@ vtk_add_test_cxx(vtkFiltersCoreCxxTests tests
TestExtractCells.cxx,NO_VALID
TestExtractCellsAlongPolyLine.cxx,NO_VALID
TestFeatureEdges.cxx,NO_VALID
TestFieldDataToDataSetAttribute.cxx,NO_VALID
TestFlyingEdges.cxx
TestGlyph3D.cxx
TestGlyph3DFollowCamera.cxx,NO_VALID
......
/*=========================================================================
Program: Visualization Toolkit
Module: TestTestFieldDataToAttributeData.cxx
Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
All rights reserved.
See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
This software is distributed WITHOUT ANY WARRANTY; without even
the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
PURPOSE. See the above copyright notice for more information.
=========================================================================*/
#include "vtkFieldDataToDataSetAttribute.h"
#include "vtkCellData.h"
#include "vtkDoubleArray.h"
#include "vtkFieldData.h"
#include "vtkImageData.h"
#include "vtkIntArray.h"
#include "vtkLogger.h"
#include "vtkMolecule.h"
#include "vtkMultiBlockDataSet.h"
#include "vtkNew.h"
#include "vtkPointData.h"
#include "vtkStringArray.h"
#include "vtkTable.h"
static const std::string FIRST_NAME = "firstArray";
static const std::string SECOND_NAME = "secondArray";
static const int FIRST_VALUE = 13;
static const double SECOND_VALUE = -3.7;
static const vtkIdType TESTED_INDEX = 42;
namespace FieldDataToAttributeDataUtils
{
void AddFieldDataArrays(vtkDataObject* obj, double shift = 0)
{
vtkNew<vtkIntArray> array1;
array1->SetName(FIRST_NAME.c_str());
vtkNew<vtkDoubleArray> array2;
array2->SetName(SECOND_NAME.c_str());
array1->InsertNextValue(FIRST_VALUE + shift);
array2->InsertNextValue(SECOND_VALUE + shift);
vtkFieldData* fieldData = obj->GetFieldData();
fieldData->AddArray(array1);
fieldData->AddArray(array2);
}
bool CheckOutput(vtkDataObject* output, vtkDataObject::AttributeTypes attributeType, int size,
const std::string& name, double value)
{
vtkDataSetAttributes* outAttribute = output->GetAttributes(attributeType);
if (!outAttribute)
{
vtkLog(ERROR,
"Cannot find attribute type " << vtkDataObject::GetAssociationTypeAsString(attributeType));
return false;
}
if (outAttribute->GetNumberOfArrays() != size)
{
vtkLog(ERROR,
"Wrong number of attribute arrays for type " << attributeType << ". Has "
<< outAttribute->GetNumberOfArrays());
return false;
}
vtkDataArray* outArray = outAttribute->GetArray(name.c_str());
if (!outArray)
{
vtkLog(ERROR, "Cannot find array in output with name '" << name << "'");
return false;
}
if (outArray->GetNumberOfTuples() != outAttribute->GetNumberOfTuples())
{
vtkLog(ERROR,
"Wrong array size : " << outArray->GetNumberOfTuples() << ". Expected "
<< outAttribute->GetNumberOfTuples());
return false;
}
if (outArray->GetTuple1(TESTED_INDEX) != value)
{
vtkLog(ERROR,
"Wrong value for array: has " << outArray->GetTuple1(TESTED_INDEX) << "instead of " << value);
return false;
}
return true;
}
bool TestDataObject(vtkDataObject* obj, vtkDataObject::AttributeTypes attributeType)
{
AddFieldDataArrays(obj);
vtkNew<vtkFieldDataToDataSetAttribute> forwarder;
forwarder->SetInputData(obj);
forwarder->SetOutputFieldType(attributeType);
forwarder->Update();
bool ret = true;
vtkDataObject* output = forwarder->GetOutput();
// Some data object can have default data arrays. For instance molecules have default AtomData
// array (atomic number) and BondData array (bond order).
int numberOfArrays = obj->GetAttributes(attributeType)->GetNumberOfArrays();
ret = ret && CheckOutput(output, attributeType, numberOfArrays + 2, FIRST_NAME, FIRST_VALUE);
ret = ret && CheckOutput(output, attributeType, numberOfArrays + 2, SECOND_NAME, SECOND_VALUE);
forwarder->ProcessAllArraysOff();
forwarder->AddFieldDataArray(SECOND_NAME.c_str());
forwarder->Update();
ret = ret && CheckOutput(output, attributeType, numberOfArrays + 1, SECOND_NAME, SECOND_VALUE);
if (!ret)
{
vtkLog(ERROR, "Test fails for " << vtkDataObject::GetAssociationTypeAsString(attributeType));
}
return ret;
}
void TestPointCellData()
{
vtkNew<vtkImageData> image;
// create more than TESTED_INDEX elements
image->SetDimensions(10, 10, 10);
TestDataObject(image, vtkDataObject::CELL);
TestDataObject(image, vtkDataObject::POINT);
}
void TestRowData()
{
vtkNew<vtkTable> table;
// create more than TESTED_INDEX elements
table->SetNumberOfRows(2 * TESTED_INDEX);
TestDataObject(table, vtkDataObject::ROW);
}
void TestVertexEdgeData()
{
vtkNew<vtkMolecule> molecule;
// create more than TESTED_INDEX elements
for (int idx = 0; idx < TESTED_INDEX * 2; idx++)
{
molecule->AppendAtom();
}
TestDataObject(molecule, vtkDataObject::VERTEX);
TestDataObject(molecule, vtkDataObject::EDGE);
}
void TestMultiBlock()
{
vtkNew<vtkImageData> image;
image->SetDimensions(10, 10, 10);
AddFieldDataArrays(image);
vtkNew<vtkImageData> image2;
image->SetDimensions(10, 10, 10);
const double shiftValue = 1;
AddFieldDataArrays(image2, shiftValue);
vtkNew<vtkMultiBlockDataSet> mblock;
mblock->SetBlock(0, image);
mblock->SetBlock(1, image2);
const auto attributeType = vtkDataObject::POINT;
vtkNew<vtkFieldDataToDataSetAttribute> forwarder;
forwarder->SetInputData(mblock);
forwarder->SetOutputFieldType(attributeType);
forwarder->Update();
vtkDataObject* output = forwarder->GetOutput();
auto outMB = vtkMultiBlockDataSet::SafeDownCast(output);
const int numberOfOutArrays = 2;
bool ret =
CheckOutput(outMB->GetBlock(0), attributeType, numberOfOutArrays, FIRST_NAME, FIRST_VALUE);
ret = ret &&
CheckOutput(outMB->GetBlock(0), attributeType, numberOfOutArrays, SECOND_NAME, SECOND_VALUE);
ret = ret &&
CheckOutput(
outMB->GetBlock(1), attributeType, numberOfOutArrays, FIRST_NAME, FIRST_VALUE + shiftValue);
ret = ret &&
CheckOutput(
outMB->GetBlock(1), attributeType, numberOfOutArrays, SECOND_NAME, SECOND_VALUE + shiftValue);
if (!ret)
{
vtkLog(ERROR, "Test fails for vtkMultiBlockDataSet");
}
}
};
int TestFieldDataToDataSetAttribute(int vtkNotUsed(argc), char* vtkNotUsed(argv)[])
{
FieldDataToAttributeDataUtils::TestPointCellData();
FieldDataToAttributeDataUtils::TestRowData();
FieldDataToAttributeDataUtils::TestVertexEdgeData();
FieldDataToAttributeDataUtils::TestMultiBlock();
return EXIT_SUCCESS;
}
......@@ -15,6 +15,7 @@ DEPENDS
VTK::CommonExecutionModel
VTK::CommonMisc
PRIVATE_DEPENDS
VTK::CommonImplicitArrays
VTK::CommonMath
VTK::CommonSystem
VTK::CommonTransforms
......
/*=========================================================================
Program: Visualization Toolkit
Module: vtkFieldDataToDataSetAttribute.cxx
Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
All rights reserved.
See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
This software is distributed WITHOUT ANY WARRANTY; without even
the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
PURPOSE. See the above copyright notice for more information.
=========================================================================*/
#include "vtkFieldDataToDataSetAttribute.h"
#include "vtkArrayDispatch.h"
#include "vtkArrayListTemplate.h"
#include "vtkConstantArray.h"
#include "vtkDataArray.h"
#include "vtkDataArrayAccessor.h"
#include "vtkDataArrayMeta.h"
#include "vtkFieldData.h"
#include "vtkInformation.h"
#include "vtkInformationVector.h"
#include "vtkObjectFactory.h"
#include "vtkStringArray.h"
#include <set>
vtkStandardNewMacro(vtkFieldDataToDataSetAttribute);
//------------------------------------------------------------------------------
int vtkFieldDataToDataSetAttribute::FillInputPortInformation(int port, vtkInformation* info)
{
if (port == 0)
{
// Skip composite data sets so that executives will handle this automatically.
info->Remove(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE());
info->Append(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkDataSet");
info->Append(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkGenericDataSet");
info->Append(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkGraph");
info->Append(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkHyperTreeGrid");
info->Append(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkTable");
}
return 1;
}
//------------------------------------------------------------------------------
void vtkFieldDataToDataSetAttribute::AddFieldDataArray(const char* name)
{
if (!name)
{
vtkErrorMacro("name cannot be null.");
return;
}
this->FieldDataArrays.insert(std::string(name));
this->Modified();
}
//------------------------------------------------------------------------------
void vtkFieldDataToDataSetAttribute::RemoveFieldDataArray(const char* name)
{
if (!name)
{
vtkErrorMacro("name cannot be null.");
return;
}
this->FieldDataArrays.erase(name);
this->Modified();
}
//------------------------------------------------------------------------------
void vtkFieldDataToDataSetAttribute::ClearFieldDataArrays()
{
if (!this->FieldDataArrays.empty())
{
this->Modified();
}
this->FieldDataArrays.clear();
}
//------------------------------------------------------------------------------
const std::set<std::string>& vtkFieldDataToDataSetAttribute::GetFieldDataArrays()
{
return this->FieldDataArrays;
}
//------------------------------------------------------------------------------
struct ArrayForwarder
{
vtkSmartPointer<vtkDataArray> Result;
vtkIdType NumberOfTuples;
template <typename ArrayType>
void operator()(ArrayType* array)
{
using SourceT = vtk::GetAPIType<ArrayType>;
auto output = vtkSmartPointer<vtkConstantArray<SourceT>>::New();
vtkDataArrayAccessor<ArrayType> access(array);
output->ConstructBackend(access.Get(0, 0));
output->SetNumberOfComponents(1);
output->SetNumberOfTuples(this->NumberOfTuples);
output->SetName(array->GetName());
this->Result = output;
}
};
//------------------------------------------------------------------------------
int vtkFieldDataToDataSetAttribute::RequestData(
vtkInformation*, vtkInformationVector** inputVector, vtkInformationVector* outputVector)
{
vtkInformation* info = outputVector->GetInformationObject(0);
vtkDataObject* output = info->Get(vtkDataObject::DATA_OBJECT());
vtkInformation* inInfo = inputVector[0]->GetInformationObject(0);
vtkDataObject* input = inInfo->Get(vtkDataObject::DATA_OBJECT());
output->ShallowCopy(input);
vtkFieldData* inputFieldData = input->GetFieldData();
vtkDataSetAttributes* outAttribute = output->GetAttributes(this->OutputFieldType);
vtkSmartPointer<vtkFieldData> originFieldData;
if (!this->ProcessAllArrays)
{
// create a vtkFieldData containing just the selected arrays.
originFieldData = vtkSmartPointer<vtkFieldData>::New();
for (const auto& name : this->FieldDataArrays)
{
vtkAbstractArray* arr = inputFieldData->GetAbstractArray(name.c_str());
if (arr == nullptr)
{
vtkWarningMacro("field data array not found: " << name);
continue;
}
originFieldData->AddArray(arr);
}
}
else
{
originFieldData = inputFieldData;
}
ArrayForwarder forwarder;
forwarder.NumberOfTuples = outAttribute->GetNumberOfTuples();
int nbOfArrays = originFieldData->GetNumberOfArrays();
for (int idx = 0; idx < nbOfArrays; idx++)
{
vtkDataArray* dataArray = originFieldData->GetArray(idx);
if (!dataArray)
{
continue;
}
if (vtkArrayDispatch::Dispatch::Execute(dataArray, forwarder))
{
outAttribute->AddArray(forwarder.Result);
}
else if (auto stringArray =
vtkStringArray::SafeDownCast(originFieldData->GetAbstractArray(idx)))
{
vtkWarningMacro("string arrays are not supported, skipping " << stringArray->GetName());
}
}
return 1;
}
//------------------------------------------------------------------------------
void vtkFieldDataToDataSetAttribute::PrintSelf(ostream& os, vtkIndent indent)
{
this->Superclass::PrintSelf(os, indent);
os << indent
<< "OutputFieldType: " << vtkDataObject::GetAssociationTypeAsString(this->OutputFieldType)
<< "\n";
os << indent << "ProcessAllArrays" << (this->ProcessAllArrays ? "On\n" : "Off\n");
os << indent << "FieldDataArrays: \n";
for (const auto& array : this->FieldDataArrays)
{
os << indent << array << "\n";
}
}
VTK_ABI_NAMESPACE_END
/*=========================================================================
Program: Visualization Toolkit
Module: vtkFieldDataToDataSetAttribute.h
Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
All rights reserved.
See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
This software is distributed WITHOUT ANY WARRANTY; without even
the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
PURPOSE. See the above copyright notice for more information.
=========================================================================*/
/**
* @class vtkFieldDataToDataSetAttribute
* @brief map field data to other attribute data
*
* vtkFieldDataToDataSetAttribute is a filter that copies field data arrays into
* another attribute data arrays.
*
* This is done at very low memory cost by using the Implicit Array infrastructure.
*
* NOTE: It copies only the first component of the first tuple into a vtkConstantArray.
* vtkStringArray are not supported.
*
* @sa
* vtkFieldData vtkCellData
*/
#ifndef vtkFieldDataToDataSetAttribute_h
#define vtkFieldDataToDataSetAttribute_h
#include "vtkFiltersCoreModule.h" // For export macro
#include "vtkPassInputTypeAlgorithm.h"
#include "vtkDataObject.h" // for vtkDataObject::AttributeType enum
#include <set>
#include <string>
VTK_ABI_NAMESPACE_BEGIN
class VTKFILTERSCORE_EXPORT vtkFieldDataToDataSetAttribute : public vtkPassInputTypeAlgorithm
{
public:
static vtkFieldDataToDataSetAttribute* New();
vtkTypeMacro(vtkFieldDataToDataSetAttribute, vtkPassInputTypeAlgorithm);
void PrintSelf(ostream& os, vtkIndent indent) override;
///@{
/**
* Activate whether to process all input arrays or only the selected ones.
* If false, only arrays selected by the user will be considered by this filter.
* The default is true.
*/
vtkSetMacro(ProcessAllArrays, bool);
vtkGetMacro(ProcessAllArrays, bool);
vtkBooleanMacro(ProcessAllArrays, bool);
///@}
///@{
/**
* Set/Get the output attribute type.
*/
vtkSetMacro(OutputFieldType, int);
vtkGetMacro(OutputFieldType, int);
///@}
/**
* Arrays to be processed.
*/
///@{
/**
* Adds an array to be processed.
* This only has an effect if ProcessAllArrays is off.
* If the name is already present, nothing happens.
*/
virtual void AddFieldDataArray(const char* name);
/**
* Removes an array to be processed.
* This only has an effect if ProcessAllArrays is off.
* If the name is not present, nothing happens.
*/
virtual void RemoveFieldDataArray(const char* name);
/**
* Removes all arrays to be processed from the list.
* This only has an effect if ProcessAllArrays is off.
*/
virtual void ClearFieldDataArrays();
/**
* Get the names of the arrays to process.
*/
virtual const std::set<std::string>& GetFieldDataArrays();
///@}
protected:
vtkFieldDataToDataSetAttribute() = default;
~vtkFieldDataToDataSetAttribute() override = default;
/**
* Reimplemented to remove Composite support. This filter
* relies on the Executive and handles composite block per block.
*/
int FillInputPortInformation(int port, vtkInformation* info) override;
/**
* Reimplemented create data arrays as required.
*/
int RequestData(vtkInformation* request, vtkInformationVector** inputVector,
vtkInformationVector* outputVector) override;
private:
vtkFieldDataToDataSetAttribute(const vtkFieldDataToDataSetAttribute&) = delete;
void operator=(const vtkFieldDataToDataSetAttribute&) = delete;
bool ProcessAllArrays = true;
/// see vtkDataObject::AttributeTypes
int OutputFieldType = vtkDataObject::POINT;
std::set<std::string> FieldDataArrays;
};
VTK_ABI_NAMESPACE_END
#endif
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment