ensightOutputVolFieldTemplates.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) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2016-2019 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "ensightOutputVolField.H"
30 #include "ensightMesh.H"
31 
32 #include "fvMesh.H"
33 #include "globalIndex.H"
34 #include "volPointInterpolation.H"
35 #include "interpolation.H"
36 #include "linear.H"
37 #include "processorFvPatch.H"
39 
40 // * * * * * * * * * * * * * * * * Detail * * * * * * * * * * * * * * * * * //
41 
42 template<class Type>
44 (
46  const ensightMesh& ensMesh,
47  ensightFile& os
48 )
49 {
50  constexpr bool parallel = true;
51 
52  const fvMesh& mesh = ensMesh.mesh();
53  const ensightCells& meshCells = ensMesh.meshCells();
54  const Map<word>& patchLookup = ensMesh.patches();
55  const HashTable<ensightFaces>& patchFaces = ensMesh.boundaryPatchFaces();
56  const HashTable<ensightFaces>& zoneFaces = ensMesh.faceZoneFaces();
57 
58  //
59  // write internalMesh, unless patch-selection was requested
60  //
61  if (ensMesh.useInternalMesh())
62  {
63  Detail::writeCellField(vf, meshCells, os, parallel);
64  }
65 
66  //
67  // write patches
68  // use sortedToc for extra safety
69  //
70  const labelList patchIds = patchLookup.sortedToc();
71  for (const label patchId : patchIds)
72  {
73  const word& patchName = patchLookup[patchId];
74  const ensightFaces& part = patchFaces[patchName];
75 
77  (
78  vf.boundaryField()[patchId],
79  part,
80  os,
81  parallel
82  );
83  }
84 
85 
86  //
87  // write faceZones, if requested
88  // use sortedToc for extra safety
89  //
90  const wordList zoneNames = zoneFaces.sortedToc();
91  if (!zoneNames.empty())
92  {
93  // Interpolates cell values to faces - needed only when exporting
94  // faceZones...
96  (
98  );
99 
100  // flat boundary field
101  // similar to volPointInterpolation::flatBoundaryField()
102 
103  Field<Type> flat(mesh.nBoundaryFaces(), Zero);
104 
105  const fvBoundaryMesh& bm = mesh.boundary();
106  forAll(vf.boundaryField(), patchi)
107  {
108  const polyPatch& pp = bm[patchi].patch();
109  const auto& bf = vf.boundaryField()[patchi];
110 
111  if (isA<processorFvPatch>(bm[patchi]))
112  {
113  // Use average value for processor faces
114  // own cell value = patchInternalField
115  // nei cell value = evaluated boundary values
117  (
118  flat,
119  bf.size(),
120  pp.offset()
121  ) = (0.5 * (bf.patchInternalField() + bf));
122  }
123  else if (!isA<emptyFvPatch>(bm[patchi]))
124  {
126  (
127  flat,
128  bf.size(),
129  pp.offset()
130  ) = bf;
131  }
132  }
133 
134  for (const word& zoneName : zoneNames)
135  {
136  const ensightFaces& part = zoneFaces[zoneName];
137 
138  // Field (local size)
139  Field<Type> values(part.size());
140 
141  // Loop over face ids to store the needed field values
142  // - internal faces use linear interpolation
143  // - boundary faces use the corresponding patch value
144  forAll(part, i)
145  {
146  const label faceId = part[i];
147  values[i] =
148  (
149  mesh.isInternalFace(faceId)
150  ? sf[faceId]
151  : flat[faceId - mesh.nInternalFaces()]
152  );
153  }
154 
155  // The field is already copied in the proper order
156  // - just need its corresponding sub-fields
157  Detail::writeFaceSubField(values, part, os, parallel);
158  }
159  }
160 
161  return true;
162 }
163 
164 
165 template<class Type>
167 (
169  const ensightMesh& ensMesh,
170  ensightFile& os
171 )
172 {
173  constexpr bool parallel = true;
174 
175  const fvMesh& mesh = ensMesh.mesh();
176  const Map<word>& patchLookup = ensMesh.patches();
177 
178  const HashTable<ensightFaces>& patchFaces = ensMesh.boundaryPatchFaces();
179  const HashTable<ensightFaces>& zoneFaces = ensMesh.faceZoneFaces();
180 
181  //
182  // write internalMesh, unless patch-selection was requested
183  //
184  if (ensMesh.useInternalMesh())
185  {
186  if (Pstream::master())
187  {
188  os.beginPart(0); // 0 = internalMesh
189  }
190 
192  (
193  "coordinates",
194  Field<Type>(pf.internalField(), ensMesh.uniquePointMap()),
195  os,
196  parallel
197  );
198  }
199 
200  //
201  // write patches
202  // use sortedToc for extra safety
203  //
204  const labelList patchIds = patchLookup.sortedToc();
205  for (const label patchId : patchIds)
206  {
207  const word& patchName = patchLookup[patchId];
208  const ensightFaces& part = patchFaces[patchName];
209 
210  const fvPatch& p = mesh.boundary()[patchId];
211 
212  // Renumber the patch points/faces into unique points
213  labelList pointToGlobal;
214  labelList uniqueMeshPointLabels;
215  autoPtr<globalIndex> globalPointsPtr =
216  mesh.globalData().mergePoints
217  (
218  p.patch().meshPoints(),
219  p.patch().meshPointMap(),
220  pointToGlobal,
221  uniqueMeshPointLabels
222  );
223 
224  if (Pstream::master())
225  {
226  os.beginPart(part.index());
227  }
228 
230  (
231  "coordinates",
232  Field<Type>(pf.internalField(), uniqueMeshPointLabels),
233  os,
234  parallel
235  );
236  }
237 
238  //
239  // write faceZones, if requested
240  //
241  const wordList zoneNames = zoneFaces.sortedToc();
242  for (const word& zoneName : zoneNames)
243  {
244  const ensightFaces& part = zoneFaces[zoneName];
245 
247  (
249  (
250  mesh.faces(),
251  part.faceIds()
252  ),
253  mesh.points()
254  );
255 
256  // Renumber the patch points/faces into unique points
257  labelList pointToGlobal;
258  labelList uniqueMeshPointLabels;
259  autoPtr<globalIndex> globalPointsPtr =
260  mesh.globalData().mergePoints
261  (
262  p.meshPoints(),
263  p.meshPointMap(),
264  pointToGlobal,
265  uniqueMeshPointLabels
266  );
267 
268  if (Pstream::master())
269  {
270  os.beginPart(part.index());
271  }
272 
274  (
275  "coordinates",
276  Field<Type>(pf.internalField(), uniqueMeshPointLabels),
277  os,
278  parallel
279  );
280  }
281 
282  return true;
283 }
284 
285 
286 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
287 
288 template<class Type>
290 (
292  const ensightMesh& ensMesh,
293  ensightFile& os,
294  const bool nodeValues
295 )
296 {
297  if (nodeValues)
298  {
300  (
302  );
303  pfld.ref().checkOut();
304  pfld.ref().rename(vf.name());
305 
306  return Detail::writePointField<Type>(pfld, ensMesh, os);
307  }
308 
309  return Detail::writeVolField<Type>(vf, ensMesh, os);
310 }
311 
312 
313 // * * * * * * * * * * * * * * * * Serial * * * * * * * * * * * * * * * * * //
314 
315 template<class Type>
317 (
319  const ensightPartFaces& part,
320  ensightFile& os
321 )
322 {
323  constexpr bool parallel = false; // serial
324 
325  const label patchi = part.patchIndex();
326 
327  if (patchi >= 0 && patchi < vf.boundaryField().size())
328  {
330  (
331  vf.boundaryField()[patchi],
332  part,
333  os,
334  parallel
335  );
336  }
337 
338  return false;
339 }
340 
341 
342 template<class Type>
344 (
346  const ensightPartCells& part,
347  ensightFile& os
348 )
349 {
350  constexpr bool parallel = false; // serial
351 
353  (
354  vf.internalField(),
355  part,
356  os,
357  parallel
358  );
359 }
360 
361 
362 template<class Type>
364 (
366  const ensightParts& list,
367  ensightFile& os
368 )
369 {
370  for (const ensightPart& part : list)
371  {
372  const ensightPartFaces* fptr = isA<ensightPartFaces>(part);
373 
374  if (fptr)
375  {
376  Serial::writeVolField(vf, *fptr, os);
377  continue;
378  }
379 
380  const ensightPartCells* cptr = isA<ensightPartCells>(part);
381 
382  if (cptr)
383  {
384  Serial::writeVolField(vf, *cptr, os);
385  continue;
386  }
387  }
388 
389  return true;
390 }
391 
392 
393 // ************************************************************************* //
Foam::ensightOutput::Serial::writeVolField
bool writeVolField(const GeometricField< Type, fvPatchField, volMesh > &vf, const ensightPartFaces &part, ensightFile &os)
Write volume field component-wise for specified faces.
Definition: ensightOutputVolFieldTemplates.C:317
Foam::ensightParts
A collection of several ensightPart elements.
Definition: ensightParts.H:54
Foam::ensightMesh
Encapsulation of volume meshes for writing in ensight format.
Definition: ensightMesh.H:68
p
volScalarField & p
Definition: createFieldRefs.H:8
processorFvPatch.H
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
uindirectPrimitivePatch.H
Foam::ensightOutput::Detail::writeVolField
bool writeVolField(const GeometricField< Type, fvPatchField, volMesh > &vf, const ensightMesh &ensMesh, ensightFile &os)
Write volume field component-wise.
Definition: ensightOutputVolFieldTemplates.C:44
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::HashTableOps::values
List< T > values(const HashTable< T, Key, Hash > &tbl, const bool doSort=false)
List of values from HashTable, optionally sorted.
Definition: HashOps.H:149
Foam::ensightMesh::useInternalMesh
bool useInternalMesh() const
Using internal?
Definition: ensightMeshI.H:48
Foam::ensightMesh::faceZoneFaces
const HashTable< ensightFaces > & faceZoneFaces() const
Face elements per selected faceZone.
Definition: ensightMeshI.H:80
Foam::SubList
A List obtained as a section of another List.
Definition: SubList.H:53
Foam::ensightMesh::mesh
const fvMesh & mesh() const
Reference to the underlying fvMesh.
Definition: ensightMeshI.H:30
globalIndex.H
Foam::Map
A HashTable to objects of type <T> with a label key.
Definition: HashTableFwd.H:46
Foam::polyPatch::offset
label offset() const
The offset where this patch starts in the boundary face list.
Definition: polyPatch.C:299
interpolation.H
Foam::ensightPartFaces::patchIndex
label patchIndex() const
Return the patch index, -1 when not in use.
Definition: ensightPartFaces.H:171
Foam::GeometricField::internalField
const Internal & internalField() const
Return a const-reference to the dimensioned internal field.
Definition: GeometricFieldI.H:43
Foam::ensightPartCells
An implementation of ensightPart to hold volume mesh cells.
Definition: ensightPartCells.H:53
Foam::ensightOutput::Detail::writeFieldComponents
bool writeFieldComponents(const char *key, const FieldContainer< Type > &fld, ensightFile &os, bool parallel)
Write field content (component-wise) for the given ensight element type.
Definition: ensightOutputTemplates.C:35
Foam::ensightFaces
Sorting/classification of faces (2D) into corresponding ensight types.
Definition: ensightFaces.H:51
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::ensightOutput::Detail::writeCellField
bool writeCellField(const Field< Type > &fld, const ensightCells &part, ensightFile &os, bool parallel)
Definition: ensightOutputTemplates.C:275
Foam::ensightOutput::Detail::writePointField
bool writePointField(const GeometricField< Type, pointPatchField, pointMesh > &pf, const ensightMesh &ensMesh, ensightFile &os)
Write point field component-wise.
Definition: ensightOutputVolFieldTemplates.C:167
Foam::fvBoundaryMesh
Foam::fvBoundaryMesh.
Definition: fvBoundaryMesh.H:57
ensightMesh.H
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:258
patchIds
labelList patchIds
Definition: convertProcessorPatches.H:67
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:62
Foam::Field
Generic templated field type.
Definition: Field.H:63
Foam::ensightOutput::writeVolField
bool writeVolField(const GeometricField< Type, fvPatchField, volMesh > &, const ensightMesh &ensMesh, ensightFile &os, const bool nodeValues=false)
Write volume field component-wise.
Definition: ensightOutputVolFieldTemplates.C:290
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
faceId
label faceId(-1)
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
Foam::ensightFaces::faceIds
const labelUList faceIds(const enum elemType) const
Return the (local) face ids of the specified element type.
Definition: ensightFacesI.H:81
Foam::ensightCells
Sorting/classification of cells (3D) into corresponding ensight element types.
Definition: ensightCells.H:53
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
ensightOutputVolField.H
A collection of functions for writing volField content in ensight format.
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:84
fvMesh.H
Foam::ensightFile
Ensight output with specialized write() for strings, integers and floats. Correctly handles binary wr...
Definition: ensightFile.H:54
Foam::ensightMesh::boundaryPatchFaces
const HashTable< ensightFaces > & boundaryPatchFaces() const
Face elements per selected patch.
Definition: ensightMeshI.H:73
Foam::ensightFaces::size
label size() const
The processor local size of all elements.
Definition: ensightFacesI.H:50
Foam::HashTable::sortedToc
List< Key > sortedToc() const
The table of contents (the keys) in sorted order.
Definition: HashTable.C:136
Foam::HashTable
A HashTable similar to std::unordered_map.
Definition: HashTable.H:105
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
volPointInterpolation.H
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::ensightMesh::patches
const Map< word > & patches() const
The list of patches to be output.
Definition: ensightMeshI.H:66
Foam::ensightFaces::index
label index() const
The index in a list.
Definition: ensightFacesI.H:38
Foam::linearInterpolate
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > linearInterpolate(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: linear.H:112
Foam::List< label >
Foam::ensightPart
Base class for ensightPartCells and ensightPartFaces.
Definition: ensightPart.H:56
Foam::ensightFile::beginPart
void beginPart(const label index)
Begin a part (0-based index internally).
Definition: ensightFile.C:319
Foam::ensightOutput::Detail::writeFaceField
bool writeFaceField(const Field< Type > &fld, const ensightFaces &part, ensightFile &os, bool parallel)
Definition: ensightOutputTemplates.C:119
Foam::ensightPartFaces
An implementation of ensightPart to hold mesh faces.
Definition: ensightPartFaces.H:53
Foam::UIndirectList
A List with indirect addressing.
Definition: fvMatrix.H:109
patchId
label patchId(-1)
Foam::GeometricField< Type, fvPatchField, volMesh >
Foam::fac::interpolate
static tmp< GeometricField< Type, faePatchField, edgeMesh > > interpolate(const GeometricField< Type, faPatchField, areaMesh > &tvf, const edgeScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
Foam::ensightMesh::meshCells
const ensightCells & meshCells() const
The volume cells (internalMesh)
Definition: ensightMeshI.H:60
Foam::ensightMesh::uniquePointMap
const labelList & uniquePointMap() const
Local points that are unique.
Definition: ensightMesh.H:362
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatch.H:90
Foam::ensightOutput::Detail::writeFaceSubField
bool writeFaceSubField(const Field< Type > &fld, const ensightFaces &part, ensightFile &os, bool parallel)
Definition: ensightOutputTemplates.C:171