File: vtkHyperTreeGridGradient.h

package info (click to toggle)
vtk9 9.5.2%2Bdfsg3-4
  • links: PTS, VCS
  • area: main
  • in suites: sid
  • size: 205,916 kB
  • sloc: cpp: 2,336,565; ansic: 327,116; python: 111,200; yacc: 4,104; java: 3,977; sh: 3,032; xml: 2,771; perl: 2,189; lex: 1,787; makefile: 178; javascript: 165; objc: 153; tcl: 59
file content (213 lines) | stat: -rw-r--r-- 5,634 bytes parent folder | download | duplicates (4)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
// SPDX-License-Identifier: BSD-3-Clause
/**
 * @class   vtkHyperTreeGridGradient
 * @brief   Compute the gradient of a scalar field
 * on a Hyper Tree Grid.
 *
 * This filter compute the gradient of a given cell scalars array on a
 * Hyper Tree Grid. This result in a new array attached to the original input.
 * This filter does not support masks.
 * In practice the mask is ignored during the processing of this filters and re-attached to the
 * output, leading to masked cell being taken into account for the gradient computation of visible
 * cells. This leads to the gradient being influenced by masked cells. This should only impact cells
 * on the boundary, where gradient is already ill-defined.
 *
 * @sa
 * vtkHyperTreeGrid vtkHyperTreeGridAlgorithm vtkGradientFilter
 *
 * @par Thanks:
 * This class was created by Charles Gueunet, 2022
 * This work was supported by Commissariat a l'Energie Atomique
 * CEA, DAM, DIF, F-91297 Arpajon, France.
 */

#ifndef vtkHyperTreeGridGradient_h
#define vtkHyperTreeGridGradient_h

#include "vtkFiltersHyperTreeModule.h" // For export macro

#include "vtkHyperTreeGridAlgorithm.h"
#include "vtkNew.h"          // for internal fields
#include "vtkSmartPointer.h" // for internal fields

#include <cstring> // for strdup, to initialize char*

VTK_ABI_NAMESPACE_BEGIN

class vtkHyperTreeGridNonOrientedCursor;
class vtkBitArray;
class vtkDoubleArray;
class vtkUnsignedCharArray;

class VTKFILTERSHYPERTREE_EXPORT vtkHyperTreeGridGradient : public vtkHyperTreeGridAlgorithm
{
public:
  enum ComputeMode
  {
    UNLIMITED = 0,
    UNSTRUCTURED
  };

  static vtkHyperTreeGridGradient* New();
  vtkTypeMacro(vtkHyperTreeGridGradient, vtkHyperTreeGridAlgorithm);
  void PrintSelf(ostream& os, vtkIndent indent) override;

  ///@{
  /**
   * Enable / disable gradient computation.
   * default is On.
   */
  vtkSetMacro(ComputeGradient, bool);
  vtkGetMacro(ComputeGradient, bool);
  vtkBooleanMacro(ComputeGradient, bool);
  ///@}

  ///@{
  /**
   * Set/Get the name of gradient vector array.
   */
  vtkSetStringMacro(GradientArrayName);
  vtkGetStringMacro(GradientArrayName);
  ///@}

  ///@{
  /**
   * Set/Get the gradient computation method to use:
   * * Unlimited: virtually reffine neighbors
   * * Unstructured: compute gradient like in UG
   * Dfault is UNLIMITED
   */
  vtkSetClampMacro(Mode, int, UNLIMITED, UNSTRUCTURED);
  vtkGetMacro(Mode, int);
  ///@}

  ///@{
  /**
   * Do we apply ratio in unlimited mode for the gradient computation ?
   * No effect in Unstructured mode
   * Default is false (intensive computation)
   */
  vtkSetMacro(ExtensiveComputation, bool);
  vtkGetMacro(ExtensiveComputation, bool);
  vtkBooleanMacro(ExtensiveComputation, bool);
  ///@}

  ///@{
  /**
   * Enable / disable divergence computation.
   * default is Off.
   */
  vtkSetMacro(ComputeDivergence, bool);
  vtkGetMacro(ComputeDivergence, bool);
  vtkBooleanMacro(ComputeDivergence, bool);
  ///@}

  ///@{
  /**
   * Set/Get the name of divergence vector array.
   */
  vtkSetStringMacro(DivergenceArrayName);
  vtkGetStringMacro(DivergenceArrayName);
  ///@}

  ///@{
  /**
   * Enable / disable vorticity computation.
   * default is Off.
   */
  vtkSetMacro(ComputeVorticity, bool);
  vtkGetMacro(ComputeVorticity, bool);
  vtkBooleanMacro(ComputeVorticity, bool);
  ///@}

  ///@{
  /**
   * Set/Get the name of vorticity array.
   */
  vtkSetStringMacro(VorticityArrayName);
  vtkGetStringMacro(VorticityArrayName);
  ///@}

  ///@{
  /**
   * Enable / disable q-criterion computation.
   * default is Off.
   */
  vtkSetMacro(ComputeQCriterion, bool);
  vtkGetMacro(ComputeQCriterion, bool);
  vtkBooleanMacro(ComputeQCriterion, bool);
  ///@}

  ///@{
  /**
   * Set/Get the name of q-criterion array.
   */
  vtkSetStringMacro(QCriterionArrayName);
  vtkGetStringMacro(QCriterionArrayName);
  ///@}

protected:
  vtkHyperTreeGridGradient();
  ~vtkHyperTreeGridGradient() override;

  /**
   * Main routine to generate gradient of hyper tree grid.
   */
  int ProcessTrees(vtkHyperTreeGrid*, vtkDataObject*) override;

  /**
   * Recursively descend into tree down to leaves to compute gradient
   * Uses a heavy supercursor
   */
  template <class Cursor, class GradWorker>
  void RecursivelyProcessGradientTree(Cursor*, GradWorker&);

  /**
   * Compute Vorticity, Divergence and QCriterion upon request,
   * from the Gradient cell array
   */
  template <class FieldsWorker>
  void ProcessFields(FieldsWorker&);

  // Fields
  // ------

  // Gradient
  bool ComputeGradient = true;
  vtkNew<vtkDoubleArray> OutGradArray;
  char* GradientArrayName = strdup("Gradient");
  int Mode = ComputeMode::UNLIMITED;
  bool ExtensiveComputation = false;

  // Divergence
  bool ComputeDivergence = false;
  vtkNew<vtkDoubleArray> OutDivArray;
  char* DivergenceArrayName = strdup("Divergence");

  // Vorticity
  bool ComputeVorticity = false;
  vtkNew<vtkDoubleArray> OutVortArray;
  char* VorticityArrayName = strdup("Vorticity");

  // QCriterion
  bool ComputeQCriterion = false;
  vtkNew<vtkDoubleArray> OutQCritArray;
  char* QCriterionArrayName = strdup("QCriterion");

  /**
   * Keep track of selected input scalars / vectors
   */
  vtkSmartPointer<vtkDataArray> InArray;

  // shortcut to HTG fields
  vtkBitArray* InMask = nullptr;
  vtkUnsignedCharArray* InGhostArray = nullptr;

private:
  vtkHyperTreeGridGradient(const vtkHyperTreeGridGradient&) = delete;
  void operator=(const vtkHyperTreeGridGradient&) = delete;
};

VTK_ABI_NAMESPACE_END
#endif // vtkHyperTreeGridGradient_h