writeFunsTemplates.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 | Copyright (C) 2011-2014 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 "writeFuns.H"
27 #include "interpolatePointToCell.H"
28 
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
30 
31 // Store List in dest
32 template<class Type>
34 (
35  const List<Type>& source,
37 )
38 {
39  forAll(source, i)
40  {
41  insert(source[i], dest);
42  }
43 }
44 
45 
47 //template<class Type>
48 //void Foam::writeFuns::insert
49 //(
50 // const labelList& map,
51 // const List<Type>& source,
52 // DynamicList<doubleScalar>& dest
53 //)
54 //{
55 // forAll(map, i)
56 // {
57 // insert(source[map[i]], dest);
58 // }
59 //}
60 
61 
62 template<class Type>
64 (
65  std::ostream& os,
66  const bool binary,
68  const vtkMesh& vMesh
69 )
70 {
71  const fvMesh& mesh = vMesh.mesh();
72 
73  const labelList& superCells = vMesh.topo().superCells();
74 
75  label nValues = mesh.nCells() + superCells.size();
76 
77  os << vvf.name() << ' ' << pTraits<Type>::nComponents << ' '
78  << nValues << " double" << std::endl;
79 
81 
82  insert(vvf, fField);
83 
84  forAll(superCells, superCellI)
85  {
86  label origCellI = superCells[superCellI];
87 
88  insert(vvf[origCellI], fField);
89  }
90  write(os, binary, fField);
91 }
92 
93 
94 template<class Type>
96 (
97  std::ostream& os,
98  const bool binary,
100  const vtkMesh& vMesh
101 )
102 {
103  const fvMesh& mesh = vMesh.mesh();
104  const vtkTopo& topo = vMesh.topo();
105 
106  const labelList& addPointCellLabels = topo.addPointCellLabels();
107  const label nTotPoints = mesh.nPoints() + addPointCellLabels.size();
108 
109  os << pvf.name() << ' ' << pTraits<Type>::nComponents << ' '
110  << nTotPoints << " double" << std::endl;
111 
113 
114  insert(pvf, fField);
115 
116  forAll(addPointCellLabels, api)
117  {
118  label origCellI = addPointCellLabels[api];
119 
120  insert(interpolatePointToCell(pvf, origCellI), fField);
121  }
122  write(os, binary, fField);
123 }
124 
125 
126 template<class Type>
128 (
129  std::ostream& os,
130  const bool binary,
133  const vtkMesh& vMesh
134 )
135 {
136  const fvMesh& mesh = vMesh.mesh();
137  const vtkTopo& topo = vMesh.topo();
138 
139  const labelList& addPointCellLabels = topo.addPointCellLabels();
140  const label nTotPoints = mesh.nPoints() + addPointCellLabels.size();
141 
142  os << vvf.name() << ' ' << pTraits<Type>::nComponents << ' '
143  << nTotPoints << " double" << std::endl;
144 
146 
147  insert(pvf, fField);
148 
149  forAll(addPointCellLabels, api)
150  {
151  label origCellI = addPointCellLabels[api];
152 
153  insert(vvf[origCellI], fField);
154  }
155  write(os, binary, fField);
156 }
157 
158 
159 template<class Type, template<class> class PatchField, class GeoMesh>
161 (
162  std::ostream& os,
163  const bool binary,
165  const vtkMesh& vMesh
166 )
167 {
168  forAll(flds, i)
169  {
170  write(os, binary, flds[i].dimensionedInternalField(), vMesh);
171  }
172 }
173 
174 
175 template<class Type>
177 (
178  std::ostream& os,
179  const bool binary,
181  const vtkMesh& vMesh
182 )
183 {
184  forAll(flds, i)
185  {
186  write(os, binary, flds[i], vMesh);
187  }
188 }
189 
190 
191 template<class Type>
193 (
194  std::ostream& os,
195  const bool binary,
196  const volPointInterpolation& pInterp,
198  const vtkMesh& vMesh
199 )
200 {
201  forAll(flds, i)
202  {
203  write(os, binary, flds[i], pInterp.interpolate(flds[i])(), vMesh);
204  }
205 }
206 
207 
208 // ************************************************************************* //
Foam::writeFuns::write
static void write(std::ostream &, const bool binary, const GeometricField< Type, pointPatchField, pointMesh > &, const vtkMesh &)
Write pointField on all mesh points. Interpolate to cell centre.
Definition: writeFunsTemplates.C:96
writeFuns.H
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::DynamicList
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:56
Foam::interpolatePointToCell
Type interpolatePointToCell(const GeometricField< Type, pointPatchField, pointMesh > &ptf, const label celli)
Definition: interpolatePointToCell.C:32
dimensionedInternalField
rDeltaT dimensionedInternalField()
interpolatePointToCell.H
Interpolates (averages) the vertex values to the cell center.
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::vtkMesh
Encapsulation of VTK mesh data. Holds mesh or meshsubset and polyhedral-cell decomposition on it.
Definition: vtkMesh.H:52
Foam::vtkMesh::topo
const vtkTopo & topo() const
topology
Definition: vtkMesh.H:110
Foam::label
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
Foam::vtkTopo::superCells
const labelList & superCells() const
Definition: vtkTopo.H:123
Foam::writeFuns::write
static void write(std::ostream &, const bool, DynamicList< floatScalar > &)
Write floats ascii or binary.
Definition: writeFuns.C:107
Foam::PtrList
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:61
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::GeoMesh
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:46
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam::vtkTopo
Polyhedral cell decomposition for VTK.
Definition: vtkTopo.H:51
Foam::volPointInterpolation::interpolate
tmp< GeometricField< Type, pointPatchField, pointMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &) const
Interpolate volField using inverse distance weighting.
Definition: volPointInterpolate.C:502
Foam::vtkMesh::mesh
const fvMesh & mesh() const
Access either mesh or submesh.
Definition: vtkMesh.H:119
Foam::vtkTopo::addPointCellLabels
const labelList & addPointCellLabels() const
Definition: vtkTopo.H:118
Foam::List< Type >
Foam::pTraits
Traits class for primitives.
Definition: pTraits.H:50
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::volPointInterpolation
Interpolate from cell centres to points (vertices) using inverse distance weighting.
Definition: volPointInterpolation.H:56
Foam::writeFuns::insert
static void insert(const point &, DynamicList< floatScalar > &dest)
Append point to given DynamicList.
Definition: writeFuns.C:170
write
Tcoeff write()
insert
timeIndices insert(timeIndex, timeDirs[timeI].value())
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:51