mirror of
https://gitlab.kitware.com/vtk/vtk-m
synced 2024-09-08 13:23:51 +00:00
Output CoordinateSystemTransforms results in Coordinates
This commit changes how CoordinateSystemTransforms write their result. Before theoy would write their result in a DataSet in which the new Coords where stored in a field with the name of: - cylindricalCoordinateSystemTransform - sphericalCoordinateSystemTransform Now, they write their results as a DataSet in which its first Coords are the transformed Coords. Previous Coordinates are appended Signed-off-by: Vicente Adolfo Bolea Sanchez <vicente.bolea@kitware.com>
This commit is contained in:
parent
fd9c21c0d4
commit
e4aa20594a
@ -28,7 +28,7 @@ public:
|
||||
using SupportedTypes = vtkm::TypeListFieldVec3;
|
||||
|
||||
VTKM_CONT
|
||||
CylindricalCoordinateTransform();
|
||||
CylindricalCoordinateTransform() = default;
|
||||
|
||||
VTKM_CONT void SetCartesianToCylindrical() { Worklet.SetCartesianToCylindrical(); }
|
||||
VTKM_CONT void SetCylindricalToCartesian() { Worklet.SetCylindricalToCartesian(); }
|
||||
@ -49,17 +49,16 @@ public:
|
||||
using SupportedTypes = vtkm::TypeListFieldVec3;
|
||||
|
||||
VTKM_CONT
|
||||
SphericalCoordinateTransform();
|
||||
SphericalCoordinateTransform() = default;
|
||||
|
||||
VTKM_CONT void SetCartesianToSpherical() { Worklet.SetCartesianToSpherical(); }
|
||||
VTKM_CONT void SetSphericalToCartesian() { Worklet.SetSphericalToCartesian(); }
|
||||
|
||||
template <typename T, typename StorageType, typename DerivedPolicy>
|
||||
VTKM_CONT vtkm::cont::DataSet DoExecute(
|
||||
const vtkm::cont::DataSet& input,
|
||||
const vtkm::cont::ArrayHandle<T, StorageType>& field,
|
||||
const vtkm::filter::FieldMetadata& fieldMeta,
|
||||
const vtkm::filter::PolicyBase<DerivedPolicy>& policy) const;
|
||||
VTKM_CONT vtkm::cont::DataSet DoExecute(const vtkm::cont::DataSet&,
|
||||
const vtkm::cont::ArrayHandle<T, StorageType>&,
|
||||
const vtkm::filter::FieldMetadata&,
|
||||
const vtkm::filter::PolicyBase<DerivedPolicy>& policy);
|
||||
|
||||
private:
|
||||
vtkm::worklet::SphericalCoordinateTransform Worklet;
|
||||
|
@ -16,31 +16,31 @@ namespace vtkm
|
||||
namespace filter
|
||||
{
|
||||
|
||||
//-----------------------------------------------------------------------------
|
||||
inline VTKM_CONT CylindricalCoordinateTransform::CylindricalCoordinateTransform()
|
||||
: Worklet()
|
||||
{
|
||||
this->SetOutputFieldName("cylindricalCoordinateSystemTransform");
|
||||
}
|
||||
|
||||
//-----------------------------------------------------------------------------
|
||||
template <typename T, typename StorageType, typename DerivedPolicy>
|
||||
inline VTKM_CONT vtkm::cont::DataSet CylindricalCoordinateTransform::DoExecute(
|
||||
const vtkm::cont::DataSet& inDataSet,
|
||||
const vtkm::cont::ArrayHandle<T, StorageType>& field,
|
||||
const vtkm::filter::FieldMetadata& fieldMetadata,
|
||||
const vtkm::filter::FieldMetadata& vtkmNotUsed(fieldMetadata),
|
||||
const vtkm::filter::PolicyBase<DerivedPolicy>&)
|
||||
{
|
||||
vtkm::cont::ArrayHandle<T> outArray;
|
||||
this->Worklet.Run(field, outArray);
|
||||
return CreateResult(inDataSet, outArray, this->GetOutputFieldName(), fieldMetadata);
|
||||
}
|
||||
vtkm::cont::DataSet outDataSet;
|
||||
|
||||
//-----------------------------------------------------------------------------
|
||||
inline VTKM_CONT SphericalCoordinateTransform::SphericalCoordinateTransform()
|
||||
: Worklet()
|
||||
{
|
||||
this->SetOutputFieldName("sphericalCoordinateSystemTransform");
|
||||
this->Worklet.Run(field, outArray);
|
||||
|
||||
// We first add the result coords to keep them at the first position
|
||||
// of the resulting dataset.
|
||||
outDataSet.AddCoordinateSystem(vtkm::cont::CoordinateSystem("coordinates", outArray));
|
||||
|
||||
for (int i = 0; i < inDataSet.GetNumberOfCoordinateSystems(); i++)
|
||||
{
|
||||
outDataSet.AddCoordinateSystem(inDataSet.GetCoordinateSystem(i));
|
||||
}
|
||||
|
||||
outDataSet.SetCellSet(inDataSet.GetCellSet());
|
||||
|
||||
return outDataSet;
|
||||
}
|
||||
|
||||
//-----------------------------------------------------------------------------
|
||||
@ -48,12 +48,26 @@ template <typename T, typename StorageType, typename DerivedPolicy>
|
||||
inline VTKM_CONT vtkm::cont::DataSet SphericalCoordinateTransform::DoExecute(
|
||||
const vtkm::cont::DataSet& inDataSet,
|
||||
const vtkm::cont::ArrayHandle<T, StorageType>& field,
|
||||
const vtkm::filter::FieldMetadata& fieldMetadata,
|
||||
const vtkm::filter::PolicyBase<DerivedPolicy>&) const
|
||||
const vtkm::filter::FieldMetadata& vtkmNotUsed(fieldMetadata),
|
||||
const vtkm::filter::PolicyBase<DerivedPolicy>&)
|
||||
{
|
||||
vtkm::cont::ArrayHandle<T> outArray;
|
||||
Worklet.Run(field, outArray);
|
||||
return CreateResult(inDataSet, outArray, this->GetOutputFieldName(), fieldMetadata);
|
||||
vtkm::cont::DataSet outDataSet;
|
||||
|
||||
this->Worklet.Run(field, outArray);
|
||||
|
||||
// We first add the result coords to keep them at the first position
|
||||
// of the resulting dataset.
|
||||
outDataSet.AddCoordinateSystem(vtkm::cont::CoordinateSystem("coordinates", outArray));
|
||||
|
||||
for (int i = 0; i < inDataSet.GetNumberOfCoordinateSystems(); i++)
|
||||
{
|
||||
outDataSet.AddCoordinateSystem(inDataSet.GetCoordinateSystem(i));
|
||||
}
|
||||
|
||||
outDataSet.SetCellSet(inDataSet.GetCellSet());
|
||||
|
||||
return outDataSet;
|
||||
}
|
||||
}
|
||||
} // namespace vtkm::filter
|
||||
|
@ -66,10 +66,10 @@ vtkm::cont::DataSet MakeTestDataSet(const CoordinateType& cType)
|
||||
vtkm::FloatDefault R = 1.0f;
|
||||
vtkm::FloatDefault eps = vtkm::Epsilon<float>();
|
||||
std::vector<vtkm::FloatDefault> Thetas = {
|
||||
eps, vtkm::Pif() / 4, vtkm::Pif() / 3, vtkm::Pif() / 2, vtkm::Pif() - eps
|
||||
eps, vtkm::Pif() / 4.0f, vtkm::Pif() / 3.0f, vtkm::Pif() / 2.0f, vtkm::Pif() - eps
|
||||
};
|
||||
std::vector<vtkm::FloatDefault> Phis = {
|
||||
eps, vtkm::TwoPif() / 4, vtkm::TwoPif() / 3, vtkm::TwoPif() / 2, vtkm::TwoPif() - eps
|
||||
eps, vtkm::TwoPif() / 4.0f, vtkm::TwoPif() / 3.0f, vtkm::TwoPif() / 2.0f, vtkm::TwoPif() - eps
|
||||
};
|
||||
for (std::size_t i = 0; i < Thetas.size(); i++)
|
||||
for (std::size_t j = 0; j < Phis.size(); j++)
|
||||
@ -136,14 +136,12 @@ void TestCoordinateSystemTransform()
|
||||
vtkm::cont::DataSet dsCart = MakeTestDataSet(CART);
|
||||
vtkm::filter::CylindricalCoordinateTransform cylTrn;
|
||||
|
||||
cylTrn.SetOutputFieldName("cylindricalCoords");
|
||||
cylTrn.SetUseCoordinateSystemAsField(true);
|
||||
cylTrn.SetCartesianToCylindrical();
|
||||
cylTrn.SetUseCoordinateSystemAsField(true);
|
||||
vtkm::cont::DataSet carToCylDataSet = cylTrn.Execute(dsCart);
|
||||
|
||||
cylTrn.SetCylindricalToCartesian();
|
||||
cylTrn.SetUseCoordinateSystemAsField(true);
|
||||
cylTrn.SetOutputFieldName("cartesianCoords");
|
||||
vtkm::cont::DataSet cylToCarDataSet = cylTrn.Execute(carToCylDataSet);
|
||||
ValidateCoordTransform(dsCart, cylToCarDataSet, { false, false, false });
|
||||
|
||||
@ -151,24 +149,20 @@ void TestCoordinateSystemTransform()
|
||||
vtkm::cont::DataSet dsCyl = MakeTestDataSet(CYL);
|
||||
cylTrn.SetCylindricalToCartesian();
|
||||
cylTrn.SetUseCoordinateSystemAsField(true);
|
||||
cylTrn.SetOutputFieldName("cartesianCoords");
|
||||
cylToCarDataSet = cylTrn.Execute(dsCyl);
|
||||
|
||||
cylTrn.SetCartesianToCylindrical();
|
||||
cylTrn.SetUseCoordinateSystemAsField(true);
|
||||
cylTrn.SetOutputFieldName("cylindricalCoords");
|
||||
carToCylDataSet = cylTrn.Execute(cylToCarDataSet);
|
||||
ValidateCoordTransform(dsCyl, carToCylDataSet, { false, true, false });
|
||||
|
||||
std::cout << "Testing SphericalCoordinateTransform Filter" << std::endl;
|
||||
|
||||
vtkm::filter::SphericalCoordinateTransform sphTrn;
|
||||
sphTrn.SetOutputFieldName("sphericalCoords");
|
||||
sphTrn.SetUseCoordinateSystemAsField(true);
|
||||
sphTrn.SetCartesianToSpherical();
|
||||
vtkm::cont::DataSet carToSphDataSet = sphTrn.Execute(dsCart);
|
||||
|
||||
sphTrn.SetOutputFieldName("cartesianCoords");
|
||||
sphTrn.SetUseCoordinateSystemAsField(true);
|
||||
sphTrn.SetSphericalToCartesian();
|
||||
vtkm::cont::DataSet sphToCarDataSet = sphTrn.Execute(carToSphDataSet);
|
||||
@ -177,13 +171,11 @@ void TestCoordinateSystemTransform()
|
||||
vtkm::cont::DataSet dsSph = MakeTestDataSet(SPH);
|
||||
sphTrn.SetSphericalToCartesian();
|
||||
sphTrn.SetUseCoordinateSystemAsField(true);
|
||||
sphTrn.SetOutputFieldName("sphericalCoords");
|
||||
sphToCarDataSet = sphTrn.Execute(dsSph);
|
||||
|
||||
sphTrn.SetCartesianToSpherical();
|
||||
sphTrn.SetUseCoordinateSystemAsField(true);
|
||||
sphTrn.SetOutputFieldName("sphericalCoords");
|
||||
carToSphDataSet = cylTrn.Execute(sphToCarDataSet);
|
||||
carToSphDataSet = sphTrn.Execute(sphToCarDataSet);
|
||||
ValidateCoordTransform(dsSph, carToSphDataSet, { false, true, true });
|
||||
}
|
||||
|
||||
|
Loading…
Reference in New Issue
Block a user