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