vtkbone
vtkboneNodeSetsByGeometry.h
Go to the documentation of this file.
1/*=========================================================================
2
3 Copyright 2010-2016, Numerics88 Solutions Ltd.
4 http://www.numerics88.com/
5
6 Copyright (c) Eric Nodwell and Steven K. Boyd
7 See Copyright.txt for details.
8
9 This software is distributed WITHOUT ANY WARRANTY; without even
10 the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
11 PURPOSE. See the above copyright notice for more information.
12=========================================================================*/
13
29#ifndef __vtkboneNodeSetsByGeometry_h
30#define __vtkboneNodeSetsByGeometry_h
31
32#include "vtkObject.h"
33#include "vtkboneWin32Header.h"
34
35// Forward declarations
38class vtkDataSet;
39class vtkIdTypeArray;
40class vtkPolyData;
41class vtkIdList;
42class vtkPoints;
43
44class VTKBONE_EXPORT vtkboneNodeSetsByGeometry : public vtkObject
45{
46public:
49 void PrintSelf(ostream& os, vtkIndent indent) override;
50
52
55 double bounds[6],
56 int specificMaterial = -1);
58
60
64 vtkIdTypeArray *input_ids,
65 vtkDataSet *data,
66 int targetCellScalar);
68
70
76 static int FindNodesOnPlane(
77 int axis,
78 float val,
79 vtkIdTypeArray *ids,
81 int specificMaterial = -1);
83
85
86 static int AddNodesOnPlane(
87 int axis,
88 float val,
89 const char* name,
91 int specificMaterial = -1);
93
95
97 int axis,
98 float val,
99 const char* name,
101 int specificMaterial = -1);
103
105
110 int axis1,
111 float val1,
112 int axis2,
113 float val2,
114 vtkIdTypeArray* ids,
116 int specificMaterial = -1);
118
120
123 int axis1,
124 float val1,
125 int axis2,
126 float val2,
127 const char* name,
129 int specificMaterial = -1);
131
133
136 int axis1,
137 float val1,
138 int axis2,
139 float val2,
140 const char* name,
142 int specificMaterial = -1);
144
146
151 int axisA,
152 float valA,
153 int axisB,
154 float valB,
155 int axisC,
156 float valC,
157 vtkIdTypeArray *ids,
159 int specificMaterial = -1);
161
163
166 int axisA,
167 float valA,
168 int axisB,
169 float valB,
170 int axisC,
171 float valC,
172 const char* name,
174 int specificMaterial = -1);
176
178
181 int axisA,
182 float valA,
183 int axisB,
184 float valB,
185 int axisC,
186 float valC,
187 const char* name,
189 int specificMaterial = -1);
191
193
201 vtkIdTypeArray *visibleNodesIds,
203 double normalVector[3],
204 int specificMaterial = -1);
206
208
210 const char* name,
212 double normalVector[3],
213 int specificMaterial = -1);
215
217
219 const char* name,
221 double normalVector[3],
222 int specificMaterial = -1);
224
225protected:
228
229private:
230 vtkboneNodeSetsByGeometry(const vtkboneNodeSetsByGeometry&); // Not implemented
231 void operator=(const vtkboneNodeSetsByGeometry&); // Not implemented
232};
233
234#endif
void operator=(const vtkObjectBase &)
data model for finite element meshes
various algorithms to select nodes sets from a mesh by geometry.
static int AddNodesIntersectingTwoPlanes(int axis1, float val1, int axis2, float val2, const char *name, vtkboneFiniteElementModel *model, int specificMaterial=-1)
static int AddNodesIntersectingThreePlanes(int axisA, float valA, int axisB, float valB, int axisC, float valC, const char *name, vtkboneFiniteElementModel *model, int specificMaterial=-1)
static int AddNodesAndElementsIntersectingTwoPlanes(int axis1, float val1, int axis2, float val2, const char *name, vtkboneFiniteElementModel *model, int specificMaterial=-1)
static void DetermineMaterialBounds(vtkUnstructuredGrid *geometry, double bounds[6], int specificMaterial=-1)
static int AddNodesAndElementsOnVisibleSurface(const char *name, vtkboneFiniteElementModel *model, double normalVector[3], int specificMaterial=-1)
static int FindNodesOnVisibleSurface(vtkIdTypeArray *visibleNodesIds, vtkUnstructuredGrid *ug, double normalVector[3], int specificMaterial=-1)
static int AddNodesAndElementsIntersectingThreePlanes(int axisA, float valA, int axisB, float valB, int axisC, float valC, const char *name, vtkboneFiniteElementModel *model, int specificMaterial=-1)
static int FindNodesOnPlane(int axis, float val, vtkIdTypeArray *ids, vtkUnstructuredGrid *ug, int specificMaterial=-1)
static int FindNodesIntersectingTwoPlanes(int axis1, float val1, int axis2, float val2, vtkIdTypeArray *ids, vtkUnstructuredGrid *ug, int specificMaterial=-1)
void PrintSelf(ostream &os, vtkIndent indent) override
static int AddNodesOnPlane(int axis, float val, const char *name, vtkboneFiniteElementModel *model, int specificMaterial=-1)
static int FilterPointListByCellScalar(vtkIdTypeArray *output_ids, vtkIdTypeArray *input_ids, vtkDataSet *data, int targetCellScalar)
static int FindNodesIntersectingThreePlanes(int axisA, float valA, int axisB, float valB, int axisC, float valC, vtkIdTypeArray *ids, vtkUnstructuredGrid *ug, int specificMaterial=-1)
static int AddNodesOnVisibleSurface(const char *name, vtkboneFiniteElementModel *model, double normalVector[3], int specificMaterial=-1)
static vtkboneNodeSetsByGeometry * New()
static int AddNodesAndElementsOnPlane(int axis, float val, const char *name, vtkboneFiniteElementModel *model, int specificMaterial=-1)