6#include <vtkCallbackCommand.h>
7#include <vtkCellData.h>
8#include <vtkCellTypes.h>
10#include <vtkImageData.h>
11#include <vtkPointData.h>
12#include <vtkPolyData.h>
13#include <vtkUnstructuredGrid.h>
14#include <vtkVersionMacros.h>
16static vtkCellArray *GetCells(vtkDataSet *dataSet) {
17 switch(dataSet->GetDataObjectType()) {
18 case VTK_UNSTRUCTURED_GRID: {
19 auto dataSetAsUG =
static_cast<vtkUnstructuredGrid *
>(dataSet);
20 return dataSetAsUG->GetCells();
23 auto dataSetAsPD =
static_cast<vtkPolyData *
>(dataSet);
24 return dataSetAsPD->GetNumberOfPolys() > 0 ? dataSetAsPD->GetPolys()
25 : dataSetAsPD->GetNumberOfLines() > 0 ? dataSetAsPD->GetLines()
26 : dataSetAsPD->GetVerts();
32static int checkCellTypes(vtkPointSet *
object) {
36#if VTK_VERSION_NUMBER >= VTK_VERSION_CHECK(9, 2, 0)
37 if(object->GetDataObjectType() == VTK_UNSTRUCTURED_GRID) {
38 auto objectAsUG = vtkUnstructuredGrid::SafeDownCast(
object);
39 auto distinctCellTypes = objectAsUG->GetDistinctCellTypesArray();
40 nTypes = distinctCellTypes->GetNumberOfTuples();
45 object->GetCellTypes(cellTypes);
46 nTypes = cellTypes->GetNumberOfTypes();
59 const auto &cellType =
object->GetCellType(0);
60 if(cellType != VTK_VERTEX && cellType != VTK_LINE
61 && cellType != VTK_TRIANGLE && cellType != VTK_TETRA)
77 void Init(vtkDataSet *dataSet) {
80 if(dataSet->IsA(
"vtkPointSet"))
81 this->observee =
static_cast<vtkObject *
>(GetCells(dataSet));
83 this->observee =
static_cast<vtkObject *
>(dataSet);
85 this->observee->AddObserver(vtkCommand::DeleteEvent,
this, 1);
92 this->observee->RemoveObserver(
this);
96 if(instance->registry.empty()) {
100 auto it = instance->registry.find(this->key);
101 if(it != instance->registry.end()) {
102 instance->registry.erase(it);
104 instance->printMsg(
"# Registered Triangulations: "
105 + std::to_string(instance->registry.size()),
114 auto cells = GetCells(dataSet);
118 if(dataSet->IsA(
"vtkImageData")) {
119 auto image =
static_cast<vtkImageData *
>(dataSet);
120 image->GetExtent(this->
extent);
121 image->GetOrigin(this->
origin);
122 image->GetSpacing(this->
spacing);
127 onDelete->Init(dataSet);
131 auto cells = GetCells(dataSet);
135 if(dataSet->IsA(
"vtkImageData")) {
136 auto image =
static_cast<vtkImageData *
>(dataSet);
143 image->GetExtent(extent_);
144 image->GetOrigin(origin_);
145 image->GetSpacing(spacing_);
146 image->GetDimensions(dimensions_);
153 for(
int i = 0; i < 6; i++)
154 if(this->
extent[i] != extent_[i])
156 for(
int i = 0; i < 3; i++)
157 if(this->
origin[i] != origin_[i] || this->
spacing[i] != spacing_[i]
167ttkTriangulationFactory::ttkTriangulationFactory() {
172 ttkTriangulationFactory::CreateImplicitTriangulation(vtkImageData *image) {
174 this->
printMsg(
"Initializing Implicit Triangulation", 0, 0,
177 auto triangulation = std::make_unique<ttk::Triangulation>();
180 image->GetExtent(extent);
183 image->GetOrigin(origin);
186 image->GetSpacing(spacing);
189 if(!spacing[1] && !spacing[2])
195 image->GetDimensions(dimensions);
197 double firstPoint[3];
198 firstPoint[0] = origin[0] + extent[0] * spacing[0];
199 firstPoint[1] = origin[1] + extent[2] * spacing[1];
200 firstPoint[2] = origin[2] + extent[4] * spacing[2];
202 triangulation->setInputGrid(firstPoint[0], firstPoint[1], firstPoint[2],
203 spacing[0], spacing[1], spacing[2], dimensions[0],
204 dimensions[1], dimensions[2]);
206 this->
printMsg(
"Initializing Implicit Triangulation", 1,
210 return triangulation;
214 ttkTriangulationFactory::CreateExplicitTriangulation(vtkPointSet *pointSet) {
217 auto points = pointSet->GetPoints();
219 this->
printErr(
"DataSet has uninitialized `vtkPoints`.");
223 auto cells = GetCells(pointSet);
225 this->
printErr(
"DataSet has uninitialized `vtkCellArray`.");
229 auto triangulation = std::make_unique<ttk::Triangulation>();
230 int const hasIndexArray
234 this->
printMsg(
"Initializing Compact Triangulation", 0, 0,
237 this->
printMsg(
"Initializing Explicit Triangulation", 0, 0,
243 auto pointDataType = points->GetDataType();
244 if(pointDataType != VTK_FLOAT && pointDataType != VTK_DOUBLE) {
245 this->
printErr(
"Unable to initialize 'ttk::Triangulation' for point "
246 "precision other than 'float' or 'double'.");
252 vtkAbstractArray *indexArray = pointSet->GetPointData()->GetAbstractArray(
254 triangulation->setStellarInputPoints(
255 points->GetNumberOfPoints(), pointDataArray,
256 (
int *)indexArray->GetVoidPointer(0), pointDataType == VTK_DOUBLE);
258 triangulation->setInputPoints(points->GetNumberOfPoints(), pointDataArray,
259 pointDataType == VTK_DOUBLE);
264 int const cellTypeStatus = checkCellTypes(pointSet);
265 if(cellTypeStatus == -1) {
266 this->
printWrn(
"Inhomogeneous cell dimensions detected.");
268 "Consider using `ttkExtract` to extract cells of a given dimension.");
270 }
else if(cellTypeStatus == -2) {
271 this->
printWrn(
"Cells are not simplices.");
272 this->
printWrn(
"Consider using `vtkTetrahedralize` in pre-processing.");
277 int const nCells = cells->GetNumberOfCells();
279 if(!cells->IsStorage64Bit()) {
280 if(cells->CanConvertTo64BitStorage()) {
281 this->
printWrn(
"Converting the cell array to 64-bit storage");
282 bool const success = cells->ConvertTo64BitStorage();
285 "Error converting the provided cell array to 64-bit storage");
290 "Cannot convert the provided cell array to 64-bit storage");
294 auto connectivity =
static_cast<vtkIdType *
>(
296 auto offsets =
static_cast<vtkIdType *
>(
302 = triangulation->setStellarInputCells(nCells, connectivity, offsets);
304 status = triangulation->setInputCells(nCells, connectivity, offsets);
309 "Run the `vtkTetrahedralize` filter to resolve the issue.");
315 this->
printMsg(
"Initializing Compact Triangulation", 1,
319 this->
printMsg(
"Initializing Explicit Triangulation", 1,
324 return triangulation;
328 ttkTriangulationFactory::CreateTriangulation(vtkDataSet *dataSet) {
329 switch(dataSet->GetDataObjectType()) {
330 case VTK_UNSTRUCTURED_GRID:
331 case VTK_POLY_DATA: {
332 return this->CreateExplicitTriangulation(
333 static_cast<vtkPointSet *
>(dataSet));
335 case VTK_IMAGE_DATA: {
336 return this->CreateImplicitTriangulation((vtkImageData *)dataSet);
339 this->
printErr(
"Unable to triangulate `"
340 + std::string(dataSet->GetClassName()) +
"`");
348 int debugLevel,
float cacheRatio, vtkDataSet *
object) {
350 instance->setDebugLevel(debugLevel);
355 auto it = instance->registry.find(key);
356 if(it != instance->registry.end()) {
358 if(it->second.isValid(
object)) {
361 triangulation = it->second.triangulation.get();
365 instance->registry.erase(key);
369 if(!triangulation && object->IsA(
"vtkImageData")) {
370 instance->FindImplicitTriangulation(
371 triangulation,
static_cast<vtkImageData *
>(
object));
373 instance->printMsg(
"Retrieving Equivalent Implicit-Triangulation",
378 triangulation = instance->CreateTriangulation(
object).release();
380 instance->registry.emplace(std::piecewise_construct,
381 std::forward_as_tuple(key),
382 std::forward_as_tuple(
object, triangulation));
387 "# Registered Triangulations: " + std::to_string(instance->registry.size()),
391 triangulation->setDebugLevel(debugLevel);
392 triangulation->setCacheSize(cacheRatio);
395 return triangulation;
398int ttkTriangulationFactory::FindImplicitTriangulation(
401 for(
const auto &it : this->
registry) {
402 if(it.second.owner->IsA(
"vtkImageData")) {
403 if(it.second.isValid(image)) {
404 triangulation = it.second.triangulation.get();
414 switch(dataSet->GetDataObjectType()) {
415 case VTK_IMAGE_DATA: {
419 auto cells = GetCells(dataSet);
#define ttkNotUsed(x)
Mark function/method parameters that are not used in the function body at all.
static RegistryKey GetKey(vtkDataSet *dataSet)
static ttk::Triangulation * GetTriangulation(int debugLevel, float cacheRatio, vtkDataSet *object)
static ttkTriangulationFactory Instance
static void * GetVoidPointer(vtkDataArray *array, vtkIdType start=0)
int printWrn(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
void setDebugMsgPrefix(const std::string &prefix)
int printMsg(const std::string &msg, const debug::Priority &priority=debug::Priority::INFO, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cout) const
int printErr(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
Triangulation is a class that provides time and memory efficient traversal methods on triangulations ...
const char compactTriangulationIndex[]
bool isValid(vtkDataSet *dataSet) const
RegistryValue(vtkDataSet *dataSet, ttk::Triangulation *triangulation_)
RegistryTriangulation triangulation
void Execute(vtkObject *, unsigned long ttkNotUsed(eventId), void *ttkNotUsed(callData)) override
static ttkOnDeleteCommand * New()
void Init(vtkDataSet *dataSet)
std::unique_ptr< ttk::Triangulation > RegistryTriangulation
printMsg(debug::output::BOLD+" | | | | | . \\ | | (__| | / __/| |_| / __/| (_) |"+debug::output::ENDCOLOR, debug::Priority::PERFORMANCE, debug::LineMode::NEW, stream)