VTK  9.1.0
vtkMultiThreshold.h
Go to the documentation of this file.
1/*=========================================================================
2
3 Program: Visualization Toolkit
4 Module: vtkMultiThreshold.h
5
6 Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
7 All rights reserved.
8 See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
9
10 This software is distributed WITHOUT ANY WARRANTY; without even
11 the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
12 PURPOSE. See the above copyright notice for more information.
13
14=========================================================================*/
15/*----------------------------------------------------------------------------
16 Copyright (c) Sandia Corporation
17 See Copyright.txt or http://www.paraview.org/HTML/Copyright.html for details.
18----------------------------------------------------------------------------*/
19
112#ifndef vtkMultiThreshold_h
113#define vtkMultiThreshold_h
114
115#include "vtkFiltersGeneralModule.h" // For export macro
116#include "vtkMath.h" // for Inf() and NegInf()
118
119#include <map> // for IntervalRules map
120#include <set> // for UpdateDependents()
121#include <string> // for holding array names in NormKey
122#include <vector> // for lists of threshold rules
123
124class vtkCell;
125class vtkCellData;
126class vtkDataArray;
127class vtkGenericCell;
128class vtkPointSet;
130
131class VTKFILTERSGENERAL_EXPORT vtkMultiThreshold : public vtkMultiBlockDataSetAlgorithm
132{
133public:
136 void PrintSelf(ostream& os, vtkIndent indent) override;
137
140 {
141 OPEN = 0,
142 CLOSED = 1
143 };
145 enum Norm
146 {
147 LINFINITY_NORM = -3,
148 L2_NORM = -2,
149 L1_NORM = -1
150 };
154 {
159 NAND
160 };
161
163
215 int AddIntervalSet(double xmin, double xmax, int omin, int omax, int assoc, const char* arrayName,
216 int component, int allScalars);
217 int AddIntervalSet(double xmin, double xmax, int omin, int omax, int assoc, int attribType,
218 int component, int allScalars);
220
222
229 int AddLowpassIntervalSet(
230 double xmax, int assoc, const char* arrayName, int component, int allScalars);
231 int AddHighpassIntervalSet(
232 double xmin, int assoc, const char* arrayName, int component, int allScalars);
233 int AddBandpassIntervalSet(
234 double xmin, double xmax, int assoc, const char* arrayName, int component, int allScalars);
235 int AddNotchIntervalSet(
236 double xlo, double xhi, int assoc, const char* arrayName, int component, int allScalars);
238
242 int AddBooleanSet(int operation, int numInputs, int* inputs);
243
247 int OutputSet(int setId);
248
252 void Reset();
253
256 typedef double (*TupleNorm)(vtkDataArray* arr, vtkIdType tuple, int component);
257
258 // NormKey must be able to use TupleNorm typedef:
259 class NormKey;
260
261 // Interval must be able to use NormKey typedef:
262 class Interval;
263
264 // Set needs to refer to boolean set pointers
265 class BooleanSet;
266
269 {
270 public:
271 int Association; // vtkDataObject::FIELD_ASSOCIATION_POINTS or
272 // vtkDataObject::FIELD_ASSOCIATION_CELLS
273 int Type; // -1 => use Name, otherwise: vtkDataSetAttributes::{SCALARS, VECTORS, TENSORS,
274 // NORMALS, TCOORDS, GLOBALIDS}
275 std::string Name; // Either empty or (when ArrayType == -1) an input array name
276 int Component; // LINFINITY_NORM, L1_NORM, L2_NORM or an integer component number
277 int AllScalars; // For Association == vtkDataObject::FIELD_ASSOCIATION_POINTS, must all points
278 // be in the interval?
279 int InputArrayIndex; // The number passed to vtkAlgorithm::SetInputArrayToProcess()
280 TupleNorm NormFunction; // A function pointer to compute the norm (or fetcht the correct
281 // component) of a tuple.
282
286 vtkIdType cellId, vtkCell* cell, vtkDataArray* array, double cellNorm[2]) const;
287
290 bool operator<(const NormKey& other) const
291 {
292 if (this->Association < other.Association)
293 return true;
294 else if (this->Association > other.Association)
295 return false;
296
297 if (this->Component < other.Component)
298 return true;
299 else if (this->Component > other.Component)
300 return false;
301
302 if ((!this->AllScalars) && other.AllScalars)
303 return true;
304 else if (this->AllScalars && (!other.AllScalars))
305 return false;
306
307 if (this->Type == -1)
308 {
309 if (other.Type == -1)
310 return this->Name < other.Name;
311 return true;
312 }
313 else
314 {
315 return this->Type < other.Type;
316 }
317 }
318 };
319
324 class Set
325 {
326 public:
327 int Id;
330
333 Set() { this->OutputId = -1; }
335 virtual ~Set() = default;
337 virtual void PrintNodeName(ostream& os);
339 virtual void PrintNode(ostream& os) = 0;
341 virtual BooleanSet* GetBooleanSetPointer();
342 virtual Interval* GetIntervalPointer();
343 };
344
346 class Interval : public Set
347 {
348 public:
350 double EndpointValues[2];
352 int EndpointClosures[2];
355
361 int Match(double cellNorm[2]);
362
363 ~Interval() override = default;
364 void PrintNode(ostream& os) override;
365 Interval* GetIntervalPointer() override;
366 };
367
369 class BooleanSet : public Set
370 {
371 public:
375 std::vector<int> Inputs;
376
378 BooleanSet(int sId, int op, int* inBegin, int* inEnd)
379 : Inputs(inBegin, inEnd)
380 {
381 this->Id = sId;
382 this->Operator = op;
383 }
384 ~BooleanSet() override = default;
385 void PrintNode(ostream& os) override;
386 BooleanSet* GetBooleanSetPointer() override;
387 };
388
389protected:
392
408 {
409 INCONCLUSIVE = -1,
410 INCLUDE = -2,
411 EXCLUDE = -3
412 };
413
418
425
432
437
439 typedef std::vector<Interval*> IntervalList;
441 typedef std::map<NormKey, IntervalList> RuleMap;
442
443 typedef std::vector<int> TruthTreeValues;
444 typedef std::vector<TruthTreeValues> TruthTree;
445
450
456 std::vector<Set*> Sets;
457
465
469 void UpdateDependents(int id, std::set<int>& unresolvedOutputs, TruthTreeValues& setStates,
470 vtkCellData* inCellData, vtkIdType inCell, vtkGenericCell* cell,
471 std::vector<vtkUnstructuredGrid*>& outv);
472
476 int AddIntervalSet(NormKey& nk, double xmin, double xmax, int omin, int omax);
477
481 void PrintGraph(ostream& os);
482
484 void operator=(const vtkMultiThreshold&) = delete;
485};
486
488 double xmax, int assoc, const char* arrayName, int component, int allScalars)
489{
490 return this->AddIntervalSet(
491 vtkMath::NegInf(), xmax, CLOSED, CLOSED, assoc, arrayName, component, allScalars);
492}
493
495 double xmin, int assoc, const char* arrayName, int component, int allScalars)
496{
497 return this->AddIntervalSet(
498 xmin, vtkMath::Inf(), CLOSED, CLOSED, assoc, arrayName, component, allScalars);
499}
500
502 double xmin, double xmax, int assoc, const char* arrayName, int component, int allScalars)
503{
504 return this->AddIntervalSet(xmin, xmax, CLOSED, CLOSED, assoc, arrayName, component, allScalars);
505}
506
508 double xlo, double xhi, int assoc, const char* arrayName, int component, int allScalars)
509{
510 int band =
511 this->AddIntervalSet(xlo, xhi, CLOSED, CLOSED, assoc, arrayName, component, allScalars);
512 if (band < 0)
513 {
514 return -1;
515 }
516 return this->AddBooleanSet(NAND, 1, &band);
517}
518
520{
521 return nullptr;
522}
523
525{
526 return nullptr;
527}
528
530{
531 return this;
532}
533
535{
536 return this;
537}
538
539#endif // vtkMultiThreshold_h
represent and manipulate cell attribute data
Definition: vtkCellData.h:42
abstract class to specify cell behavior
Definition: vtkCell.h:67
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:59
provides thread-safe access to cells
a simple class to control print indentation
Definition: vtkIndent.h:43
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
static double NegInf()
Special IEEE-754 number used to represent negative infinity.
static double Inf()
Special IEEE-754 number used to represent positive infinity.
Superclass for algorithms that produce only vtkMultiBlockDataSet as output.
A subset of a mesh represented as a boolean set operation.
int Operator
The boolean operation that will be performed on the inputs to obtain the output.
BooleanSet(int sId, int op, int *inBegin, int *inEnd)
Construct a new set with the given ID, operator, and inputs.
void PrintNode(ostream &os) override
Print a graphviz node name for use in an edge statement.
std::vector< int > Inputs
A list of input sets. These may be IntervalSets or BooleanSets.
~BooleanSet() override=default
BooleanSet * GetBooleanSetPointer() override
Avoid dynamic_casts. Subclasses must override.
A subset of a mesh represented by a range of acceptable attribute values.
void PrintNode(ostream &os) override
Print a graphviz node name for use in an edge statement.
~Interval() override=default
int Match(double cellNorm[2])
Does the specified range fall inside the interval? For cell-centered attributes, only cellNorm[0] is ...
NormKey Norm
This contains information about the attribute over which the interval is defined.
Interval * GetIntervalPointer() override
A class with comparison operator used to index input array norms used in threshold rules.
void ComputeNorm(vtkIdType cellId, vtkCell *cell, vtkDataArray *array, double cellNorm[2]) const
Compute the norm of a cell by calling NormFunction for all its points or for its single cell-centered...
bool operator<(const NormKey &other) const
A partial ordering of NormKey objects is required for them to serve as keys in the vtkMultiThreshold:...
A base class for representing threshold sets.
int OutputId
A unique identifier for this set.
virtual Interval * GetIntervalPointer()
virtual void PrintNodeName(ostream &os)
Print a graphviz node label statement (with fancy node name and shape).
virtual void PrintNode(ostream &os)=0
Print a graphviz node name for use in an edge statement.
virtual BooleanSet * GetBooleanSetPointer()
Avoid dynamic_casts. Subclasses must override.
Set()
The index of the output mesh that will hold this set or -1 if the set is not output.
virtual ~Set()=default
Virtual destructor since we have virtual members.
Threshold cells within multiple intervals.
~vtkMultiThreshold() override
TruthTree DependentSets
A list of boolean sets whose values depend on the given set.
int NumberOfOutputs
The number of output datasets.
vtkMultiThreshold(const vtkMultiThreshold &)=delete
static vtkMultiThreshold * New()
void operator=(const vtkMultiThreshold &)=delete
Ruling
When an interval is evaluated, its value is used to update a truth table.
Norm
Norms that can be used to threshold vector attributes.
std::map< NormKey, IntervalList > RuleMap
A map describing the IntervalSets that share a common attribute and norm.
std::vector< Set * > Sets
A list of rules keyed by their unique integer ID.
std::vector< Interval * > IntervalList
A list of pointers to IntervalSets.
int NextArrayIndex
A variable used to store the next index to use when calling SetInputArrayToProcess.
RuleMap IntervalRules
A set of threshold rules sorted by the attribute+norm to which they are applied.
int AddBandpassIntervalSet(double xmin, double xmax, int assoc, const char *arrayName, int component, int allScalars)
These convenience members make it easy to insert closed intervals.
std::vector< int > TruthTreeValues
void Reset()
Remove all the intervals currently defined.
int AddLowpassIntervalSet(double xmax, int assoc, const char *arrayName, int component, int allScalars)
These convenience members make it easy to insert closed intervals.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
void PrintGraph(ostream &os)
Print out a graphviz-formatted text description of all the sets.
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This function performs the actual thresholding.
int AddBooleanSet(int operation, int numInputs, int *inputs)
Create a new mesh subset using boolean operations on pre-existing sets.
int OutputSet(int setId)
Create an output mesh containing a boolean or interval subset of the input mesh.
void UpdateDependents(int id, std::set< int > &unresolvedOutputs, TruthTreeValues &setStates, vtkCellData *inCellData, vtkIdType inCell, vtkGenericCell *cell, std::vector< vtkUnstructuredGrid * > &outv)
Recursively update the setStates and unresolvedOutputs vectors based on this->DependentSets.
int AddIntervalSet(double xmin, double xmax, int omin, int omax, int assoc, int attribType, int component, int allScalars)
Add a mesh subset to be computed by thresholding an attribute of the input mesh.
int AddIntervalSet(double xmin, double xmax, int omin, int omax, int assoc, const char *arrayName, int component, int allScalars)
Add a mesh subset to be computed by thresholding an attribute of the input mesh.
int AddNotchIntervalSet(double xlo, double xhi, int assoc, const char *arrayName, int component, int allScalars)
These convenience members make it easy to insert closed intervals.
SetOperation
Operations that can be performed on sets to generate another set.
@ WOR
Include elements that belong to an odd number of input sets (a kind of "winding XOR")
@ XOR
Include an element if it belongs to exactly one input set.
@ NAND
Only include elements that don't belong to any input set.
@ AND
Only include an element if it belongs to all the input sets.
@ OR
Include an element if it belongs to any input set.
int AddHighpassIntervalSet(double xmin, int assoc, const char *arrayName, int component, int allScalars)
These convenience members make it easy to insert closed intervals.
Closure
Whether the endpoint value of an interval should be included or excluded.
@ CLOSED
Specify a closed interval.
int FillInputPortInformation(int port, vtkInformation *info) override
We accept any mesh that is descended from vtkPointSet.
std::vector< TruthTreeValues > TruthTree
int AddIntervalSet(NormKey &nk, double xmin, double xmax, int omin, int omax)
A utility method called by the public AddInterval members.
concrete class for storing a set of points
Definition: vtkPointSet.h:76
dataset represents arbitrary combinations of all possible cell types
@ component
Definition: vtkX3D.h:181
@ info
Definition: vtkX3D.h:382
@ port
Definition: vtkX3D.h:453
@ string
Definition: vtkX3D.h:496
int vtkIdType
Definition: vtkType.h:332