File: vtkCellGridSidesQuery.h

package info (click to toggle)
vtk9 9.5.2%2Bdfsg3-6
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 205,984 kB
  • sloc: cpp: 2,336,570; 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: 181; javascript: 165; objc: 153; tcl: 59
file content (239 lines) | stat: -rw-r--r-- 10,172 bytes parent folder | download
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
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
// SPDX-License-Identifier: BSD-3-Clause
#ifndef vtkCellGridSidesQuery_h
#define vtkCellGridSidesQuery_h

#include "vtkCellGridQuery.h"

#include "vtkStringToken.h" // For API.

#include <functional>
#include <map>
#include <set>
#include <unordered_map>
#include <vector>

VTK_ABI_NAMESPACE_BEGIN

class vtkCellGridSidesCache;
class vtkIdTypeArray;

/**\brief A cell-grid query for enumerating sides of cells.
 *
 * This query runs in 3 passes (see vtkCellGridSidesQuery::PassWork):
 *
 * + In the first pass, responders invoke the AddSides() method on
 *   this query, entries are added to this->Hashes storage indicating
 *   the cells which are bounded by a given shape + connectivity.
 * + In the second pass, responders mark the entries created above and
 *   add entries in this->Sides. This reorganizes the hashes into groups
 *   more amenable to output as side arrays. This pass is called
 *   "Summarization," since not every input side identified will be
 *   output.
 * + In the third and final pass, responders create new cells in
 *   the output cell-grid that correspond to the selected sides of
 *   the input.
 */
class VTKCOMMONDATAMODEL_EXPORT vtkCellGridSidesQuery : public vtkCellGridQuery
{
public:
  static vtkCellGridSidesQuery* New();
  vtkTypeMacro(vtkCellGridSidesQuery, vtkCellGridQuery);
  void PrintSelf(ostream& os, vtkIndent indent) override;

  /// An enum specifying which side(s) each responder should generate.
  enum SideFlags : int
  {
    // Individual bits
    VerticesOfEdges = 0x01,    //!< Input edges should produce endpoint vertices.
    VerticesOfSurfaces = 0x02, //!< Input surfaces should produce corner vertices.
    EdgesOfSurfaces = 0x04,    //!< Input surfaces should produce bounding edges.
    VerticesOfVolumes = 0x08,  //!< Input volumes should produce corner vertices.
    EdgesOfVolumes = 0x10,     //!< Input volumes should produce edges bounding faces.
    SurfacesOfVolumes = 0x20,  //!< Input volumes should produce bounding surfaces.

    // Useful (but not exhaustive) combinations
    SurfacesOfInputs = 0x20,   //!< Produce surfaces (i.e., omit sides of renderable geometry).
    EdgesOfInputs = 0x14,      //!< Produce edges of inputs of higher dimension.
    VerticesOfInputs = 0x0b,   //!< Produce vertices of inputs of higher dimension.
    AllSides = 0x3f,           //!< Produce all sides of inputs that have any sides.
    NextLowestDimension = 0x25 //!< Given an input of dimension D, produce sides of dimension D-1.
  };

  /// An enum specifying the work responders should perform for each pass.
  enum PassWork : int
  {
    /// Responders should call AddSide() on each cell's sides to insert
    /// entries into this->SideCache.
    HashSides = 0,
    /// Responders should claim entries in this->SideCache by transcribing them
    /// to this->Sides and then deleting the entry in this->SideCache to prevent
    /// multiple responders from processing them.
    Summarize = 1,
    /// Responders should insert new side sets into their parent cell-grid
    /// by examining this->Sides.
    GenerateSideSets = 2
  };

  /// An enum specifying the strategy by which input hashes are summarized
  /// into output Sides entries.
  enum SummaryStrategy
  {
    /// The number of hash entries for a given side should be used to
    /// decide whether to include a hash's side in the output.
    /// If a side occurs an odd number of times, it should be included
    /// in the output.
    Winding,
    /// If a hash entry entry exists, a single side should be included
    /// in the output.
    AnyOccurrence,
    /// Hashes for shapes (a) of dimension (d-1) that occur an odd number of
    /// times and (b) of dimension < (d-1) that occur once or more should
    /// be included in the output. (Dimension d is the dimension of the
    /// input shapes, whether they are cells or sides.
    Boundary
  };

  /// Indicate how output should be generated or marked so selection works as expected.
  enum SelectionMode
  {
    Input, //!< Input shapes should be selected when output sides are picked.
    Output //!< Output sides should be selected when they are picked.
  };

  /// A structure created by the GetSideSetArrays() method for responders to use.
  struct SideSetArray
  {
    /// The type of parent cell which created the sides.
    vtkStringToken CellType;
    /// The shape of all the sides in the \a Sides array.
    vtkStringToken SideShape;
    /// An array of tuples of (cell-id, side-id) specifying sides.
    vtkSmartPointer<vtkIdTypeArray> Sides;
  };

  /// Set/get whether renderable cells should be included in the output
  /// or the output should strictly contain sides of cells.
  ///
  /// A cell is renderable if it is of dimension 2 or less (i.e., surfaces,
  /// edges, and vertices are all renderable; volumetric cells are not).
  ///
  /// The default is false.
  vtkSetMacro(PreserveRenderableInputs, vtkTypeBool);
  vtkGetMacro(PreserveRenderableInputs, vtkTypeBool);
  vtkBooleanMacro(PreserveRenderableInputs, vtkTypeBool);

  /// Set/get whether to omit computation of sides for renderable cells.
  ///
  /// A cell is renderable if it is of dimension 2 or less (i.e., surfaces,
  /// edges, and vertices are all renderable; volumetric cells are not).
  /// This setting, when used in combination with PreserveRenderableInputs,
  /// allows the filter to behave similar to vtkPolyData surface extraction
  /// filters; volumetric cells will have sides computed but others will be
  /// passed through from the input unaltered.
  ///
  /// The default is false.
  vtkSetMacro(OmitSidesForRenderableInputs, vtkTypeBool);
  vtkGetMacro(OmitSidesForRenderableInputs, vtkTypeBool);
  vtkBooleanMacro(OmitSidesForRenderableInputs, vtkTypeBool);

  /// Set/get which sides to generate given input cells/sides.
  ///
  /// OutputDimensionControl is a bit-vector taking values from the SideFlags
  /// enumeration. It determines which sides of the input should be generated.
  /// The default is SideFlags::SurfacesOfInputs, which will only emit surfaces
  /// of volumetric cells.
  vtkSetMacro(OutputDimensionControl, int);
  vtkGetMacro(OutputDimensionControl, int);
  vtkBooleanMacro(OutputDimensionControl, int);

  /// Set/get the strategy responders should use to generate entries in
  /// Sides from entries in SideCache.
  ///
  /// The default is BoundaryStrategy.
  vtkSetEnumMacro(Strategy, SummaryStrategy);
  vtkGetEnumMacro(Strategy, SummaryStrategy);
  void SetStrategyToWinding() { this->SetStrategy(SummaryStrategy::Winding); }
  void SetStrategyToAnyOccurrence() { this->SetStrategy(SummaryStrategy::AnyOccurrence); }
  void SetStrategyToBoundary() { this->SetStrategy(SummaryStrategy::Boundary); }
  /// This method exists for ParaView to set the strategy.
  virtual void SetStrategy(int strategy)
  {
    this->SetStrategy(static_cast<SummaryStrategy>(strategy));
  }

  /// Set/get whether the extracted sides should be marked as selectable
  /// or whether their originating data should be selectable.
  ///
  /// Responders should use this to *either*:
  /// (a) mark the output to indicate what shapes should be selected upon being picked; or
  /// (b) output different shapes so that picking implicitly results in the proper shape
  ///     being picked.
  ///
  /// The default SelectionMode::Input indicates the input data should be selected.
  /// Other values indicate the generated output sides should be selected.
  vtkSetEnumMacro(SelectionType, SelectionMode);
  vtkGetEnumMacro(SelectionType, SelectionMode);
  /// This method exists for ParaView to set the selection mode.
  virtual void SetSelectionType(int selnType)
  {
    this->SetSelectionType(static_cast<SelectionMode>(selnType));
  }

  bool Initialize() override;
  void StartPass() override;
  bool IsAnotherPassRequired() override;
  bool Finalize() override;

  std::unordered_map<vtkStringToken,
    std::unordered_map<vtkStringToken, std::unordered_map<vtkIdType, std::set<int>>>>&
  GetSides()
  {
    return this->Sides;
  }

  /// Return arrays of cell+side IDs for the given \a cellType.
  std::vector<SideSetArray> GetSideSetArrays(vtkStringToken cellType);

  /// Return a string-token with the given selection mode or vice-versa.
  static vtkStringToken SelectionModeToLabel(SelectionMode mode);
  static SelectionMode SelectionModeFromLabel(vtkStringToken token);

  /// Return a string-token with the given summarization strategy or vice-versa.
  static vtkStringToken SummaryStrategyToLabel(SummaryStrategy strategy);
  static SummaryStrategy SummaryStrategyFromLabel(vtkStringToken token);

  /// Set/get cached hashtable of sides.
  ///
  /// The idea is that vtkCellGridSidesCache is generic enough to accommodate
  /// a wide variety of cell types and that many of them will be capable of
  /// having sides that are conformal to cells of different types that may
  /// reside in the same cell-grid. Filters may own this cache or they may
  /// attach it to a collection of cell-grid objects that participate by inserting
  /// their cells' sides into the cache. (For example, all the cell-grids within
  /// a partitioned dataset collection may wish to insert sides in the same cache.)
  vtkGetObjectMacro(SideCache, vtkCellGridSidesCache);
  virtual void SetSideCache(vtkCellGridSidesCache* cache);

protected:
  vtkCellGridSidesQuery() = default;
  ~vtkCellGridSidesQuery() override;

  vtkTypeBool PreserveRenderableInputs{ false };
  vtkTypeBool OmitSidesForRenderableInputs{ false };
  int OutputDimensionControl{ SideFlags::SurfacesOfInputs };
  SelectionMode SelectionType{ SelectionMode::Input };
  SummaryStrategy Strategy{ SummaryStrategy::Boundary };
  vtkCellGridSidesCache* SideCache{ nullptr };
  bool TemporarySideCache{ false };
  std::unordered_map<vtkStringToken,
    std::unordered_map<vtkStringToken, std::unordered_map<vtkIdType, std::set<int>>>>
    Sides;

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

VTK_ABI_NAMESPACE_END
#endif // vtkCellGridSidesQuery_h