vtkWriteFieldOpsTemplates.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2020 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "vtkWriteFieldOps.H"
27 #include "interpolatePointToCell.H"
28 
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
30 
31 template<class Type>
33 (
34  std::ostream& os,
35  const bool binary,
36  const DimensionedField<Type, volMesh>& df,
37  const vtkMesh& vMesh
38 )
39 {
40  const fvMesh& mesh = vMesh.mesh();
41 
42  const labelList& superCells = vMesh.topo().superCells();
43 
44  label nValues = mesh.nCells() + superCells.size();
45 
46  os << df.name() << ' ' << pTraits<Type>::nComponents << ' '
47  << nValues << " float" << std::endl;
48 
49  DynamicList<floatScalar> fField(pTraits<Type>::nComponents*nValues);
50 
51  insert(df, fField);
52 
53  forAll(superCells, superCelli)
54  {
55  label origCelli = superCells[superCelli];
56 
57  insert(df[origCelli], fField);
58  }
59  write(os, binary, fField);
60 }
61 
62 
63 template<class Type>
65 (
66  std::ostream& os,
67  const bool binary,
68  const GeometricField<Type, pointPatchField, pointMesh>& pvf,
69  const vtkMesh& vMesh
70 )
71 {
72  const fvMesh& mesh = vMesh.mesh();
73  const vtkTopo& topo = vMesh.topo();
74 
75  const labelList& addPointCellLabels = topo.addPointCellLabels();
76  const label nTotPoints = mesh.nPoints() + addPointCellLabels.size();
77 
78  os << pvf.name() << ' ' << pTraits<Type>::nComponents << ' '
79  << nTotPoints << " float" << std::endl;
80 
81  DynamicList<floatScalar> fField(pTraits<Type>::nComponents*nTotPoints);
82 
83  insert(pvf, fField);
84 
85  forAll(addPointCellLabels, api)
86  {
87  label origCelli = addPointCellLabels[api];
88 
89  insert(interpolatePointToCell(pvf, origCelli), fField);
90  }
91  write(os, binary, fField);
92 }
93 
94 
95 template<class Type>
97 (
98  std::ostream& os,
99  const bool binary,
100  const GeometricField<Type, fvPatchField, volMesh>& vvf,
101  const GeometricField<Type, pointPatchField, pointMesh>& pvf,
102  const vtkMesh& vMesh
103 )
104 {
105  const fvMesh& mesh = vMesh.mesh();
106  const vtkTopo& topo = vMesh.topo();
107 
108  const labelList& addPointCellLabels = topo.addPointCellLabels();
109  const label nTotPoints = mesh.nPoints() + addPointCellLabels.size();
110 
111  os << vvf.name() << ' ' << pTraits<Type>::nComponents << ' '
112  << nTotPoints << " float" << std::endl;
113 
114  DynamicList<floatScalar> fField(pTraits<Type>::nComponents*nTotPoints);
115 
116  insert(pvf, fField);
117 
118  forAll(addPointCellLabels, api)
119  {
120  label origCelli = addPointCellLabels[api];
121 
122  insert(vvf[origCelli], fField);
123  }
124  write(os, binary, fField);
125 }
126 
127 
128 template<class Type, template<class> class PatchField, class GeoMesh>
130 (
131  std::ostream& os,
132  const bool binary,
133  const PtrList<GeometricField<Type, PatchField, GeoMesh>>& flds,
134  const vtkMesh& vMesh
135 )
136 {
137  forAll(flds, i)
138  {
139  write(os, binary, flds[i], vMesh);
140  }
141 }
142 
143 
144 template<class Type>
146 (
147  std::ostream& os,
148  const bool binary,
149  const volPointInterpolation& pInterp,
150  const PtrList<GeometricField<Type, fvPatchField, volMesh>>& flds,
151  const vtkMesh& vMesh
152 )
153 {
154  forAll(flds, i)
155  {
156  write(os, binary, flds[i], pInterp.interpolate(flds[i])(), vMesh);
157  }
158 }
159 
160 
161 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
void insert(const scalar, DynamicList< floatScalar > &)
Append scalar to given DynamicList.
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
List< label > labelList
A List of labels.
Definition: labelList.H:56
Type interpolatePointToCell(const GeometricField< Type, pointPatchField, pointMesh > &ptf, const label celli)
Interpolates (averages) the vertex values to the cell center.