VTK
/Users/kitware/Dashboards/MyTests/VTK_BLD_Release_docs/Utilities/Doxygen/dox/Filters/General/vtkMultiThreshold.h
Go to the documentation of this file.
00001 /*=========================================================================
00002 
00003   Program:   Visualization Toolkit
00004   Module:    vtkMultiThreshold.h
00005 
00006   Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
00007   All rights reserved.
00008   See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
00009 
00010      This software is distributed WITHOUT ANY WARRANTY; without even
00011      the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
00012      PURPOSE.  See the above copyright notice for more information.
00013 
00014 =========================================================================*/
00015 /*----------------------------------------------------------------------------
00016  Copyright (c) Sandia Corporation
00017  See Copyright.txt or http://www.paraview.org/HTML/Copyright.html for details.
00018 ----------------------------------------------------------------------------*/
00019 
00103 #ifndef vtkMultiThreshold_h
00104 #define vtkMultiThreshold_h
00105 
00106 #include "vtkFiltersGeneralModule.h" // For export macro
00107 #include "vtkMultiBlockDataSetAlgorithm.h"
00108 #include "vtkMath.h" // for Inf() and NegInf()
00109 
00110 #include <vector> // for lists of threshold rules
00111 #include <map> // for IntervalRules map
00112 #include <set> // for UpdateDependents()
00113 #include <string> // for holding array names in NormKey
00114 
00115 class vtkCell;
00116 class vtkCellData;
00117 class vtkDataArray;
00118 class vtkGenericCell;
00119 class vtkPointSet;
00120 class vtkUnstructuredGrid;
00121 
00122 class VTKFILTERSGENERAL_EXPORT vtkMultiThreshold : public vtkMultiBlockDataSetAlgorithm
00123 {
00124 public:
00125   vtkTypeMacro(vtkMultiThreshold,vtkMultiBlockDataSetAlgorithm);
00126   static vtkMultiThreshold* New();
00127   virtual void PrintSelf( ostream& os, vtkIndent indent );
00128 
00129   //BTX
00131   enum Closure {
00132     OPEN=0,   
00133     CLOSED=1  
00134   };
00136   enum Norm {
00137     LINFINITY_NORM=-3, 
00138     L2_NORM=-2,        
00139     L1_NORM=-1         
00140   };
00142   enum SetOperation {
00143     AND, 
00144     OR, 
00145     XOR, 
00146     WOR, 
00147     NAND 
00148   };
00149   //ETX
00150 
00152 
00195   int AddIntervalSet( double xmin, double xmax, int omin, int omax,
00196     int assoc, const char* arrayName, int component, int allScalars );
00197   int AddIntervalSet( double xmin, double xmax, int omin, int omax,
00198     int assoc, int attribType, int component, int allScalars );
00200 
00202 
00208   int AddLowpassIntervalSet( double xmax, int assoc, const char* arrayName, int component, int allScalars );
00209   int AddHighpassIntervalSet( double xmin, int assoc, const char* arrayName, int component, int allScalars );
00210   int AddBandpassIntervalSet( double xmin, double xmax, int assoc, const char* arrayName, int component, int allScalars );
00211   int AddNotchIntervalSet( double xlo, double xhi, int assoc, const char* arrayName, int component, int allScalars );
00213 
00216   int AddBooleanSet( int operation, int numInputs, int* inputs );
00217 
00220   int OutputSet( int setId );
00221 
00223   void Reset();
00224 
00225   //BTX
00227   typedef double (*TupleNorm)( vtkDataArray* arr, vtkIdType tuple, int component );
00228 
00229   // NormKey must be able to use TupleNorm typedef:
00230   class NormKey;
00231 
00232   // Interval must be able to use NormKey typedef:
00233   class Interval;
00234 
00235   // Set needs to refer to boolean set pointers
00236   class BooleanSet;
00237 
00239   class NormKey {
00240   public:
00241     int Association; // vtkDataObject::FIELD_ASSOCIATION_POINTS or vtkDataObject::FIELD_ASSOCIATION_CELLS
00242     int Type; // -1 => use Name, otherwise: vtkDataSetAttributes::{SCALARS, VECTORS, TENSORS, NORMALS, TCOORDS, GLOBALIDS}
00243     std::string Name; // Either empty or (when ArrayType == -1) an input array name
00244     int Component; // LINFINITY_NORM, L1_NORM, L2_NORM or an integer component number
00245     int AllScalars; // For Association == vtkDataObject::FIELD_ASSOCIATION_POINTS, must all points be in the interval?
00246     int InputArrayIndex; // The number passed to vtkAlgorithm::SetInputArrayToProcess()
00247     TupleNorm NormFunction; // A function pointer to compute the norm (or fetcht the correct component) of a tuple.
00248 
00250     void ComputeNorm( vtkIdType cellId, vtkCell* cell, vtkDataArray* array, double cellNorm[2] ) const;
00251 
00253     bool operator < ( const NormKey& other ) const {
00254       if ( this->Association < other.Association )
00255         return true;
00256       else if ( this->Association > other.Association )
00257         return false;
00258 
00259       if ( this->Component < other.Component )
00260         return true;
00261       else if ( this->Component > other.Component )
00262         return false;
00263 
00264       if ( (! this->AllScalars) && other.AllScalars )
00265         return true;
00266       else if ( this->AllScalars && (! other.AllScalars) )
00267         return false;
00268 
00269       if ( this->Type == -1 )
00270         {
00271         if ( other.Type == -1 )
00272           return this->Name < other.Name;
00273         return true;
00274         }
00275       else
00276         {
00277         return this->Type < other.Type;
00278         }
00279     }
00280   };
00281 
00286   class Set {
00287   public:
00288     int Id; 
00289     int OutputId; 
00290 
00292     Set() {
00293       this->OutputId = -1;
00294     }
00296     virtual ~Set() { }
00298     virtual void PrintNodeName( ostream& os );
00300     virtual void PrintNode( ostream& os ) = 0;
00302     virtual BooleanSet* GetBooleanSetPointer();
00303     virtual Interval* GetIntervalPointer();
00304   };
00305 
00307   class Interval : public Set {
00308   public:
00310     double EndpointValues[2];
00312     int EndpointClosures[2];
00314     NormKey Norm;
00315 
00320     int Match( double cellNorm[2] );
00321 
00322     virtual ~Interval() { }
00323     virtual void PrintNode( ostream& os );
00324     virtual Interval* GetIntervalPointer();
00325   };
00326 
00328   class BooleanSet : public Set {
00329   public:
00331     int Operator;
00333     std::vector<int> Inputs;
00334 
00336     BooleanSet( int sId, int op, int* inBegin, int* inEnd ) : Inputs( inBegin, inEnd ) {
00337       this->Id = sId;
00338       this->Operator = op;
00339     }
00340     virtual ~BooleanSet() { }
00341     virtual void PrintNode( ostream& os );
00342     virtual BooleanSet* GetBooleanSetPointer();
00343   };
00344   //ETX
00345 
00346 protected:
00347 
00348   vtkMultiThreshold();
00349   virtual ~vtkMultiThreshold();
00350 
00351   //BTX
00353 
00364   enum Ruling {
00365     INCONCLUSIVE=-1,
00366     INCLUDE=-2,
00367     EXCLUDE=-3
00368   };
00369   //ETX
00371 
00373   virtual int RequestData( vtkInformation*, vtkInformationVector**, vtkInformationVector* );
00374 
00378   virtual int FillInputPortInformation( int port, vtkInformation* info );
00379 
00384   int NextArrayIndex;
00385 
00387   int NumberOfOutputs;
00388 
00389   //BTX
00391   typedef std::vector<Interval*> IntervalList;
00393   typedef std::map<NormKey,IntervalList> RuleMap;
00394 
00395   typedef std::vector<int> TruthTreeValues;
00396   typedef std::vector<TruthTreeValues> TruthTree;
00397 
00400   RuleMap IntervalRules;
00401 
00405   std::vector<Set*> Sets;
00406 
00412   TruthTree DependentSets;
00413 
00415 
00417   void UpdateDependents(
00418     int id, std::set<int>& unresolvedOutputs, TruthTreeValues& setStates,
00419     vtkCellData* inCellData, vtkIdType cellId, vtkGenericCell* cell, std::vector<vtkUnstructuredGrid*>& outv );
00421 
00423   int AddIntervalSet( NormKey& nk, double xmin, double xmax, int omin, int omax );
00424 
00425   //ETX
00426 
00428   void PrintGraph( ostream& os );
00429 
00430   vtkMultiThreshold( const vtkMultiThreshold& ); // Not implemented.
00431   void operator = ( const vtkMultiThreshold& ); // Not implemented.
00432 };
00433 
00434 inline int vtkMultiThreshold::AddLowpassIntervalSet( double xmax, int assoc, const char* arrayName, int component, int allScalars )
00435 {
00436   return this->AddIntervalSet( vtkMath::NegInf(), xmax, CLOSED, CLOSED, assoc, arrayName, component, allScalars );
00437 }
00438 
00439 inline int vtkMultiThreshold::AddHighpassIntervalSet( double xmin, int assoc, const char* arrayName, int component, int allScalars )
00440 {
00441   return this->AddIntervalSet( xmin, vtkMath::Inf(), CLOSED, CLOSED, assoc, arrayName, component, allScalars );
00442 }
00443 
00444 inline int vtkMultiThreshold::AddBandpassIntervalSet(
00445   double xmin, double xmax, int assoc, const char* arrayName, int component, int allScalars )
00446 {
00447   return this->AddIntervalSet( xmin, xmax, CLOSED, CLOSED, assoc, arrayName, component, allScalars );
00448 }
00449 
00450 inline int vtkMultiThreshold::AddNotchIntervalSet(
00451   double xlo, double xhi, int assoc, const char* arrayName, int component, int allScalars )
00452 {
00453   int band = this->AddIntervalSet( xlo, xhi, CLOSED, CLOSED, assoc, arrayName, component, allScalars );
00454   if ( band < 0 )
00455     {
00456     return -1;
00457     }
00458   return this->AddBooleanSet( NAND, 1, &band );
00459 }
00460 
00461 inline vtkMultiThreshold::Interval* vtkMultiThreshold::Set::GetIntervalPointer()
00462 {
00463   return 0;
00464 }
00465 
00466 inline vtkMultiThreshold::BooleanSet* vtkMultiThreshold::Set::GetBooleanSetPointer()
00467 {
00468   return 0;
00469 }
00470 
00471 inline vtkMultiThreshold::Interval* vtkMultiThreshold::Interval::GetIntervalPointer()
00472 {
00473   return this;
00474 }
00475 
00476 inline vtkMultiThreshold::BooleanSet* vtkMultiThreshold::BooleanSet::GetBooleanSetPointer()
00477 {
00478   return this;
00479 }
00480 
00481 #endif // vtkMultiThreshold_h