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
132#ifndef vtkMultiThreshold_h
133#define vtkMultiThreshold_h
134
135#include "vtkFiltersGeneralModule.h" // For export macro
136#include "vtkMath.h" // for Inf() and NegInf()
138
139#include <map> // for IntervalRules map
140#include <set> // for UpdateDependents()
141#include <string> // for holding array names in NormKey
142#include <vector> // for lists of threshold rules
143
144class vtkCell;
145class vtkCellData;
146class vtkDataArray;
147class vtkGenericCell;
148class vtkPointSet;
150
151class VTKFILTERSGENERAL_EXPORT vtkMultiThreshold : public vtkMultiBlockDataSetAlgorithm
152{
153public:
156 void PrintSelf(ostream& os, vtkIndent indent) override;
157
160 {
161 OPEN = 0,
162 CLOSED = 1
163 };
165 enum Norm
166 {
167 LINFINITY_NORM = -3,
168 L2_NORM = -2,
169 L1_NORM = -1
170 };
174 {
179 NAND
180 };
181
183
235 int AddIntervalSet(double xmin, double xmax, int omin, int omax, int assoc, const char* arrayName,
236 int component, int allScalars);
237 int AddIntervalSet(double xmin, double xmax, int omin, int omax, int assoc, int attribType,
238 int component, int allScalars);
240
242
249 int AddLowpassIntervalSet(
250 double xmax, int assoc, const char* arrayName, int component, int allScalars);
251 int AddHighpassIntervalSet(
252 double xmin, int assoc, const char* arrayName, int component, int allScalars);
253 int AddBandpassIntervalSet(
254 double xmin, double xmax, int assoc, const char* arrayName, int component, int allScalars);
255 int AddNotchIntervalSet(
256 double xlo, double xhi, int assoc, const char* arrayName, int component, int allScalars);
258
262 int AddBooleanSet(int operation, int numInputs, int* inputs);
263
267 int OutputSet(int setId);
268
272 void Reset();
273
276 typedef double (*TupleNorm)(vtkDataArray* arr, vtkIdType tuple, int component);
277
278 // NormKey must be able to use TupleNorm typedef:
279 class NormKey;
280
281 // Interval must be able to use NormKey typedef:
282 class Interval;
283
284 // Set needs to refer to boolean set pointers
285 class BooleanSet;
286
289 {
290 public:
291 int Association; // vtkDataObject::FIELD_ASSOCIATION_POINTS or
292 // vtkDataObject::FIELD_ASSOCIATION_CELLS
293 int Type; // -1 => use Name, otherwise: vtkDataSetAttributes::{SCALARS, VECTORS, TENSORS,
294 // NORMALS, TCOORDS, GLOBALIDS}
295 std::string Name; // Either empty or (when ArrayType == -1) an input array name
296 int Component; // LINFINITY_NORM, L1_NORM, L2_NORM or an integer component number
297 int AllScalars; // For Association == vtkDataObject::FIELD_ASSOCIATION_POINTS, must all points
298 // be in the interval?
299 int InputArrayIndex; // The number passed to vtkAlgorithm::SetInputArrayToProcess()
300 TupleNorm NormFunction; // A function pointer to compute the norm (or fetcht the correct
301 // component) of a tuple.
302
306 vtkIdType cellId, vtkCell* cell, vtkDataArray* array, double cellNorm[2]) const;
307
310 bool operator<(const NormKey& other) const
311 {
312 if (this->Association < other.Association)
313 return true;
314 else if (this->Association > other.Association)
315 return false;
316
317 if (this->Component < other.Component)
318 return true;
319 else if (this->Component > other.Component)
320 return false;
321
322 if ((!this->AllScalars) && other.AllScalars)
323 return true;
324 else if (this->AllScalars && (!other.AllScalars))
325 return false;
326
327 if (this->Type == -1)
328 {
329 if (other.Type == -1)
330 return this->Name < other.Name;
331 return true;
332 }
333 else
334 {
335 return this->Type < other.Type;
336 }
337 }
338 };
339
344 class Set
345 {
346 public:
347 int Id;
350
353 Set() { this->OutputId = -1; }
355 virtual ~Set() = default;
357 virtual void PrintNodeName(ostream& os);
359 virtual void PrintNode(ostream& os) = 0;
361 virtual BooleanSet* GetBooleanSetPointer();
362 virtual Interval* GetIntervalPointer();
363 };
364
366 class Interval : public Set
367 {
368 public:
370 double EndpointValues[2];
372 int EndpointClosures[2];
375
381 int Match(double cellNorm[2]);
382
383 ~Interval() override = default;
384 void PrintNode(ostream& os) override;
385 Interval* GetIntervalPointer() override;
386 };
387
389 class BooleanSet : public Set
390 {
391 public:
395 std::vector<int> Inputs;
396
398 BooleanSet(int sId, int op, int* inBegin, int* inEnd)
399 : Inputs(inBegin, inEnd)
400 {
401 this->Id = sId;
402 this->Operator = op;
403 }
404 ~BooleanSet() override = default;
405 void PrintNode(ostream& os) override;
406 BooleanSet* GetBooleanSetPointer() override;
407 };
408
409protected:
412
428 {
429 INCONCLUSIVE = -1,
430 INCLUDE = -2,
431 EXCLUDE = -3
432 };
433
438
445
452
457
459 typedef std::vector<Interval*> IntervalList;
461 typedef std::map<NormKey, IntervalList> RuleMap;
462
463 typedef std::vector<int> TruthTreeValues;
464 typedef std::vector<TruthTreeValues> TruthTree;
465
470
476 std::vector<Set*> Sets;
477
485
489 void UpdateDependents(int id, std::set<int>& unresolvedOutputs, TruthTreeValues& setStates,
490 vtkCellData* inCellData, vtkIdType inCell, vtkGenericCell* cell,
491 std::vector<vtkUnstructuredGrid*>& outv);
492
496 int AddIntervalSet(NormKey& nk, double xmin, double xmax, int omin, int omax);
497
501 void PrintGraph(ostream& os);
502
504 void operator=(const vtkMultiThreshold&) = delete;
505};
506
508 double xmax, int assoc, const char* arrayName, int component, int allScalars)
509{
510 return this->AddIntervalSet(
511 vtkMath::NegInf(), xmax, CLOSED, CLOSED, assoc, arrayName, component, allScalars);
512}
513
515 double xmin, int assoc, const char* arrayName, int component, int allScalars)
516{
517 return this->AddIntervalSet(
518 xmin, vtkMath::Inf(), CLOSED, CLOSED, assoc, arrayName, component, allScalars);
519}
520
522 double xmin, double xmax, int assoc, const char* arrayName, int component, int allScalars)
523{
524 return this->AddIntervalSet(xmin, xmax, CLOSED, CLOSED, assoc, arrayName, component, allScalars);
525}
526
528 double xlo, double xhi, int assoc, const char* arrayName, int component, int allScalars)
529{
530 int band =
531 this->AddIntervalSet(xlo, xhi, CLOSED, CLOSED, assoc, arrayName, component, allScalars);
532 if (band < 0)
533 {
534 return -1;
535 }
536 return this->AddBooleanSet(NAND, 1, &band);
537}
538
540{
541 return nullptr;
542}
543
545{
546 return nullptr;
547}
548
550{
551 return this;
552}
553
555{
556 return this;
557}
558
559#endif // vtkMultiThreshold_h
represent and manipulate cell attribute data
Definition: vtkCellData.h:142
abstract class to specify cell behavior
Definition: vtkCell.h:147
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:159
provides thread-safe access to cells
a simple class to control print indentation
Definition: vtkIndent.h:113
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:106
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