Loading...
Searching...
No Matches
foamVtkToolsI.H
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2017-2025 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
29
31(
32 vtkUnsignedCharArray* array,
33 const label size
34)
35{
36 array->SetNumberOfComponents(1);
37 array->SetNumberOfTuples(size);
38
39 return UList<uint8_t>(array->WritePointer(0, size), size);
40}
41
42
44(
45 vtkIdTypeArray* array,
46 const label size
47)
48{
49 array->SetNumberOfComponents(1);
50 array->SetNumberOfTuples(size);
51
52 return UList<vtkIdType>(array->WritePointer(0, size), size);
53}
54
55
56inline vtkSmartPointer<vtkPoints>
58(
59 vtkIdType count
60)
61{
62 auto vtkpoints = vtkSmartPointer<vtkPoints>::New();
63 // if (narrowing || (sizeof(float) == sizeof(Foam::scalar)))
64 // {
65 // vtkpoints->SetDataTypeToFloat();
66 // }
67 // else
68 // {
69 // vtkpoints->SetDataTypeToDouble();
70 // }
71 vtkpoints->SetNumberOfPoints(count);
72
73 return vtkpoints;
74}
75
76
77inline vtkSmartPointer<vtkPoints>
79(
80 const UList<point>& pts
81)
82{
83 auto vtkpoints = vtk::Tools::NewPoints(pts.size());
84
85 vtkIdType pointId = 0;
86 for (const auto& p : pts)
87 {
88 vtkpoints->SetPoint(pointId++, p.cdata());
89 }
90
91 return vtkpoints;
92}
93
94
95inline vtkSmartPointer<vtkPoints>
97(
98 const UList<point>& pts,
99 const labelUList& addr
100)
101{
102 auto vtkpoints = vtk::Tools::NewPoints(addr.size());
103
104 vtkIdType pointId = 0;
105 for (const label pointi : addr)
106 {
107 const auto& p = pts[pointi];
108 vtkpoints->SetPoint(pointId++, p.cdata());
110
111 return vtkpoints;
112}
113
114
115inline vtkSmartPointer<vtkPolyData>
117(
118 const UList<point>& pts
119)
120{
121 auto vtkmesh = vtkSmartPointer<vtkPolyData>::New();
122
123 vtkmesh->SetPoints(Tools::Points(pts));
124 vtkmesh->SetVerts(Tools::identityVertices(pts.size()));
125
126 return vtkmesh;
127}
128
129
130inline vtkSmartPointer<vtkPolyData>
132(
133 const UList<point>& pts,
134 const labelUList& addr
135)
136{
137 auto vtkmesh = vtkSmartPointer<vtkPolyData>::New();
138
139 vtkmesh->SetPoints(Tools::Points(pts, addr));
140 vtkmesh->SetVerts(Tools::identityVertices(addr.size()));
141
142 return vtkmesh;
143}
144
145
146inline Foam::scalarMinMax Foam::vtk::Tools::rangeOf(vtkDataArray* data)
147{
148 double range[2]{GREAT, -GREAT};
149
150 if (data)
151 {
152 if (data->GetNumberOfComponents() == 1)
153 {
154 data->GetRange(range);
155 }
156 else
157 {
158 // Mag
159 data->GetRange(range, -1);
160 }
162
163 return scalarMinMax(range[0], range[1]);
164}
165
166
167inline vtkSmartPointer<vtkCellArray>
169(
170 const vtkIdType size
171)
172{
173 auto cells = vtkSmartPointer<vtkCellArray>::New();
174
175 #ifdef VTK_CELL_ARRAY_V2
176
177 // Offsets
178 // [0, n1, n1+n2, n1+n2+n3... ]
179
180 auto offsets = vtkSmartPointer<vtkIdTypeArray>::New();
181
182 {
183 const vtkIdType nOffsets(size+1);
184
185 offsets->SetNumberOfTuples(nOffsets);
186
187 vtkIdType* iter = offsets->WritePointer(0, nOffsets);
188
189 std::iota(iter, (iter + nOffsets), vtkIdType(0));
190 }
191
192
193 auto connect = vtkSmartPointer<vtkIdTypeArray>::New();
194
195 // Connectivity
196 {
197 const vtkIdType nConnect(size);
198
199 connect->SetNumberOfTuples(nConnect);
200
201 vtkIdType* iter = connect->WritePointer(0, nConnect);
202
203 std::iota(iter, (iter + nConnect), vtkIdType(0));
204 }
205
206 // Move into a vtkCellArray
207
208 cells->SetData(offsets, connect);
209
210 #else
211
212 // In VTK-8.2.0 and older,
213 // sizes are interwoven (prefixed) in the connectivity
214
215 // Connectivity size, with prefixed size information
216
217 // Cell connectivity for vertex
218 // [size, ids.., size, ids...] -> therefore [1, id, 1, id, ...]
219
220 const vtkIdType nElem(size);
221 const vtkIdType nConnect(2*size);
222
223 {
224 cells->GetData()->SetNumberOfTuples(nConnect);
225
226 vtkIdType* iter = cells->WritePointer(nElem, nConnect);
227
228 // Fill in the connectivity array, with prefixed size information
229
230 for (vtkIdType id = 0; id < nElem; ++id)
231 {
232 *(iter++) = 1;
233 *(iter++) = id;
234 }
235 }
236
237 #endif
238
239 return cells;
240}
241
242
243// ************************************************************************* //
scalar range
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition UList.H:89
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
volScalarField & p
const cellShapeList & cells
UList< uint8_t > asUList(vtkUnsignedCharArray *array, const label size)
Wrap vtkUnsignedCharArray as a UList.
scalarMinMax rangeOf(vtkDataArray *data)
Min/Max of scalar, or mag() of non-scalars. Includes nullptr check.
vtkSmartPointer< vtkPoints > Points(const UList< point > &pts)
Return a list of points as vtkPoints.
vtkSmartPointer< vtkCellArray > identityVertices(const vtkIdType size)
An identity list of VTK_VERTEX.
vtkSmartPointer< vtkPoints > NewPoints(vtkIdType count)
Allocate vtkPoints.
vtkSmartPointer< vtkPolyData > Vertices(const UList< point > &pts)
Return vtkPolyData of vertices for each point.
MinMax< scalar > scalarMinMax
A scalar min/max range.
Definition MinMax.H:97
UList< label > labelUList
A UList of labels.
Definition UList.H:75
const pointField & pts