Skip to content
Snippets Groups Projects
Commit 2379018f authored by Julien Fausty's avatar Julien Fausty Committed by Kitware Robot
Browse files

Merge topic 'add-pass-through-cells-HTGGeometry'


5e1968a7 HTGGeometry: add PassThroughCellIds changelog
7782cbe8 HTGGeometry: add testing for PassThroughCellIds option
192daf2b HTGGeometry: add PassThroughCellIds option for hardware picking

Acked-by: default avatarKitware Robot <kwrobot@kitware.com>
Acked-by: default avatarbuildbot <buildbot@kitware.com>
Reviewed-by: default avatarMathieu Westphal <mathieu.westphal@kitware.com>
Reviewed-by: default avatarThomas Galland <thomas.galland@kitware.com>
Merge-request: !10131
parents a197a461 5e1968a7
No related branches found
No related tags found
No related merge requests found
## Addition of `PassThroughCellIds` property in `vtkHyperTreeGridGeometry` filter
The `vtkHyperTreeGridGeometry` filter now provides and option `PassThroughCellIds` (default `false`) to pass through original cell IDs from the input `vtkHyperTreeGrid` to the output `vtkPolyData`. One can also name the array using the `OriginalCellIdArrayName` string proerty.
This option, often found in geometry filters, allows one to make the resulting `vtkPolyData` compatible with the `vtkHardwareSelector` for making selections on the original `vtkHyperTreeGrid`. It can also help for debugging purposes as well as general additional visualization data when mapping the surface to color.
......@@ -38,6 +38,7 @@ set(test_sources
TestHyperTreeGridBinaryClipPlanes.cxx
TestHyperTreeGridBinaryEllipseMaterial.cxx
TestHyperTreeGridBinaryHyperbolicParaboloidMaterial.cxx
TestHyperTreeGridGeometryPassCellIds.cxx
TestHyperTreeGridTernary2D.cxx
TestHyperTreeGridTernary2DBiMaterial.cxx
TestHyperTreeGridTernary2DFullMaterialBits.cxx
......
/*=========================================================================
Program: Visualization Toolkit
Module: TestHyperTreeGridGeometryPassCellIds.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.
=========================================================================*/
// Funded by CEA, DAM, DIF, F-91297 Arpajon, France
#include "vtkHyperTreeGrid.h"
#include "vtkHyperTreeGridGeometry.h"
#include "vtkHyperTreeGridPreConfiguredSource.h"
#include "vtkActor.h"
#include "vtkCamera.h"
#include "vtkDataSet.h"
#include "vtkLookupTable.h"
#include "vtkNew.h"
#include "vtkPolyDataMapper.h"
#include "vtkProperty.h"
#include "vtkRegressionTestImage.h"
#include "vtkRenderWindow.h"
#include "vtkRenderWindowInteractor.h"
#include "vtkRenderer.h"
int TestHyperTreeGridGeometryPassCellIds(int argc, char* argv[])
{
vtkNew<vtkHyperTreeGridPreConfiguredSource> htgSource;
htgSource->SetHTGMode(vtkHyperTreeGridPreConfiguredSource::CUSTOM);
htgSource->SetCustomArchitecture(vtkHyperTreeGridPreConfiguredSource::UNBALANCED);
htgSource->SetCustomDim(3);
htgSource->SetCustomFactor(3);
htgSource->SetCustomDepth(4);
std::vector<unsigned int> subdivs = { 3, 3, 3 };
std::vector<double> extent = { -1, 1, -1, 1, -1, 1 };
htgSource->SetCustomSubdivisions(subdivs.data());
htgSource->SetCustomExtent(extent.data());
vtkNew<vtkHyperTreeGridGeometry> geom;
geom->SetPassThroughCellIds(true);
geom->SetInputConnection(htgSource->GetOutputPort());
geom->Update();
vtkNew<vtkPolyDataMapper> mapper;
mapper->SetInputConnection(geom->GetOutputPort());
auto ds = vtkDataSet::SafeDownCast(geom->GetOutput());
vtkNew<vtkLookupTable> lut;
lut->SetNumberOfTableValues(100);
lut->SetTableRange(0, ds->GetNumberOfCells());
mapper->ScalarVisibilityOn();
mapper->SetLookupTable(lut);
mapper->UseLookupTableScalarRangeOn();
mapper->SetScalarModeToUseCellFieldData();
mapper->ColorByArrayComponent("vtkOriginalCellIds", 0);
mapper->InterpolateScalarsBeforeMappingOn();
vtkNew<vtkActor> actor;
actor->SetMapper(mapper);
actor->GetProperty()->SetRepresentationToSurface();
actor->GetProperty()->EdgeVisibilityOn();
vtkNew<vtkRenderer> renderer;
renderer->AddActor(actor);
vtkNew<vtkRenderWindow> renWin;
renWin->AddRenderer(renderer.Get());
vtkCamera* camera = renderer->GetActiveCamera();
camera->SetPosition(-1.5, -1.5, -1.5);
renderer->ResetCamera();
renWin->Render();
return !vtkRegressionTester::Test(argc, argv, renWin, 10);
}
6f4aadcc67f2a4af91c697d3650c346e34ec973df964f7410b6126047a4792fcca910ebfbf80fcb583426c13fc83bda46608521116e9a73a498a2a2e4ebd1663
......@@ -43,6 +43,23 @@
#include <set>
#include <vector>
namespace
{
bool PassCellId(vtkDataArray* cellIds, vtkIdType inId, vtkIdType outId)
{
auto typedCellIds = vtkIdTypeArray::SafeDownCast(cellIds);
if (!typedCellIds)
{
vtkErrorWithObjectMacro(nullptr, "Pass cell ids array has wrong type.");
return false;
}
typedCellIds->InsertValue(outId, inId);
return true;
}
}
VTK_ABI_NAMESPACE_BEGIN
static constexpr unsigned int VonNeumannCursors3D[] = { 0, 1, 2, 4, 5, 6 };
static constexpr unsigned int VonNeumannOrientations3D[] = { 2, 1, 0, 0, 1, 2 };
......@@ -292,6 +309,14 @@ int vtkHyperTreeGridGeometry::ProcessTrees(vtkHyperTreeGrid* input, vtkDataObjec
this->OutData = output->GetCellData();
this->OutData->CopyAllocate(this->InData);
if (this->PassThroughCellIds)
{
vtkNew<vtkIdTypeArray> originalCellIds;
originalCellIds->SetName(this->OriginalCellIdArrayName.c_str());
originalCellIds->SetNumberOfComponents(1);
this->OutData->AddArray(originalCellIds);
}
// Retrieve material mask
this->Mask = input->HasMask() ? input->GetMask() : nullptr;
......@@ -495,6 +520,10 @@ void vtkHyperTreeGridGeometry::ProcessLeaf1D(vtkHyperTreeGridNonOrientedGeometry
// Copy edge data from that of the cell from which it comes
this->OutData->CopyData(this->InData, inId, outId);
if (this->PassThroughCellIds)
{
::PassCellId(this->OutData->GetArray(this->OriginalCellIdArrayName.c_str()), inId, outId);
}
}
//------------------------------------------------------------------------------
......@@ -797,6 +826,10 @@ void vtkHyperTreeGridGeometry::ProcessLeaf3D(
// Copy face data from that of the cell from which it comes
this->OutData->CopyData(this->InData, inId, outId);
if (this->PassThroughCellIds)
{
::PassCellId(this->OutData->GetArray(this->OriginalCellIdArrayName.c_str()), inId, outId);
}
} // if ( nA > 0 )
// Create face B when its vertices are present
......@@ -840,6 +873,10 @@ void vtkHyperTreeGridGeometry::ProcessLeaf3D(
// Copy face data from that of the cell from which it comes
this->OutData->CopyData(this->InData, inId, outId);
if (this->PassThroughCellIds)
{
::PassCellId(this->OutData->GetArray(this->OriginalCellIdArrayName.c_str()), inId, outId);
}
} // if ( nB > 0 )
} // if ( this->HasInterface )
}
......@@ -942,6 +979,10 @@ void vtkHyperTreeGridGeometry::AddFace(vtkIdType useId, const double* origin, co
// Copy face data from that of the cell from which it comes
this->OutData->CopyData(this->InData, useId, outId);
if (this->PassThroughCellIds)
{
::PassCellId(this->OutData->GetArray(this->OriginalCellIdArrayName.c_str()), useId, outId);
}
}
//------------------------------------------------------------------------------
void vtkHyperTreeGridGeometry::AddFace2(vtkIdType inId, vtkIdType useId, const double* origin,
......@@ -1273,6 +1314,10 @@ void vtkHyperTreeGridGeometry::AddFace2(vtkIdType inId, vtkIdType useId, const d
// Copy face data from that of the cell from which it comes
this->OutData->CopyData(this->InData, useId, outId);
if (this->PassThroughCellIds)
{
::PassCellId(this->OutData->GetArray(this->OriginalCellIdArrayName.c_str()), inId, outId);
}
} // if ( create )
}
VTK_ABI_NAMESPACE_END
......@@ -64,6 +64,32 @@ public:
vtkGetMacro(Merging, bool);
///@}
//@{
/**
* Set/Get for the PassThroughCellIds boolean.
*
* When set to true this boolean ensures an array named with whatever is
* in `OriginalCellIdArrayName` gets created in the output holding the
* original cell ids
*
* default is false
*/
vtkSetMacro(PassThroughCellIds, bool);
vtkGetMacro(PassThroughCellIds, bool);
vtkBooleanMacro(PassThroughCellIds, bool);
/**
* Set/Get the OriginalCellIdArrayName string.
*
* When PassThroughCellIds is set to true, the name of the generated
* array is whatever is set in this variable.
*
* default to vtkOriginalCellIds
*/
vtkSetMacro(OriginalCellIdArrayName, std::string);
vtkGetMacro(OriginalCellIdArrayName, std::string);
//@}
protected:
vtkHyperTreeGridGeometry();
~vtkHyperTreeGridGeometry() override;
......@@ -143,6 +169,18 @@ protected:
*/
vtkCellArray* Cells;
/**
* Boolean for passing cell ids to poly data
*
* default is false
*/
bool PassThroughCellIds = false;
/**
* Name of the array holding original cell ids in output if PassThroughCellIds is true
*/
std::string OriginalCellIdArrayName = "vtkOriginalCellIds";
/**
*JB Un locator est utilise afin de produire un maillage avec moins
*JB de points. Le gain en 3D est de l'ordre d'un facteur 4 !
......
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