VTK
vtkExtractCTHPart.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkExtractCTHPart.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 =========================================================================*/
43 #ifndef vtkExtractCTHPart_h
44 #define vtkExtractCTHPart_h
45 
46 #include "vtkFiltersParallelModule.h" // For export macro
48 
49 class vtkAppendPolyData;
50 class vtkContourFilter;
51 class vtkDataArray;
52 class vtkDataSet;
54 class vtkDoubleArray;
55 class vtkExtractCTHPartInternal;
56 class vtkImageData;
59 class vtkPlane;
60 class vtkPolyData;
61 class vtkRectilinearGrid;
62 class vtkUniformGrid;
64 class vtkExtractCTHPartFragments;
65 
66 //#define EXTRACT_USE_IMAGE_DATA 1
67 
68 class VTKFILTERSPARALLEL_EXPORT vtkExtractCTHPart : public vtkMultiBlockDataSetAlgorithm
69 {
70 public:
73  void PrintSelf(ostream& os, vtkIndent indent);
74 
76 
79  void AddVolumeArrayName(const char*);
82  const char* GetVolumeArrayName(int idx);
84 
86 
92  vtkGetObjectMacro(Controller,vtkMultiProcessController);
94 
96 
99  vtkSetMacro(Capping, bool);
100  vtkGetMacro(Capping, bool);
101  vtkBooleanMacro(Capping, bool);
103 
105 
109  vtkSetMacro(GenerateTriangles, bool);
110  vtkGetMacro(GenerateTriangles, bool);
111  vtkBooleanMacro(GenerateTriangles, bool);
113 
115 
120  vtkSetMacro(RemoveGhostCells, bool);
121  vtkGetMacro(RemoveGhostCells, bool);
122  vtkBooleanMacro(RemoveGhostCells, bool);
124 
126 
129  void SetClipPlane(vtkPlane *clipPlane);
130  vtkGetObjectMacro(ClipPlane, vtkPlane);
132 
137 
139 
143  vtkSetClampMacro(VolumeFractionSurfaceValue, double, 0.0, 1.0);
144  vtkGetMacro(VolumeFractionSurfaceValue, double);
146 
147 protected:
150 
152  virtual int RequestData(
154 
160 
166  vtkPolyData* output, vtkCompositeDataSet* input, const char*arrayName);
167 
169  vtkPolyData *output,
170  int maxFlag,
171  int originExtents[3],
172  int ext[6],
173  int aAxis,
174  int bAxis,
175  int cAxis);
176 
182  int IsGhostFace(int axis0,
183  int maxFlag,
184  int dims[3],
185  vtkUnsignedCharArray *ghostArray);
186 
187  void TriggerProgressEvent(double val);
188 
193  bool Capping;
197 private:
198  vtkExtractCTHPart(const vtkExtractCTHPart&) VTK_DELETE_FUNCTION;
199  void operator=(const vtkExtractCTHPart&) VTK_DELETE_FUNCTION;
200 
201  class VectorOfFragments;
202 
207  template <class T>
208  bool ExtractClippedContourOnBlock(
209  vtkExtractCTHPart::VectorOfFragments& fragments, T* input, const char* arrayName);
210 
215  template <class T>
216  bool ExtractContourOnBlock(
217  vtkExtractCTHPart::VectorOfFragments& fragments, T* input, const char* arrayName);
218 
223  template <class T>
224  void ExtractExteriorSurface(
225  vtkExtractCTHPart::VectorOfFragments& fragments, T* input);
226 
230  void ExecuteCellDataToPointData(
231  vtkDataArray *cellVolumeFraction, vtkDoubleArray *pointVolumeFraction, const int *dims);
232 
233  double ProgressShift;
234  double ProgressScale;
235 
236  class ScaledProgress;
237  friend class ScaledProgress;
238  vtkExtractCTHPartInternal* Internals;
239 };
240 #endif
appends one or more polygonal datasets together
abstract superclass for composite (multi-block or AMR) datasets
generate isosurfaces/isolines from scalar values
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:55
Extracts outer (polygonal) surface.
abstract class to specify dataset behavior
Definition: vtkDataSet.h:63
dynamic, self-adjusting array of double
Generates surface of a CTH volume fraction.
bool ExtractContour(vtkPolyData *output, vtkCompositeDataSet *input, const char *arrayName)
Extract contour for a particular array over the entire input dataset.
vtkMultiProcessController * Controller
virtual int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *)
This is called by the superclass.
void SetClipPlane(vtkPlane *clipPlane)
Set, get or manipulate the implicit clipping plane.
int GetNumberOfVolumeArrayNames()
static vtkExtractCTHPart * New()
void RemoveVolumeArrayNames()
int IsGhostFace(int axis0, int maxFlag, int dims[3], vtkUnsignedCharArray *ghostArray)
Is block face on axis0 (either min or max depending on the maxFlag) composed of only ghost cells?
vtkMTimeType GetMTime()
Look at clip plane to compute MTime.
virtual int FillInputPortInformation(int port, vtkInformation *info)
Fill the input port information objects for this algorithm.
void ExecuteFaceQuads(vtkDataSet *input, vtkPolyData *output, int maxFlag, int originExtents[3], int ext[6], int aAxis, int bAxis, int cAxis)
void PrintSelf(ostream &os, vtkIndent indent)
Methods invoked by print to print information about the object including superclasses.
void SetController(vtkMultiProcessController *controller)
Get/Set the parallel controller.
bool ComputeGlobalBounds(vtkCompositeDataSet *input)
Compute the bounds over the composite dataset, some sub-dataset can be on other processors.
const char * GetVolumeArrayName(int idx)
void TriggerProgressEvent(double val)
double VolumeFractionSurfaceValueInternal
void AddVolumeArrayName(const char *)
Select cell-data arrays (volume-fraction arrays) to contour with.
topologically and geometrically regular array of data
Definition: vtkImageData.h:46
a simple class to control print indentation
Definition: vtkIndent.h:40
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
Superclass for algorithms that produce only vtkMultiBlockDataSet as output.
Multiprocessing communication superclass.
perform various plane computations
Definition: vtkPlane.h:38
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition: vtkPolyData.h:86
a dataset that is topologically regular with variable spacing in the three coordinate directions
image data with blanking
dynamic, self-adjusting array of unsigned char
@ info
Definition: vtkX3D.h:376
@ port
Definition: vtkX3D.h:447
vtkSetMacro(IgnoreDriverBugs, bool)
Updates the extensions string.
vtkBooleanMacro(IgnoreDriverBugs, bool)
Updates the extensions string.
vtkTypeUInt64 vtkMTimeType
Definition: vtkType.h:248