foamVtkToolsTemplates.C
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-2020 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
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 InNamespace
27  Foam::vtk::Tools
28 
29 \*---------------------------------------------------------------------------*/
30 
31 // OpenFOAM includes
32 #include "error.H"
33 
34 // VTK includes
35 #include "vtkFloatArray.h"
36 #include "vtkCellData.h"
37 #include "vtkPointData.h"
38 #include "vtkSmartPointer.h"
39 
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 
42 template<class Face>
43 vtkSmartPointer<vtkCellArray>
45 {
47 
48  #ifdef VTK_CELL_ARRAY_V2
49 
50  // Offsets
51  // [0, n1, n1+n2, n1+n2+n3... ]
52 
53  const vtkIdType nOffsets(faces.size()+1);
54 
55  auto offsets = vtkSmartPointer<vtkIdTypeArray>::New();
56 
57  vtkIdType nConnect(0);
58  {
59  offsets->SetNumberOfTuples(nOffsets);
60 
61  vtkIdType* iter = offsets->WritePointer(0, nOffsets);
62 
63  // Assign offsets, determine overall connectivity size
64 
65  *iter = 0;
66  for (const auto& f : faces)
67  {
68  nConnect += f.size();
69 
70  *(++iter) = nConnect;
71  }
72  }
73 
74 
75  // Cell connectivity for polygons
76  // [verts..., verts... ]
77 
78  auto connect = vtkSmartPointer<vtkIdTypeArray>::New();
79 
80  {
81  connect->SetNumberOfTuples(nConnect);
82 
83  vtkIdType* iter = connect->WritePointer(0, nConnect);
84 
85  // Fill in the connectivity array
86 
87  for (const auto& f : faces)
88  {
89  for (const label verti : f)
90  {
91  *(iter++) = verti;
92  }
93  }
94  }
95 
96  // Move into a vtkCellArray
97 
98  cells->SetData(offsets, connect);
99 
100  #else
101 
102  // In VTK-8.2.0 and older,
103  // sizes are interwoven (prefixed) in the connectivity
104 
105  // Cell connectivity for polygons
106  // [n1, verts..., n2, verts... ]
107 
108 
109  const vtkIdType nElem(faces.size());
110 
111  // Connectivity size, with prefixed size information
112  vtkIdType nConnect(faces.size());
113  for (const auto& f : faces)
114  {
115  nConnect += f.size();
116  }
117 
118  {
119  cells->GetData()->SetNumberOfTuples(nConnect);
120 
121  vtkIdType* iter = cells->WritePointer(nElem, nConnect);
122 
123  // Fill in the connectivity array, with prefixed size information
124 
125  for (const auto& f : faces)
126  {
127  *(iter++) = f.size();
128 
129  for (const label verti : f)
130  {
131  *(iter++) = verti;
132  }
133  }
134  }
135 
136  #endif
137 
138  return cells;
139 }
140 
141 
142 template<class PatchType>
143 vtkSmartPointer<vtkPoints>
145 {
146  // Local patch points to vtkPoints
147  return Tools::Points(p.localPoints());
148 }
149 
150 
151 template<class PatchType>
152 vtkSmartPointer<vtkCellArray>
154 {
155  return Tools::Faces(p.localFaces());
156 }
157 
158 
159 template<class PatchType>
160 vtkSmartPointer<vtkPolyData>
162 {
163  auto vtkmesh = vtkSmartPointer<vtkPolyData>::New();
164 
165  vtkmesh->SetPoints(points(p));
166  vtkmesh->SetPolys(faces(p));
167 
168  return vtkmesh;
169 }
170 
171 
172 template<class Face>
173 vtkSmartPointer<vtkPolyData>
175 (
176  const UList<point>& pts,
177  const UList<Face>& fcs
178 )
179 {
180  auto vtkmesh = vtkSmartPointer<vtkPolyData>::New();
181 
182  vtkmesh->SetPoints(Tools::Points(pts));
183  vtkmesh->SetPolys(Tools::Faces(fcs));
184 
185  return vtkmesh;
186 }
187 
188 
189 template<class PatchType>
190 vtkSmartPointer<vtkFloatArray>
192 {
194 
195  array->SetNumberOfComponents(3);
196  array->SetNumberOfTuples(p.size());
197 
198  // Unit normals for patch faces.
199  // Cached values if available or loop over faces (avoid triggering cache)
200 
201  vtkIdType faceId = 0;
202 
203  if (p.hasFaceNormals())
204  {
205  for (const vector& n : p.faceNormals())
206  {
207  array->SetTuple(faceId++, n.v_);
208  }
209  }
210  else
211  {
212  for (const auto& f : p)
213  {
214  const vector n(f.unitNormal(p.points()));
215  array->SetTuple(faceId++, n.v_);
216  }
217  }
218 
219  return array;
220 }
221 
222 
223 template<class PatchType>
224 vtkSmartPointer<vtkPoints>
226 {
227  auto vtkpoints = vtkSmartPointer<vtkPoints>::New();
228 
229  vtkpoints->SetNumberOfPoints(p.size());
230 
231  // Use cached values if available or loop over faces
232  // (avoid triggering cache)
233 
234  vtkIdType pointId = 0;
235 
236  if (p.hasFaceCentres())
237  {
238  for (const point& pt : p.faceCentres())
239  {
240  vtkpoints->SetPoint(pointId++, pt.v_);
241  }
242  }
243  else
244  {
245  for (const auto& f : p)
246  {
247  const point pt(f.centre(p.points()));
248  vtkpoints->SetPoint(pointId++, pt.v_);
249  }
250  }
251 
252  return vtkpoints;
253 }
254 
255 
256 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
257 //
258 // Low-Level conversions
259 //
260 
261 template<class Type>
263 (
264  vtkFloatArray* array,
265  const UList<Type>& input,
266  vtkIdType start
267 )
268 {
269  const int nComp(pTraits<Type>::nComponents);
270 
271  if (array->GetNumberOfComponents() != nComp)
272  {
274  << "vtk array '" << array->GetName()
275  << "' has mismatch in number of components for type '"
277  << "' : target array has " << array->GetNumberOfComponents()
278  << " components instead of " << nComp
279  << nl;
280  }
281 
282  const vtkIdType maxSize = array->GetNumberOfTuples();
283  const vtkIdType endPos = start + vtkIdType(input.size());
284 
285  if (!maxSize)
286  {
287  // no-op
288  return 0;
289  }
290  else if (start < 0 || vtkIdType(start) >= maxSize)
291  {
293  << "vtk array '" << array->GetName()
294  << "' copy with out-of-range [0," << long(maxSize) << ")"
295  << " starting at " << long(start)
296  << nl;
297 
298  return 0;
299  }
300  else if (endPos > maxSize)
301  {
303  << "vtk array '" << array->GetName()
304  << "' copy ends out-of-range (" << long(maxSize) << ")"
305  << " using sizing (start,size) = ("
306  << long(start) << "," << input.size() << ")"
307  << nl;
308 
309  return 0;
310  }
311 
312  float scratch[pTraits<Type>::nComponents];
313 
314  for (const Type& val : input)
315  {
316  foamToVtkTuple(scratch, val);
317  array->SetTuple(start++, scratch);
318  }
319 
320  return input.size();
321 }
322 
323 
324 template<class Type>
325 vtkSmartPointer<vtkFloatArray>
327 (
328  const word& name,
329  const label size
330 )
331 {
333 
334  data->SetName(name.c_str());
335  data->SetNumberOfComponents(static_cast<int>(pTraits<Type>::nComponents));
336  data->SetNumberOfTuples(size);
337 
338  // Fill() was not available before VTK-8
339  #if (VTK_MAJOR_VERSION < 8)
340  for (int i = 0; i < data->GetNumberOfComponents(); ++i)
341  {
342  data->FillComponent(i, 0);
343  }
344  #else
345  data->Fill(0);
346  #endif
347 
348  return data;
349 }
350 
351 
352 template<class Type>
353 vtkSmartPointer<vtkFloatArray>
355 (
356  const word& name,
357  const UList<Type>& fld
358 )
359 {
361 
362  data->SetName(name.c_str());
363  data->SetNumberOfComponents(static_cast<int>(pTraits<Type>::nComponents));
364  data->SetNumberOfTuples(fld.size());
365 
367 
368  return data;
369 }
370 
371 
372 // ************************************************************************* //
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::vtk::Tools::Faces
vtkSmartPointer< vtkCellArray > Faces(const UList< Face > &faces)
Convert a list of faces (or triFaces) to vtk polygon cells.
Definition: foamVtkToolsTemplates.C:44
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::vtk::Tools::Patch::mesh
static vtkSmartPointer< vtkPolyData > mesh(const PatchType &p)
Convert patch points/faces to vtkPolyData.
Definition: foamVtkToolsTemplates.C:161
Foam::vtk::Tools::Patch::faceNormals
static vtkSmartPointer< vtkFloatArray > faceNormals(const PatchType &p)
Convert patch face normals to vtkFloatArray.
Definition: foamVtkToolsTemplates.C:191
Foam::vtk::Tools::Patch::faceCentres
static vtkSmartPointer< vtkPoints > faceCentres(const PatchType &p)
Return patch face centres as vtkPoints.
Definition: foamVtkToolsTemplates.C:225
Foam::vtk::Tools::zeroField
vtkSmartPointer< vtkFloatArray > zeroField(const word &name, const label size)
Create named field initialized to zero.
Definition: foamVtkToolsTemplates.C:327
Foam::VectorSpace::v_
Cmpt v_[Ncmpts]
The components of this vector space.
Definition: VectorSpace.H:83
Foam::vtk::Tools::transcribeFloatData
label transcribeFloatData(vtkFloatArray *array, const UList< Type > &input, vtkIdType start=0)
Copy list to pre-allocated vtk array.
n
label n
Definition: TABSMDCalcMethod2.H:31
error.H
faceId
label faceId(-1)
fld
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Definition: gmvOutputLagrangian.H:23
Foam::vtk::Tools::convertFieldToVTK
vtkSmartPointer< vtkFloatArray > convertFieldToVTK(const word &name, const UList< Type > &fld)
Convert field data to a vtkFloatArray.
Definition: foamVtkToolsTemplates.C:355
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::nl
constexpr char nl
Definition: Ostream.H:404
f
labelList f(nPoints)
Foam::Vector< scalar >
Foam::pTraits
A traits class, which is primarily used for primitives.
Definition: pTraits.H:56
Foam::UList< Face >
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::vtk::Tools::Patch::points
static vtkSmartPointer< vtkPoints > points(const PatchType &p)
Return local patch points as vtkPoints.
Definition: foamVtkToolsTemplates.C:144
Foam::input
static Istream & input(Istream &is, IntRange< T > &range)
Definition: IntRanges.C:55
Foam::vtk::Tools::Points
vtkSmartPointer< vtkPoints > Points(const UList< point > &pts)
Return a list of points as vtkPoints.
Definition: foamVtkToolsI.H:59
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
Foam::UList::size
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
Foam::vtk::Tools::foamToVtkTuple
void foamToVtkTuple(float output[], const Type &val)
Copy/transcribe OpenFOAM data types to VTK format.
Definition: foamVtkToolsI.H:214
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:328
Foam::data
Database for solution data, solver performance and other reduced data.
Definition: data.H:55
Foam::vtk::Tools::Patch::faces
static vtkSmartPointer< vtkCellArray > faces(const PatchType &p)
Convert patch faces to vtk polygon cells.
Definition: foamVtkToolsTemplates.C:153