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) 2016-2022 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 \*---------------------------------------------------------------------------*/
27 
28 #include "ensightOutputVolField.H"
29 #include "ensightMesh.H"
30 
31 #include "fvMesh.H"
32 #include "linear.H"
33 #include "volPointInterpolation.H"
34 #include "interpolation.H"
36 #include "DynamicField.H"
37 
38 // * * * * * * * * * * * * * * * * Detail * * * * * * * * * * * * * * * * * //
39 
40 template<class Type>
42 (
44  ensightFile& os,
46  const ensightMesh& ensMesh
47 )
48 {
49  bool parallel = Pstream::parRun();
50 
51  const fvMesh& mesh = vf.mesh();
52  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
53 
54  const Map<ensightCells>& cellZoneParts = ensMesh.cellZoneParts();
55  const Map<ensightFaces>& faceZoneParts = ensMesh.faceZoneParts();
56  const Map<ensightFaces>& boundaryParts = ensMesh.boundaryParts();
57 
58  // Write internalMesh and cellZones - sorted by index
59  for (const label zoneId : cellZoneParts.sortedToc())
60  {
61  const ensightCells& part = cellZoneParts[zoneId];
62 
63  ensightOutput::writeField(scratch, os, vf, part, parallel);
64  }
65 
66 
67  // Write patches - sorted by index
68  for (const label patchId : boundaryParts.sortedToc())
69  {
70  const ensightFaces& part = boundaryParts[patchId];
71 
72  if (patchId < 0 || patchId >= pbm.size())
73  {
74  // Future handling of combined patches?
75  continue;
76  }
77 
78  const label patchStart = pbm[patchId].start();
79 
80  // Either use a flat boundary field for all patches,
81  // or patch-local face ids
82 
83  // Operate on a copy
84  ensightFaces localPart(part);
85 
86  // Change from global faceIds to patch-local
87  localPart.decrFaceIds(patchStart);
88 
90  (
91  scratch,
92  os,
93  vf.boundaryField()[patchId],
94  localPart,
95  parallel
96  );
97  }
98 
99 
100  // No face zones data
101  if (faceZoneParts.empty())
102  {
103  return true;
104  }
105 
106 
107  // Flat boundary field
108  // similar to volPointInterpolation::flatBoundaryField()
109 
110  Field<Type> flat(pbm.nFaces(), Foam::zero{});
111 
112  forAll(vf.boundaryField(), patchi)
113  {
114  const polyPatch& pp = pbm[patchi];
115  const auto& pfld = vf.boundaryField()[patchi];
116 
117  // Note: restrict transcribing to actual size of the patch field
118  // - handles "empty" patch type etc.
119 
120  SubList<Type> slice(flat, pfld.size(), pp.offset());
121 
122  if (isA<processorPolyPatch>(pp))
123  {
124  // Use average value for processor faces
125  // own cell value = patchInternalField
126  // nei cell value = evaluated boundary values
127  slice = (0.5 * (pfld.patchInternalField() + pfld));
128  }
129  else if (!isA<emptyPolyPatch>(pp))
130  {
131  slice = pfld;
132  }
133  }
134 
135 
136  // Interpolate cell values to faces
137  // - only needed when exporting faceZones...
138 
141 
142  const auto& sfld = tsfld();
143 
144 
145  // Local output buffer
146  label maxLen = 0;
147 
148  forAllConstIters(faceZoneParts, iter)
149  {
150  maxLen = max(maxLen, iter.val().size());
151  }
152 
153  DynamicField<Type> values(maxLen);
154 
155 
156  //
157  // Write requested faceZones - sorted by index
158  //
159  for (const label zoneId : faceZoneParts.sortedToc())
160  {
161  const ensightFaces& part = faceZoneParts[zoneId];
162 
163  // Loop over face ids to store the needed field values
164  // - internal faces use linear interpolation
165  // - boundary faces use the corresponding patch value
166 
167  // Local copy of the field
168  values.resize_nocopy(part.size());
169  values = Zero;
170 
171  auto valIter = values.begin();
172 
173  for (const label faceId : part.faceIds())
174  {
175  *valIter =
176  (
177  mesh.isInternalFace(faceId)
178  ? sfld[faceId]
179  : flat[faceId - mesh.nInternalFaces()]
180  );
181 
182  ++valIter;
183  }
184 
185  // The field is already in the proper element order
186  // - just need its corresponding sub-fields
188  (
189  scratch, os, values, part, parallel
190  );
191  }
193  return true;
194 }
195 
196 
197 template<class Type>
199 (
201  ensightFile& os,
203  const ensightMesh& ensMesh
204 )
205 {
206  bool parallel = Pstream::parRun();
207 
208  const polyMesh& mesh = ensMesh.mesh();
209 
210  const Map<ensightCells>& cellZoneParts = ensMesh.cellZoneParts();
211  const Map<ensightFaces>& faceZoneParts = ensMesh.faceZoneParts();
212  const Map<ensightFaces>& boundaryParts = ensMesh.boundaryParts();
213 
214  //
215  // Write internalMesh and cellZones - sorted by index
216  //
217  for (const label zoneId : cellZoneParts.sortedToc())
218  {
219  const ensightCells& part = cellZoneParts[zoneId];
220 
221  labelList uniquePointLabels;
222  part.uniqueMeshPoints(mesh, uniquePointLabels, parallel);
223 
224  if (Pstream::master())
225  {
226  os.beginPart(part.index());
227  }
228 
229  // Skip if empty
230  if
231  (
232  parallel
233  ? returnReduceAnd(uniquePointLabels.empty())
234  : uniquePointLabels.empty()
235  )
236  {
237  continue;
238  }
239 
241  (
242  scratch,
243  os,
245  UIndirectList<Type>(pf.internalField(), uniquePointLabels),
246  parallel
247  );
248  }
249 
250 
251  //
252  // Write patches - sorted by index
253  //
254  for (const label patchId : boundaryParts.sortedToc())
255  {
256  const ensightFaces& part = boundaryParts[patchId];
257 
258  labelList uniquePointLabels;
259  part.uniqueMeshPoints(mesh, uniquePointLabels, parallel);
260 
261  if (Pstream::master())
262  {
263  os.beginPart(part.index());
264  }
265 
266  // Skip if empty
267  if
268  (
269  parallel
270  ? returnReduceAnd(uniquePointLabels.empty())
271  : uniquePointLabels.empty()
272  )
273  {
274  continue;
275  }
276 
277  const auto& bfld = pf.boundaryField()[patchId];
278 
279  // Only valuePointPatchField is actually derived from Field
280  const auto* vpp = isA<Field<Type>>(bfld);
281  if (vpp)
282  {
283  // Require patch local indices, not mesh point labels.
284  // But need to use polyPatch meshPointMap() to recover the
285  // local indices since the ensight output will have jumbled
286  // the face output order
287 
288  const polyPatch& pp = mesh.boundaryMesh()[patchId];
289 
290  for (label& pointi : uniquePointLabels)
291  {
292  pointi = pp.meshPointMap()[pointi];
293  }
294 
296  (
297  scratch,
298  os,
300  UIndirectList<Type>(*vpp, uniquePointLabels),
301  parallel
302  );
303  }
304  else
305  {
307  (
308  scratch,
309  os,
311  UIndirectList<Type>(pf.internalField(), uniquePointLabels),
312  parallel
313  );
314  }
315  }
316 
317  //
318  // Write requested faceZones - sorted by index
319  //
320  for (const label zoneId : faceZoneParts.sortedToc())
321  {
322  const ensightFaces& part = faceZoneParts[zoneId];
323 
324  labelList uniquePointLabels;
325  part.uniqueMeshPoints(mesh, uniquePointLabels, parallel);
326 
327  if (Pstream::master())
328  {
329  os.beginPart(part.index());
330  }
331 
332  // Skip if empty
333  if
334  (
335  parallel
336  ? returnReduceAnd(uniquePointLabels.empty())
337  : uniquePointLabels.empty()
338  )
339  {
340  continue;
341  }
342 
343  // CAVEAT - does not properly handle valuePointPatchField,
344  // uses internalField only
345 
347  (
348  scratch,
349  os,
351  UIndirectList<Type>(pf.internalField(), uniquePointLabels),
352  parallel
353  );
354  }
355 
356  return true;
357 }
358 
359 
360 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
361 
362 template<class Type>
364 (
366  ensightFile& os,
368  const ensightMesh& ensMesh,
369  const bool nodeValues
370 )
371 {
372  if (nodeValues)
373  {
375  (
377  );
378  pfld.ref().checkOut();
379  pfld.ref().rename(vf.name());
380 
381  return ensightOutput::writePointField<Type>(scratch, os, pfld, ensMesh);
382  }
383 
384  return ensightOutput::writeVolField<Type>(scratch, os, vf, ensMesh);
385 }
386 
387 
388 // ************************************************************************* //
A variant of OFstream with specialised handling for Ensight writing of strings, integers and floats (...
Definition: ensightFile.H:47
label patchId(-1)
const polyBoundaryMesh & pbm
label uniqueMeshPoints(const polyMesh &mesh, labelList &uniqueMeshPointLabels, bool parallel) const
Globally unique mesh points. Required when writing point fields.
label faceId(-1)
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > linearInterpolate(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: linear.H:120
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition: tmpI.H:235
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
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.
const Internal & internalField() const noexcept
Return a const-reference to the dimensioned internal field.
bool empty() const noexcept
True if List is empty (ie, size() is zero)
Definition: UList.H:675
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
Sorting/classification of faces (2D) into corresponding ensight types.
Definition: ensightFaces.H:66
bool writeFaceSubField(ensightOutput::floatBufferType &scratch, ensightFile &os, const Field< Type > &fld, const ensightFaces &part, bool parallel)
Write a sub-field of faces values as an indirect list, using the sub-list sizing information from ens...
Generic GeometricField class.
A collection of functions for writing volField content in ensight format.
bool writeField(ensightOutput::floatBufferType &scratch, ensightFile &os, const Field< Type > &fld, const ensightCells &part, bool parallel)
Write a field of cell values as an indirect list, using the cell ids from ensightCells.
const polyMesh & mesh() const noexcept
Reference to the underlying polyMesh.
Definition: ensightMesh.H:191
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
List< T > values(const HashTable< T, Key, Hash > &tbl, const bool doSort=false)
List of values from HashTable, optionally sorted.
Definition: HashOps.H:164
bool returnReduceAnd(const bool value, const label comm=UPstream::worldComm)
Perform logical (and) MPI Allreduce on a copy. Uses UPstream::reduceAnd.
dynamicFvMesh & mesh
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:51
void writeFieldComponents(ensightOutput::floatBufferType &scratch, ensightFile &os, const char *key, const FieldContainer< Type > &fld, bool parallel)
Write field content (component-wise) for the given ensight element type.
Generic templated field type.
Definition: Field.H:62
label size(const elemType etype) const
Processor-local size of the specified element type.
Definition: ensightFacesI.H:60
Sorting/classification of cells (3D) into corresponding ensight element types.
Definition: ensightCells.H:50
const Map< ensightFaces > & boundaryParts() const noexcept
Face elements per selected patch, lookup by patch index.
Definition: ensightMesh.H:227
Dynamically sized Field.
Definition: DynamicField.H:45
const Map< ensightFaces > & faceZoneParts() const noexcept
Face elements per faceZone, lookup by zone index.
Definition: ensightMesh.H:217
const Map< ensightCells > & cellZoneParts() const noexcept
Face elements per selected patch, lookup by patch index.
Definition: ensightMesh.H:207
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO...
const labelList & faceIds() const noexcept
Processor-local face ids of all elements.
Definition: ensightFacesI.H:72
const Mesh & mesh() const noexcept
Return mesh.
bool empty() const noexcept
True if the hash table is empty.
Definition: HashTable.H:353
Encapsulation of volume meshes for writing in ensight format. It manages cellZones, facesZone, patches.
Definition: ensightMesh.H:78
OBJstream os(runTime.globalPath()/outputName)
void decrFaceIds(const label off)
Decrease face ids by specified offset value.
label index() const noexcept
The index in a list (0-based)
Definition: ensightPart.H:153
PtrList< coordinateSystem > coordinates(solidRegions.size())
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
List< label > sortedToc() const
The table of contents (the keys) in sorted order.
Definition: HashTable.C:163
A List with indirect addressing. Like IndirectList but does not store addressing. ...
Definition: faMatrix.H:52
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition: zero.H:57
bool writePointField(ensightOutput::floatBufferType &scratch, ensightFile &os, const GeometricField< Type, pointPatchField, pointMesh > &pf, const ensightMesh &ensMesh)
Write point field component-wise.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:75
List< label > labelList
A List of labels.
Definition: List.H:62
A class for managing temporary objects.
Definition: HashPtrTable.H:50
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:69
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
A HashTable to objects of type <T> with a label key.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
bool writeVolField(ensightOutput::floatBufferType &scratch, ensightFile &os, const GeometricField< Type, fvPatchField, volMesh > &vf, const ensightMesh &ensMesh)
Write volume field component-wise.