VTK
|
00001 /*========================================================================= 00002 00003 Program: Visualization Toolkit 00004 Module: $RCSfile: vtkConnectivityFilter.h,v $ 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 =========================================================================*/ 00052 #ifndef __vtkConnectivityFilter_h 00053 #define __vtkConnectivityFilter_h 00054 00055 #include "vtkUnstructuredGridAlgorithm.h" 00056 00057 #define VTK_EXTRACT_POINT_SEEDED_REGIONS 1 00058 #define VTK_EXTRACT_CELL_SEEDED_REGIONS 2 00059 #define VTK_EXTRACT_SPECIFIED_REGIONS 3 00060 #define VTK_EXTRACT_LARGEST_REGION 4 00061 #define VTK_EXTRACT_ALL_REGIONS 5 00062 #define VTK_EXTRACT_CLOSEST_POINT_REGION 6 00063 00064 class vtkDataArray; 00065 class vtkFloatArray; 00066 class vtkIdList; 00067 class vtkIdTypeArray; 00068 class vtkIntArray; 00069 00070 class VTK_GRAPHICS_EXPORT vtkConnectivityFilter : public vtkUnstructuredGridAlgorithm 00071 { 00072 public: 00073 vtkTypeRevisionMacro(vtkConnectivityFilter,vtkUnstructuredGridAlgorithm); 00074 void PrintSelf(ostream& os, vtkIndent indent); 00075 00077 static vtkConnectivityFilter *New(); 00078 00080 00083 vtkSetMacro(ScalarConnectivity,int); 00084 vtkGetMacro(ScalarConnectivity,int); 00085 vtkBooleanMacro(ScalarConnectivity,int); 00087 00089 00091 vtkSetVector2Macro(ScalarRange,double); 00092 vtkGetVector2Macro(ScalarRange,double); 00094 00096 00097 vtkSetClampMacro(ExtractionMode,int, 00098 VTK_EXTRACT_POINT_SEEDED_REGIONS,VTK_EXTRACT_CLOSEST_POINT_REGION); 00099 vtkGetMacro(ExtractionMode,int); 00100 void SetExtractionModeToPointSeededRegions() 00101 {this->SetExtractionMode(VTK_EXTRACT_POINT_SEEDED_REGIONS);}; 00102 void SetExtractionModeToCellSeededRegions() 00103 {this->SetExtractionMode(VTK_EXTRACT_CELL_SEEDED_REGIONS);}; 00104 void SetExtractionModeToLargestRegion() 00105 {this->SetExtractionMode(VTK_EXTRACT_LARGEST_REGION);}; 00106 void SetExtractionModeToSpecifiedRegions() 00107 {this->SetExtractionMode(VTK_EXTRACT_SPECIFIED_REGIONS);}; 00108 void SetExtractionModeToClosestPointRegion() 00109 {this->SetExtractionMode(VTK_EXTRACT_CLOSEST_POINT_REGION);}; 00110 void SetExtractionModeToAllRegions() 00111 {this->SetExtractionMode(VTK_EXTRACT_ALL_REGIONS);}; 00112 const char *GetExtractionModeAsString(); 00114 00116 void InitializeSeedList(); 00117 00119 void AddSeed(vtkIdType id); 00120 00122 void DeleteSeed(vtkIdType id); 00123 00125 void InitializeSpecifiedRegionList(); 00126 00128 void AddSpecifiedRegion(int id); 00129 00131 void DeleteSpecifiedRegion(int id); 00132 00134 00136 vtkSetVector3Macro(ClosestPoint,double); 00137 vtkGetVectorMacro(ClosestPoint,double,3); 00139 00141 int GetNumberOfExtractedRegions(); 00142 00144 00145 vtkSetMacro(ColorRegions,int); 00146 vtkGetMacro(ColorRegions,int); 00147 vtkBooleanMacro(ColorRegions,int); 00149 00150 protected: 00151 vtkConnectivityFilter(); 00152 ~vtkConnectivityFilter(); 00153 00154 // Usual data generation method 00155 virtual int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *); 00156 virtual int FillInputPortInformation(int port, vtkInformation *info); 00157 00158 int ColorRegions; //boolean turns on/off scalar gen for separate regions 00159 int ExtractionMode; //how to extract regions 00160 vtkIdList *Seeds; //id's of points or cells used to seed regions 00161 vtkIdList *SpecifiedRegionIds; //regions specified for extraction 00162 vtkIdTypeArray *RegionSizes; //size (in cells) of each region extracted 00163 00164 double ClosestPoint[3]; 00165 00166 int ScalarConnectivity; 00167 double ScalarRange[2]; 00168 00169 void TraverseAndMark(vtkDataSet *input); 00170 00171 private: 00172 // used to support algorithm execution 00173 vtkFloatArray *CellScalars; 00174 vtkIdList *NeighborCellPointIds; 00175 vtkIdType *Visited; 00176 vtkIdType *PointMap; 00177 vtkIdTypeArray *NewScalars; 00178 vtkIdTypeArray *NewCellScalars; 00179 vtkIdType RegionNumber; 00180 vtkIdType PointNumber; 00181 vtkIdType NumCellsInRegion; 00182 vtkDataArray *InScalars; 00183 vtkIdList *Wave; 00184 vtkIdList *Wave2; 00185 vtkIdList *PointIds; 00186 vtkIdList *CellIds; 00187 private: 00188 vtkConnectivityFilter(const vtkConnectivityFilter&); // Not implemented. 00189 void operator=(const vtkConnectivityFilter&); // Not implemented. 00190 }; 00191 00193 inline const char *vtkConnectivityFilter::GetExtractionModeAsString(void) 00194 { 00195 if ( this->ExtractionMode == VTK_EXTRACT_POINT_SEEDED_REGIONS ) 00196 { 00197 return "ExtractPointSeededRegions"; 00198 } 00199 else if ( this->ExtractionMode == VTK_EXTRACT_CELL_SEEDED_REGIONS ) 00200 { 00201 return "ExtractCellSeededRegions"; 00202 } 00203 else if ( this->ExtractionMode == VTK_EXTRACT_SPECIFIED_REGIONS ) 00204 { 00205 return "ExtractSpecifiedRegions"; 00206 } 00207 else if ( this->ExtractionMode == VTK_EXTRACT_ALL_REGIONS ) 00208 { 00209 return "ExtractAllRegions"; 00210 } 00211 else if ( this->ExtractionMode == VTK_EXTRACT_CLOSEST_POINT_REGION ) 00212 { 00213 return "ExtractClosestPointRegion"; 00214 } 00215 else 00216 { 00217 return "ExtractLargestRegion"; 00218 } 00219 } 00220 00221 #endif 00222 00223