VTK  9.2.6
vtkArrayCalculator.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkArrayCalculator.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 =========================================================================*/
75 #ifndef vtkArrayCalculator_h
76 #define vtkArrayCalculator_h
77 
78 #include "vtkDataObject.h" // For attribute types
79 #include "vtkFiltersCoreModule.h" // For export macro
81 #include "vtkTuple.h" // needed for vtkTuple
82 #include <vector> // needed for vector
83 
84 class vtkDataSet;
85 
86 class VTKFILTERSCORE_EXPORT vtkArrayCalculator : public vtkPassInputTypeAlgorithm
87 {
88 public:
90  void PrintSelf(ostream& os, vtkIndent indent) override;
91 
92  static vtkArrayCalculator* New();
93 
95 
98  vtkSetStringMacro(Function);
99  vtkGetStringMacro(Function);
101 
103 
113  void AddScalarArrayName(const char* arrayName, int component = 0);
114  void AddVectorArrayName(
115  const char* arrayName, int component0 = 0, int component1 = 1, int component2 = 2);
117 
119 
125  void AddScalarVariable(const char* variableName, const char* arrayName, int component = 0);
126  void AddVectorVariable(const char* variableName, const char* arrayName, int component0 = 0,
127  int component1 = 1, int component2 = 2);
129 
131 
137  void AddCoordinateScalarVariable(const char* variableName, int component = 0);
138  void AddCoordinateVectorVariable(
139  const char* variableName, int component0 = 0, int component1 = 1, int component2 = 2);
141 
143 
149  vtkSetStringMacro(ResultArrayName);
150  vtkGetStringMacro(ResultArrayName);
152 
154 
158  vtkGetMacro(ResultArrayType, int);
159  vtkSetMacro(ResultArrayType, int);
161 
163 
169  vtkGetMacro(CoordinateResults, vtkTypeBool);
170  vtkSetMacro(CoordinateResults, vtkTypeBool);
171  vtkBooleanMacro(CoordinateResults, vtkTypeBool);
173 
175 
180  vtkGetMacro(ResultNormals, bool);
181  vtkSetMacro(ResultNormals, bool);
182  vtkBooleanMacro(ResultNormals, bool);
184 
186 
191  vtkGetMacro(ResultTCoords, bool);
192  vtkSetMacro(ResultTCoords, bool);
193  vtkBooleanMacro(ResultTCoords, bool);
195 
199  const char* GetAttributeTypeAsString();
200 
201  static const int DEFAULT_ATTRIBUTE_TYPE = -1;
203 
209  vtkSetMacro(AttributeType, int);
210  vtkGetMacro(AttributeType, int);
211  void SetAttributeTypeToDefault() { this->SetAttributeType(DEFAULT_ATTRIBUTE_TYPE); }
212  void SetAttributeTypeToPointData() { this->SetAttributeType(vtkDataObject::POINT); }
213  void SetAttributeTypeToCellData() { this->SetAttributeType(vtkDataObject::CELL); }
214  void SetAttributeTypeToEdgeData() { this->SetAttributeType(vtkDataObject::EDGE); }
215  void SetAttributeTypeToVertexData() { this->SetAttributeType(vtkDataObject::VERTEX); }
216  void SetAttributeTypeToRowData() { this->SetAttributeType(vtkDataObject::ROW); }
218 
222  void RemoveAllVariables();
223 
227  virtual void RemoveScalarVariables();
228 
232  virtual void RemoveVectorVariables();
233 
237  virtual void RemoveCoordinateScalarVariables();
238 
242  virtual void RemoveCoordinateVectorVariables();
243 
245 
248  const std::vector<std::string>& GetScalarArrayNames() { return this->ScalarArrayNames; }
249  std::string GetScalarArrayName(int i);
250  const std::vector<std::string>& GetVectorArrayNames() { return this->VectorArrayNames; }
251  std::string GetVectorArrayName(int i);
252  const std::vector<std::string>& GetScalarVariableNames() { return this->ScalarVariableNames; }
253  std::string GetScalarVariableName(int i);
254  const std::vector<std::string>& GetVectorVariableNames() { return this->VectorVariableNames; }
255  std::string GetVectorVariableName(int i);
256  const std::vector<int>& GetSelectedScalarComponents() { return this->SelectedScalarComponents; }
257  int GetSelectedScalarComponent(int i);
258  const std::vector<vtkTuple<int, 3>>& GetSelectedVectorComponents()
259  {
260  return this->SelectedVectorComponents;
261  }
262  vtkTuple<int, 3> GetSelectedVectorComponents(int i);
263  int GetNumberOfScalarArrays() { return static_cast<int>(this->ScalarArrayNames.size()); }
264  int GetNumberOfVectorArrays() { return static_cast<int>(this->VectorArrayNames.size()); }
266 
268 
274  vtkSetMacro(ReplaceInvalidValues, vtkTypeBool);
275  vtkGetMacro(ReplaceInvalidValues, vtkTypeBool);
276  vtkBooleanMacro(ReplaceInvalidValues, vtkTypeBool);
277  vtkSetMacro(ReplacementValue, double);
278  vtkGetMacro(ReplacementValue, double);
280 
282 
287  vtkSetMacro(IgnoreMissingArrays, bool);
288  vtkGetMacro(IgnoreMissingArrays, bool);
289  vtkBooleanMacro(IgnoreMissingArrays, bool);
291 
296  {
297  FunctionParser, // vtkFunctionParser
298  ExprTkFunctionParser, // vtkExprTkFunctionParser
299  NumberOfFunctionParserTypes
300  };
301 
303 
307  vtkSetEnumMacro(FunctionParserType, FunctionParserTypes);
309  {
310  this->FunctionParserType = FunctionParserTypes::FunctionParser;
311  this->Modified();
312  }
314  {
315  this->FunctionParserType = FunctionParserTypes::ExprTkFunctionParser;
316  this->Modified();
317  }
318  vtkGetEnumMacro(FunctionParserType, FunctionParserTypes);
320 
325  vtkDataSet* GetDataSetOutput();
326 
327 protected:
329  ~vtkArrayCalculator() override;
330 
331  int FillInputPortInformation(int, vtkInformation*) override;
332 
334 
338  int GetAttributeTypeFromInput(vtkDataObject* input);
339 
347  static std::string CheckValidVariableName(const char* variableName);
348 
350 
351  char* Function;
353  std::vector<std::string> ScalarArrayNames;
354  std::vector<std::string> VectorArrayNames;
355  std::vector<std::string> ScalarVariableNames;
356  std::vector<std::string> VectorVariableNames;
358  std::vector<int> SelectedScalarComponents;
359  std::vector<vtkTuple<int, 3>> SelectedVectorComponents;
360 
364 
368  std::vector<std::string> CoordinateScalarVariableNames;
369  std::vector<std::string> CoordinateVectorVariableNames;
371  std::vector<vtkTuple<int, 3>> SelectedCoordinateVectorComponents;
372 
374 
375 private:
376  vtkArrayCalculator(const vtkArrayCalculator&) = delete;
377  void operator=(const vtkArrayCalculator&) = delete;
378 
379  // Do the bulk of the work
380  template <typename TFunctionParser>
381  int ProcessDataObject(vtkDataObject* input, vtkDataObject* output);
382 };
383 
384 #endif
vtkTypeBool ReplaceInvalidValues
const std::vector< vtkTuple< int, 3 > > & GetSelectedVectorComponents()
Methods to get information about the current variables.
std::vector< std::string > ScalarArrayNames
const std::vector< std::string > & GetScalarVariableNames()
Methods to get information about the current variables.
Superclass for algorithms that produce output of the same type as input.
Store vtkAlgorithm input/output information.
perform mathematical operations on data in field data arrays
abstract class to specify dataset behavior
Definition: vtkDataSet.h:62
void SetAttributeTypeToVertexData()
Control which AttributeType the filter operates on (point data or cell data for vtkDataSets).
void SetFunctionParserTypeToFunctionParser()
Set/Get the FunctionParser type that will be used.
void SetAttributeTypeToDefault()
Control which AttributeType the filter operates on (point data or cell data for vtkDataSets).
void SetAttributeTypeToPointData()
Control which AttributeType the filter operates on (point data or cell data for vtkDataSets).
const std::vector< int > & GetSelectedScalarComponents()
Methods to get information about the current variables.
void SetAttributeTypeToCellData()
Control which AttributeType the filter operates on (point data or cell data for vtkDataSets).
int GetNumberOfVectorArrays()
Methods to get information about the current variables.
std::vector< std::string > CoordinateVectorVariableNames
int vtkTypeBool
Definition: vtkABI.h:69
const std::vector< std::string > & GetVectorVariableNames()
Methods to get information about the current variables.
void SetFunctionParserTypeToExprTkFunctionParser()
Set/Get the FunctionParser type that will be used.
std::vector< std::string > VectorVariableNames
a simple class to control print indentation
Definition: vtkIndent.h:39
int FillInputPortInformation(int port, vtkInformation *info) override
Fill the input port information objects for this algorithm.
int GetNumberOfScalarArrays()
Methods to get information about the current variables.
virtual void Modified()
Update the modification time for this object.
std::vector< std::string > CoordinateScalarVariableNames
std::vector< std::string > ScalarVariableNames
FunctionParserTypes FunctionParserType
std::vector< vtkTuple< int, 3 > > SelectedVectorComponents
const std::vector< std::string > & GetVectorArrayNames()
Methods to get information about the current variables.
std::vector< vtkTuple< int, 3 > > SelectedCoordinateVectorComponents
FunctionParserTypes
Enum that includes the types of parsers that can be used.
std::vector< std::string > VectorArrayNames
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
void SetAttributeTypeToEdgeData()
Control which AttributeType the filter operates on (point data or cell data for vtkDataSets).
const std::vector< std::string > & GetScalarArrayNames()
Methods to get information about the current variables.
Store zero or more vtkInformation instances.
vtkTypeBool CoordinateResults
general representation of visualization data
Definition: vtkDataObject.h:65
virtual int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *)
This is called within ProcessRequest when a request asks the algorithm to do its work.
void SetAttributeTypeToRowData()
Control which AttributeType the filter operates on (point data or cell data for vtkDataSets).
std::vector< int > SelectedScalarComponents
static vtkPassInputTypeAlgorithm * New()
std::vector< int > SelectedCoordinateScalarComponents