convertVolumeFields.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-2023 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM, distributed under GPL-3.0-or-later.
12 
13 Description
14  Code chunk for converting volume and dimensioned fields
15  included by foamToVTK.
16 
17 \*---------------------------------------------------------------------------*/
18 
19 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
20 
21 {
22  using reportFields = foamToVtkReportFields;
23 
24  const label nVolFields =
25  (
26  (doInternal || doBoundary)
27  ? objects.count(stringListOps::foundOp<word>(fieldTypes::volume))
28  : 0
29  );
30 
31  const label nDimFields =
32  (
33  (doInternal || doBoundary)
34  ? objects.count(stringListOps::foundOp<word>(fieldTypes::internal))
35  : 0
36  );
37 
38  label nPointFields =
39  (
40  doPointValues
41  ? objects.count(stringListOps::foundOp<word>(fieldTypes::point))
42  : 0
43  );
44 
45 
46  if (doInternal || doBoundary)
47  {
49  }
50  if (doInternal)
51  {
53  }
54 
55 
56  // Setup for the vtm writer.
57  // For legacy format, the information added is simply ignored.
58 
59  fileName vtmOutputBase
60  (
61  outputDir/regionDir/vtkName + timeDesc
62  );
63 
64  // Combined internal + boundary in a vtm file
66 
67  // Collect individual boundaries into a vtm file
69 
70  // Setup the internal writer
71  autoPtr<vtk::internalWriter> internalWriter;
72 
73  // Interpolator for volume and dimensioned fields
74  autoPtr<volPointInterpolation> pInterp;
75 
76  if (doInternal)
77  {
78  if (vtuMeshCells.empty())
79  {
80  // Use the appropriate mesh (baseMesh or subMesh)
81  vtuMeshCells.reset(meshProxy.mesh());
82 
83  if (doPointValues && vtuMeshCells.manifold())
84  {
85  doPointValues = false;
86  nPointFields = 0;
87  Warning
88  << nl
89  << "Manifold cells detected (eg, AMI) - disabling PointData"
90  << nl
91  << endl;
92  }
93  }
94 
96  (
97  meshProxy.mesh(),
98  vtuMeshCells,
99  writeOpts,
100  // The output base name for internal
101  (
102  writeOpts.legacy()
103  ? vtmOutputBase
104  : (vtmOutputBase / "internal")
105  ),
106  UPstream::parRun()
107  );
108 
109  // No sub-block for internal
110  vtmWriter.append_vtu
111  (
112  "internal",
113  vtmOutputBase.name()/"internal"
114  );
115 
116  Info<< " Internal : "
117  << args.relativePath(internalWriter->output()) << nl;
118 
119  internalWriter->writeTimeValue(mesh.time().value());
120  internalWriter->writeGeometry();
121 
122  if (doPointValues)
123  {
124  pInterp.reset(new volPointInterpolation(mesh));
125  }
126  }
127 
128 
129  // Setup the patch writers
130 
131  const polyBoundaryMesh& patches = mesh.boundaryMesh();
132 
133  PtrList<vtk::patchWriter> patchWriters;
134  PtrList<PrimitivePatchInterpolation<primitivePatch>> patchInterps;
137  if (doBoundary)
138  {
139  patchIds = getSelectedPatches(patches, patchSelector);
140  }
141 
142  if (oneBoundary && patchIds.size())
143  {
145  (
146  meshProxy.mesh(),
147  patchIds,
148  writeOpts,
149  nearCellValue,
150  // Output one patch: "boundary"
151  (
152  writeOpts.legacy()
153  ?
154  (
155  outputDir/regionDir/"boundary"
156  / (meshProxy.useSubMesh() ? meshProxy.name() : "boundary")
157  + timeDesc
158  )
159  : (vtmOutputBase / "boundary")
160  ),
161  UPstream::parRun()
162  );
163 
164  // No sub-block for one-patch
165  vtmWriter.append_vtp
166  (
167  "boundary",
168  vtmOutputBase.name()/"boundary"
169  );
170 
171  Info<< " Boundaries: "
172  << args.relativePath(writer->output()) << nl;
173 
174  writer->writeTimeValue(timeValue);
175  writer->writeGeometry();
176 
177  // Transfer writer to list for later use
178  patchWriters.resize(1);
179  patchWriters.set(0, writer);
180 
181  // Avoid patchInterpolation for each sub-patch
182  patchInterps.resize(1); // == nullptr
183  }
184  else if (patchIds.size())
185  {
186  patchWriters.resize(patchIds.size());
187  if (doPointValues)
188  {
189  patchInterps.resize(patchIds.size());
190  }
191 
192  label nPatchWriters = 0;
193  label nPatchInterps = 0;
194 
195  for (const label patchId : patchIds)
196  {
197  const polyPatch& pp = patches[patchId];
198 
200  (
201  meshProxy.mesh(),
202  labelList(one{}, pp.index()),
203  writeOpts,
204  nearCellValue,
205  // Output patch: "boundary"/name
206  (
207  writeOpts.legacy()
208  ?
209  (
210  outputDir/regionDir/pp.name()
211  / (meshProxy.useSubMesh() ? meshProxy.name() : pp.name())
212  + timeDesc
213  )
214  : (vtmOutputBase / "boundary" / pp.name())
215  ),
216  UPstream::parRun()
217  );
218 
219  if (!nPatchWriters)
220  {
221  vtmWriter.beginBlock("boundary");
222  vtmBoundaries.beginBlock("boundary");
223  }
224 
225  vtmWriter.append_vtp
226  (
227  pp.name(),
228  vtmOutputBase.name()/"boundary"/pp.name()
229  );
230 
231  vtmBoundaries.append_vtp
232  (
233  pp.name(),
234  "boundary"/pp.name()
235  );
236 
237  Info<< " Boundary : "
238  << args.relativePath(writer->output()) << nl;
239 
240  writer->writeTimeValue(timeValue);
241  writer->writeGeometry();
242 
243  // Transfer writer to list for later use
245 
246  if (patchInterps.size())
247  {
248  patchInterps.set
249  (
250  nPatchInterps++,
251  new PrimitivePatchInterpolation<primitivePatch>(pp)
252  );
253  }
254  }
255 
256  if (nPatchWriters)
257  {
258  vtmWriter.endBlock("boundary");
259  }
260 
261  patchWriters.resize(nPatchWriters);
262  patchInterps.resize(nPatchInterps);
263  }
264 
265  // With point-interpolation, cache fields to avoid multiple re-reading
266  std::unique_ptr<objectRegistry> cacheFieldsPtr;
267 
268  if (doPointValues && (nVolFields || nDimFields))
269  {
271  (
272  new objectRegistry
273  (
274  IOobject
275  (
276  "foamToVTK::volume",
277  runTime.timeName(),
278  runTime,
279  IOobject::NO_REGISTER
280  )
281  )
282  );
283  }
284 
285 
286  // CellData
287  {
288  // Begin CellData
289  if (internalWriter)
290  {
291  // Optionally with (cellID, procID) fields
292  internalWriter->beginCellData
293  (
294  (withMeshIds ? 1 + (internalWriter->parallel() ? 1 : 0) : 0)
296  );
297 
298  if (withMeshIds)
299  {
300  internalWriter->writeCellIDs();
301  internalWriter->writeProcIDs(); // parallel only
302  }
303  }
304 
305  if (nVolFields || withMeshIds)
306  {
307  for (vtk::patchWriter& writer : patchWriters)
308  {
309  // Optionally with (patchID, procID) fields
310  writer.beginCellData
311  (
312  (withMeshIds ? 2 : 0)
313  + nVolFields
314  );
315 
316  if (withMeshIds)
317  {
318  writer.writePatchIDs();
319  writer.writeProcIDs();
320  }
321  }
322  }
323 
325  (
327  patchWriters,
328  meshProxy,
329  objects,
330  true, // syncPar
331  cacheFieldsPtr.get()
332  );
333 
335  (
337  meshProxy,
338  objects,
339  true, // syncPar
340  cacheFieldsPtr.get()
341  );
342 
343  // End CellData is implicit
344  }
345 
346 
347  // PointData
348  // - only construct pointMesh on request since it constructs
349  // edge addressing
350  if (doPointValues)
351  {
352  // Begin PointData
353  if (internalWriter)
354  {
355  internalWriter->beginPointData
356  (
358  + (withPointIds ? 1 : 0)
359  );
360 
361  if (withPointIds)
362  {
363  internalWriter->writePointIDs();
364  }
365  }
366 
367  forAll(patchWriters, writeri)
368  {
369  const label nPatchFields =
370  (
371  (patchInterps.test(writeri) ? nVolFields : 0)
372  + nPointFields
373  );
374 
375  if (nPatchFields)
376  {
377  patchWriters[writeri].beginPointData(nPatchFields);
378  }
379  }
380 
382  (
385  meshProxy,
386  objects,
387  true, // syncPar
388  cacheFieldsPtr.get()
389  );
390 
392  (
394  meshProxy,
395  objects,
396  true, // syncPar
397  cacheFieldsPtr.get()
398  );
399 
401  (
403  patchWriters,
404  meshProxy,
405  objects,
406  true // syncPar
407  );
408 
409  // End PointData is implicit
410  }
411 
412 
413  // Finish writers
414  if (internalWriter)
415  {
416  internalWriter->close();
417  }
418 
419  for (vtk::patchWriter& writer : patchWriters)
420  {
421  writer.close();
422  }
423 
424  pInterp.clear();
425  patchWriters.clear();
426  patchInterps.clear();
427 
428 
429  // Collective output
430 
431  if (UPstream::master())
432  {
433  // Naming for vtm, file series etc.
434  fileName outputName(vtmOutputBase);
435 
436  if (writeOpts.legacy())
437  {
438  if (doInternal)
439  {
440  // Add to file-series and emit as JSON
441 
443 
444  fileName seriesName(vtk::seriesWriter::base(outputName));
445 
446  vtk::seriesWriter& series = vtkSeries(seriesName);
447 
448  // First time?
449  // Load from file, verify against filesystem,
450  // prune time >= currentTime
451  if (series.empty())
452  {
453  series.load(seriesName, true, timeValue);
454  }
455 
456  series.append(timeValue, timeDesc);
457  series.write(seriesName);
458  }
459  }
460  else
461  {
462  if (vtmWriter.size())
463  {
464  // Emit ".vtm"
465 
466  outputName.ext(vtmWriter.ext());
467 
468  fileName seriesName(vtk::seriesWriter::base(outputName));
469 
470  vtmWriter.setTime(timeValue);
471  vtmWriter.write(outputName);
472 
473  // Add to file-series and emit as JSON
474 
475  vtk::seriesWriter& series = vtkSeries(seriesName);
476 
477  // First time?
478  // Load from file, verify against filesystem,
479  // prune time >= currentTime
480  if (series.empty())
481  {
482  series.load(seriesName, true, timeValue);
483  }
484 
485  series.append(timeValue, outputName);
486  series.write(seriesName);
487 
488  // Add to multi-region vtm
489  vtmMultiRegion.add(regionName, regionDir, vtmWriter);
490  }
491 
492  if (vtmBoundaries.size())
493  {
494  // Emit "boundary.vtm" with collection of boundaries
495 
496  // Naming for vtm
497  fileName outputName(vtmOutputBase / "boundary");
498  outputName.ext(vtmBoundaries.ext());
499 
500  vtmBoundaries.setTime(timeValue);
501  vtmBoundaries.write(outputName);
502  }
503  }
504  }
505 }
506 
507 
508 // ************************************************************************* //
label patchId(-1)
labelList patchIds
autoPtr< vtk::internalWriter > internalWriter
vtk::lineWriter writer(edgeCentres, edgeList::null(), fileName(aMesh.time().globalPath()/"finiteArea-edgesCentres"))
label writeAllVolFields(ensightCase &ensCase, const ensightMesh &ensMesh, const IOobjectList &objects, const bool nearCellValue=false)
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
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.
PtrList< vtk::patchWriter > patchWriters
const label nDimFields
label writeAllPointFields(ensightCase &ensCase, const ensightMesh &ensMesh, const IOobjectList &objects)
std::unique_ptr< objectRegistry > cacheFieldsPtr
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
word outputName("finiteArea-edges.obj")
dynamicFvMesh & mesh
label nPatchWriters
fileName vtmOutputBase(outputDir/regionDir/vtkName+timeDesc)
const polyBoundaryMesh & patches
fileName relativePath(const fileName &input, const bool caseTag=false) const
Return the input relative to the globalPath by stripping off a leading value of the globalPath...
Definition: argListI.H:87
const label nVolFields
const wordList volume
Standard volume field types (scalar, vector, tensor, etc)
PtrList< PrimitivePatchInterpolation< primitivePatch > > patchInterps
messageStream Warning
Warning stream (stdout output on master, null elsewhere), with additional &#39;FOAM Warning&#39; header text...
vector point
Point is a vector.
Definition: point.H:37
const wordList internal
Standard dimensioned field types (scalar, vector, tensor, etc)
const Foam::Enum< fileTag > fileExtension
File extension (without ".") for some vtk XML file content types.
Foam::word regionName(args.getOrDefault< word >("region", Foam::polyMesh::defaultRegion))
vtk::vtmWriter vtmBoundaries
const word & regionDir
vtk::vtmWriter vtmWriter
label nPointFields
messageStream Info
Information stream (stdout output on master, null elsewhere)
List< label > labelList
A List of labels.
Definition: List.H:62
writeAllDimFields(internalWriter, meshProxy, objects, true, cacheFieldsPtr.get())
Foam::argList args(argc, argv)
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
autoPtr< volPointInterpolation > pInterp