writeVolFields.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) 2018-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
28 
29 Description
30  Read volume fields from disk
31  and write with vtk::internalWriter and vtk::patchWriter
32 
33 SourceFiles
34  writeVolFields.H
35 
36 \*---------------------------------------------------------------------------*/
37 
38 #ifndef writeVolFields_H
39 #define writeVolFields_H
40 
41 #include "readFields.H"
42 #include "foamVtkInternalWriter.H"
43 #include "foamVtkPatchWriter.H"
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 namespace Foam
48 {
49 
50 template<class GeoField>
51 bool writeVolField
52 (
55  const tmp<GeoField>& tfield
56 )
57 {
58  if (!tfield.valid())
59  {
60  return false;
61  }
62 
63  const auto& field = tfield();
64 
65  // Internal
66  if (internalWriter)
67  {
68  internalWriter->write(field);
69  }
70 
71  // Boundary
73  {
75  }
76 
77  tfield.clear();
78  return true;
79 }
80 
81 
82 template<class GeoField>
83 bool writeVolField
84 (
87 
90 
91  const tmp<GeoField>& tfield
92 )
93 {
94  if (!tfield.valid())
95  {
96  return false;
97  }
98 
99  const auto& field = tfield();
100 
101  // Internal
102  if (internalWriter && pInterp)
103  {
104  internalWriter->write(field, *pInterp);
105  }
106 
107  // Boundary
108  label writeri = 0;
110  {
111  if (writeri < patchInterps.size() && patchInterps.set(writeri))
112  {
113  writer.write(field, patchInterps[writeri]);
114  }
115  ++writeri;
116  }
117 
118  tfield.clear();
119  return true;
120 }
121 
122 
123 template<class GeoField>
124 label writeVolFields
125 (
128  const fvMeshSubsetProxy& proxy,
129  const IOobjectList& objects,
130  const bool syncPar
131 )
132 {
133  label count = 0;
134 
135  for (const word& fieldName : objects.sortedNames<GeoField>())
136  {
137  if
138  (
139  writeVolField<GeoField>
140  (
142  patchWriters,
143  getField<GeoField>(proxy, objects, fieldName, syncPar)
144  )
145  )
146  {
147  ++count;
148  }
149  }
150 
151  return count;
152 }
153 
154 
155 template<class GeoField>
156 label writeVolFields
157 (
160 
163 
164  const fvMeshSubsetProxy& proxy,
165  const IOobjectList& objects,
166  const bool syncPar
167 )
168 {
169  label count = 0;
170 
171  for (const word& fieldName : objects.sortedNames<GeoField>())
172  {
173  if
174  (
175  writeVolField<GeoField>
176  (
179  getField<GeoField>(proxy, objects, fieldName, syncPar)
180  )
181  )
182  {
183  ++count;
184  }
185  }
186 
187  return count;
188 }
189 
190 
191 label writeAllVolFields
192 (
195 
196  const fvMeshSubsetProxy& proxy,
197  const IOobjectList& objects,
198  const bool syncPar
199 )
200 {
201  #undef foamToVtk_WRITE_FIELD
202  #define foamToVtk_WRITE_FIELD(FieldType) \
203  writeVolFields<FieldType> \
204  ( \
205  internalWriter, \
206  patchWriters, \
207  proxy, \
208  objects, \
209  syncPar \
210  )
211 
212  label count = 0;
218 
219  #undef foamToVTK_WRITE_FIELD
220  return count;
221 }
222 
223 
224 label writeAllVolFields
225 (
228 
231 
232  const fvMeshSubsetProxy& proxy,
233  const IOobjectList& objects,
234  const bool syncPar
235 )
236 {
237  #undef foamToVtk_WRITE_FIELD
238  #define foamToVtk_WRITE_FIELD(FieldType) \
239  writeVolFields<FieldType> \
240  ( \
241  internalWriter, pInterp, \
242  patchWriters, patchInterps, \
243  proxy, \
244  objects, \
245  syncPar \
246  )
247 
248 
249  label count = 0;
255 
256  #undef foamToVtk_WRITE_FIELD
257  return count;
258 }
259 
260 
261 } // End namespace Foam
262 
263 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
264 
265 #endif
266 
267 // ************************************************************************* //
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::tmp::clear
void clear() const noexcept
Definition: tmpI.H:287
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
internalWriter
autoPtr< vtk::internalWriter > internalWriter
Definition: convertProcessorPatches.H:60
Foam::writer::write
virtual void write(const coordSet &, const wordList &, const List< const Field< Type > * > &, Ostream &) const =0
General entry point for writing.
Foam::IOobjectList::sortedNames
wordList sortedNames() const
The sorted names of the IOobjects.
Definition: IOobjectList.C:345
patchWriters
PtrList< vtk::patchWriter > patchWriters
Definition: convertProcessorPatches.H:63
Foam::writeVolField
bool writeVolField(ensightCase &ensCase, const ensightMesh &ensMesh, const tmp< GeometricField< Type, fvPatchField, volMesh >> &tfield, const bool nearCellValue=false)
Definition: writeVolFields.H:42
Foam::writeAllVolFields
label writeAllVolFields(ensightCase &ensCase, const ensightMesh &ensMesh, const IOobjectList &objects, const bool nearCellValue=false)
Definition: writeVolFields.H:126
Foam::UPtrList
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition: UPtrList.H:62
field
rDeltaTY field()
Foam::writer
Base class for graphics format writing. Entry points are.
Definition: writer.H:81
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::IOobjectList
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:55
Foam::writeVolFields
label writeVolFields(ensightCase &ensCase, const ensightMesh &ensMesh, const IOobjectList &objects, const bool nearCellValue=false)
Definition: writeVolFields.H:90
foamVtkPatchWriter.H
Foam::vtk::patchWriter
Write OpenFOAM patches and patch fields in VTP or legacy vtk format.
Definition: foamVtkPatchWriter.H:68
readFields.H
Helper routines for reading a field or fields, optionally with a mesh subset (using fvMeshSubsetProxy...
Foam::fvMeshSubsetProxy
Simple proxy for holding a mesh, or mesh-subset. The subMeshes are currently limited to cellSet or ce...
Definition: fvMeshSubsetProxy.H:55
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::PrimitivePatchInterpolation
Interpolation class within a primitive patch. Allows interpolation from points to faces and vice vers...
Definition: PrimitivePatchInterpolation.H:53
Foam::BitOps::count
unsigned int count(const UList< bool > &bools, const bool val=true)
Count number of 'true' entries.
Definition: BitOps.H:77
foamToVtk_WRITE_FIELD
#define foamToVtk_WRITE_FIELD(FieldType)
foamVtkInternalWriter.H
Foam::tmp::valid
bool valid() const noexcept
Identical to good(), or bool operator.
Definition: tmp.H:292
Foam::GeometricField< scalar, fvPatchField, volMesh >
patchInterps
PtrList< PrimitivePatchInterpolation< primitivePatch > > patchInterps
Definition: convertVolumeFields.H:130
pInterp
autoPtr< volPointInterpolation > pInterp
Definition: convertVolumeFields.H:81