VTK  9.1.0
vtkReflectionFilter.h
Go to the documentation of this file.
1/*=========================================================================
2
3 Program: Visualization Toolkit
4 Module: vtkReflectionFilter.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=========================================================================*/
44#ifndef vtkReflectionFilter_h
45#define vtkReflectionFilter_h
46
48#include "vtkFiltersGeneralModule.h" // For export macro
49
51class vtkDataSet;
52
53class VTKFILTERSGENERAL_EXPORT vtkReflectionFilter : public vtkDataObjectAlgorithm
54{
55public:
57
59 void PrintSelf(ostream& os, vtkIndent indent) override;
60
62 {
63 USE_X_MIN = 0,
64 USE_Y_MIN = 1,
65 USE_Z_MIN = 2,
66 USE_X_MAX = 3,
67 USE_Y_MAX = 4,
68 USE_Z_MAX = 5,
69 USE_X = 6,
70 USE_Y = 7,
71 USE_Z = 8
72 };
73
75
78 vtkSetClampMacro(Plane, int, 0, 8);
79 vtkGetMacro(Plane, int);
80 void SetPlaneToX() { this->SetPlane(USE_X); }
81 void SetPlaneToY() { this->SetPlane(USE_Y); }
82 void SetPlaneToZ() { this->SetPlane(USE_Z); }
83 void SetPlaneToXMin() { this->SetPlane(USE_X_MIN); }
84 void SetPlaneToYMin() { this->SetPlane(USE_Y_MIN); }
85 void SetPlaneToZMin() { this->SetPlane(USE_Z_MIN); }
86 void SetPlaneToXMax() { this->SetPlane(USE_X_MAX); }
87 void SetPlaneToYMax() { this->SetPlane(USE_Y_MAX); }
88 void SetPlaneToZMax() { this->SetPlane(USE_Z_MAX); }
90
92
96 vtkSetMacro(Center, double);
97 vtkGetMacro(Center, double);
99
101
105 vtkSetMacro(CopyInput, vtkTypeBool);
106 vtkGetMacro(CopyInput, vtkTypeBool);
107 vtkBooleanMacro(CopyInput, vtkTypeBool);
109
111
118 vtkSetMacro(FlipAllInputArrays, bool);
119 vtkGetMacro(FlipAllInputArrays, bool);
120 vtkBooleanMacro(FlipAllInputArrays, bool);
122
123protected:
126
133
137 virtual int RequestDataInternal(vtkDataSet* input, vtkUnstructuredGrid* output, double bounds[6]);
138
142 virtual int ComputeBounds(vtkDataObject* input, double bounds[6]);
143
148 vtkDataSet* input, vtkUnstructuredGrid* output, vtkIdType cellId, vtkIdType numInputPoints);
149
152
153 void FlipTuple(double* tuple, int* mirrorDir, int nComp);
154
155 int Plane;
156 double Center;
159
160private:
162 void operator=(const vtkReflectionFilter&) = delete;
163};
164
165#endif
Superclass for algorithms that produce only data object as output.
general representation of visualization data
abstract class to specify dataset behavior
Definition: vtkDataSet.h:166
a simple class to control print indentation
Definition: vtkIndent.h:113
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
reflects a data set across a plane
virtual vtkIdType ReflectNon3DCell(vtkDataSet *input, vtkUnstructuredGrid *output, vtkIdType cellId, vtkIdType numInputPoints)
Generate new, non-3D cell and return the generated cells id.
void SetPlaneToY()
Set the normal of the plane to use as mirror.
void SetPlaneToZMin()
Set the normal of the plane to use as mirror.
virtual int RequestDataInternal(vtkDataSet *input, vtkUnstructuredGrid *output, double bounds[6])
Actual implementation for reflection.
void SetPlaneToX()
Set the normal of the plane to use as mirror.
void SetPlaneToYMin()
Set the normal of the plane to use as mirror.
void FlipTuple(double *tuple, int *mirrorDir, int nComp)
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
void SetPlaneToXMax()
Set the normal of the plane to use as mirror.
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
void SetPlaneToZMax()
Set the normal of the plane to use as mirror.
void SetPlaneToYMax()
Set the normal of the plane to use as mirror.
virtual int ComputeBounds(vtkDataObject *input, double bounds[6])
Internal method to compute bounds.
static vtkReflectionFilter * New()
void SetPlaneToZ()
Set the normal of the plane to use as mirror.
int RequestDataObject(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called by the superclass.
int FillInputPortInformation(int port, vtkInformation *info) override
Fill the input port information objects for this algorithm.
void SetPlaneToXMin()
Set the normal of the plane to use as mirror.
~vtkReflectionFilter() override
dataset represents arbitrary combinations of all possible cell types
@ info
Definition: vtkX3D.h:382
@ port
Definition: vtkX3D.h:453
int vtkTypeBool
Definition: vtkABI.h:69
int vtkIdType
Definition: vtkType.h:332