foamToVTK.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2024 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Application
25  foamToVTK
26 
27 Description
28  Legacy VTK file format writer.
29 
30  - Handles volFields, volFields::Internal, pointFields,
31  surfaceScalarField and surfaceVectorField.
32  - Mesh topo changes.
33  - Both ascii and binary.
34  - Single time step writing.
35  - Write subset only.
36  - Automatic decomposition of cells; polygons on boundary undecomposed since
37  handled by vtk.
38 
39 Usage
40  \b foamToVTK [OPTION]
41 
42  Options:
43  - \par -ascii
44  Write VTK data in ASCII format instead of binary.
45 
46  - \par -mesh <name>
47  Use a different mesh name (instead of -region)
48 
49  - \par -fields <fields>
50  Convert selected fields only. For example,
51  \verbatim
52  -fields "( p T U )"
53  \endverbatim
54  The quoting is required to avoid shell expansions and to pass the
55  information as a single argument.
56 
57  - \par -surfaceFields
58  Write surfaceScalarFields (e.g., phi)
59 
60  - \par -cellSet <name>
61  - \par -faceSet <name>
62 
63  - \par -pointSet <name>
64  Restrict conversion to the cellSet, faceSet or pointSet.
65 
66  - \par -nearCellValue
67  Output cell value on patches instead of patch value itself
68 
69  - \par -noInternal
70  Do not generate file for mesh, only for patches
71 
72  - \par -noPointValues
73  No pointFields
74 
75  - \par -noFaceZones
76  No faceZones
77 
78  - \par -noLinks
79  (in parallel) do not link processor files to master
80 
81  - \par -polyhedra
82  Which cell types to write as polyhedra - 'none', 'polyhedra', or 'all'
83 
84  - \par -allPatches
85  Combine all patches into a single file
86 
87  - \par -excludePatches <patchNames>
88  Specify patches (wildcards) to exclude. For example,
89  \verbatim
90  -excludePatches '( inlet_1 inlet_2 "proc.*")'
91  \endverbatim
92  The quoting is required to avoid shell expansions and to pass the
93  information as a single argument. The double quotes denote a regular
94  expression.
95 
96  - \par -useTimeName
97  use the time index in the VTK file name instead of the time index
98 
99  Note:
100  mesh subset is handled by vtkMesh. Slight inconsistency in
101  interpolation: on the internal field it interpolates the whole volField
102  to the whole-mesh pointField and then selects only those values it
103  needs for the subMesh (using the fvMeshSubset cellMap(), pointMap()
104  functions). For the patches however it uses the
105  fvMeshSubset.interpolate function to directly interpolate the
106  whole-mesh values onto the subset patch.
107 
108  \par new file format:
109  no automatic timestep recognition.
110  However can have .pvd file format which refers to time simulation
111  if XML *.vtu files are available:
112 
113  \verbatim
114  <?xml version="1.0"?>
115  <VTKFile type="Collection"
116  version="0.1"
117  byte_order="LittleEndian"
118  compressor="vtkZLibDataCompressor">
119  <Collection>
120  <DataSet timestep="50" file="pitzDaily_2.vtu"/>
121  <DataSet timestep="100" file="pitzDaily_3.vtu"/>
122  <DataSet timestep="150" file="pitzDaily_4.vtu"/>
123  <DataSet timestep="200" file="pitzDaily_5.vtu"/>
124  <DataSet timestep="250" file="pitzDaily_6.vtu"/>
125  <DataSet timestep="300" file="pitzDaily_7.vtu"/>
126  <DataSet timestep="350" file="pitzDaily_8.vtu"/>
127  <DataSet timestep="400" file="pitzDaily_9.vtu"/>
128  <DataSet timestep="450" file="pitzDaily_10.vtu"/>
129  <DataSet timestep="500" file="pitzDaily_11.vtu"/>
130  </Collection>
131  </VTKFile>
132  \endverbatim
133 
134 \*---------------------------------------------------------------------------*/
135 
136 #include "argList.H"
137 #include "timeSelector.H"
138 #include "pointMesh.H"
139 #include "volPointInterpolation.H"
140 #include "emptyPolyPatch.H"
141 #include "nonConformalPolyPatch.H"
142 #include "labelIOField.H"
143 #include "scalarIOField.H"
144 #include "sphericalTensorIOField.H"
145 #include "symmTensorIOField.H"
146 #include "tensorIOField.H"
147 #include "faceZoneList.H"
148 #include "Cloud.H"
149 #include "passiveParticle.H"
150 #include "stringListOps.H"
151 
152 #include "vtkMesh.H"
153 #include "readFields.H"
154 #include "vtkWriteOps.H"
155 
156 #include "internalWriter.H"
157 #include "patchWriter.H"
158 #include "lagrangianWriter.H"
159 
160 #include "writeFaceSet.H"
161 #include "writePointSet.H"
162 #include "surfaceMeshWriter.H"
163 #include "writeSurfFields.H"
164 
165 using namespace Foam;
166 
167 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
168 
169 template<class GeoField>
170 void print(const char* msg, Ostream& os, const PtrList<const GeoField>& flds)
171 {
172  if (flds.size())
173  {
174  os << msg;
175  forAll(flds, i)
176  {
177  os << ' ' << flds[i].name();
178  }
179  os << endl;
180  }
181 }
182 
183 
184 void print(Ostream& os, const wordList& flds)
185 {
186  forAll(flds, i)
187  {
188  os << ' ' << flds[i];
189  }
190  os << endl;
191 }
192 
193 
194 labelList getSelectedPatches
195 (
196  const polyBoundaryMesh& patches,
197  const List<wordRe>& excludePatches
198 )
199 {
200  DynamicList<label> patchIDs(patches.size());
201 
202  Info<< "Combining patches:" << endl;
203 
205  {
206  const polyPatch& pp = patches[patchi];
207 
208  if
209  (
210  isA<emptyPolyPatch>(pp)
211  || isA<nonConformalPolyPatch>(pp)
212  || (Pstream::parRun() && isType<processorPolyPatch>(pp))
213  )
214  {
215  Info<< " discarding empty/nonConformal/processor patch "
216  << patchi << " " << pp.name() << endl;
217  }
218  else if (findStrings(excludePatches, pp.name()))
219  {
220  Info<< " excluding patch " << patchi
221  << " " << pp.name() << endl;
222  }
223  else
224  {
225  patchIDs.append(patchi);
226  Info<< " patch " << patchi << " " << pp.name() << endl;
227  }
228  }
229 
230  return patchIDs.shrink();
231 }
232 
233 
234 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
235 
236 int main(int argc, char *argv[])
237 {
239  (
240  "legacy VTK file format writer"
241  );
243 
244  #include "addRegionOption.H"
246  (
247  "fields",
248  "wordList",
249  "only convert the specified fields - eg '(p T U)'"
250  );
252  (
253  "cellSet",
254  "name",
255  "convert a mesh subset corresponding to the specified cellSet"
256  );
258  (
259  "faceSet",
260  "name",
261  "restrict conversion to the specified faceSet"
262  );
264  (
265  "pointSet",
266  "name",
267  "restrict conversion to the specified pointSet"
268  );
270  (
271  "ascii",
272  "write in ASCII format instead of binary"
273  );
275  (
276  "polyhedra",
277  "types",
278  "cell types to write as polyhedra - 'none', 'polyhedra', or 'all'"
279  );
281  (
282  "surfaceFields",
283  "write surfaceScalarFields (e.g., phi)"
284  );
286  (
287  "nearCellValue",
288  "use cell value on patches instead of patch value itself"
289  );
291  (
292  "noInternal",
293  "do not generate file for mesh, only for patches"
294  );
296  (
297  "noPointValues",
298  "no pointFields"
299  );
301  (
302  "allPatches",
303  "combine all patches into a single file"
304  );
306  (
307  "excludePatches",
308  "wordReList",
309  "a list of patches to exclude - eg '( inlet \".*Wall\" )' "
310  );
312  (
313  "noFaceZones",
314  "no faceZones"
315  );
317  (
318  "noLinks",
319  "don't link processor VTK files - parallel only"
320  );
322  (
323  "useTimeName",
324  "use the time name instead of the time index when naming the files"
325  );
326 
327  #include "setRootCase.H"
328  #include "createTime.H"
329 
330  const bool doWriteInternal = !args.optionFound("noInternal");
331  const bool doFaceZones = !args.optionFound("noFaceZones");
332  const bool doLinks = !args.optionFound("noLinks");
333  bool binary = !args.optionFound("ascii");
334  const bool useTimeName = args.optionFound("useTimeName");
335  const vtkTopo::vtkPolyhedra polyhedra =
337  [
339  (
340  "polyhedra",
342  )
343  ];
344 
345  if (binary && (sizeof(floatScalar) != 4 || sizeof(label) != 4))
346  {
348  << "Using ASCII rather than binary VTK format because "
349  "floatScalar and/or label are not 4 bytes in size."
350  << nl << endl;
351  binary = false;
352  }
353 
354  const bool nearCellValue = args.optionFound("nearCellValue");
355 
356  if (nearCellValue)
357  {
359  << "Using neighbouring cell value instead of patch value"
360  << nl << endl;
361  }
362 
363  const bool noPointValues = args.optionFound("noPointValues");
364 
365  if (noPointValues)
366  {
368  << "Outputting cell values only" << nl << endl;
369  }
370 
371  const bool allPatches = args.optionFound("allPatches");
372 
373  List<wordRe> excludePatches;
374  if (args.optionFound("excludePatches"))
375  {
376  args.optionLookup("excludePatches")() >> excludePatches;
377 
378  Info<< "Not including patches " << excludePatches << nl << endl;
379  }
380 
381  word cellSetName;
382  word faceSetName;
383  word pointSetName;
384  string vtkName = runTime.caseName();
385 
386  if (args.optionReadIfPresent("cellSet", cellSetName))
387  {
388  vtkName = cellSetName;
389  }
390  else if (Pstream::parRun())
391  {
392  // Strip off leading casename, leaving just processor_DDD ending.
393  string::size_type i = vtkName.rfind("processor");
394 
395  if (i != string::npos)
396  {
397  vtkName = vtkName.substr(i);
398  }
399  }
400  args.optionReadIfPresent("faceSet", faceSetName);
401  args.optionReadIfPresent("pointSet", pointSetName);
402 
403 
404 
406 
408 
409  // VTK/ directory in the case
410  fileName fvPath(runTime.path()/"VTK");
411 
412  // Directory of mesh (region0 gets filtered out)
413  fileName regionPrefix = "";
414 
416  {
417  fvPath = fvPath/regionName;
418  regionPrefix = regionName;
419  }
420 
421  if (isDir(fvPath))
422  {
423  if
424  (
425  args.optionFound("time")
426  || args.optionFound("latestTime")
427  || cellSetName.size()
428  || faceSetName.size()
429  || pointSetName.size()
431  )
432  {
433  Info<< "Keeping old VTK files in " << fvPath << nl << endl;
434  }
435  else
436  {
437  Info<< "Deleting old VTK files in " << fvPath << nl << endl;
438 
439  rmDir(fvPath);
440  }
441  }
442 
443  mkDir(fvPath);
444 
445 
446  // Mesh wrapper; does subsetting and decomposition
447  vtkMesh vMesh(mesh, polyhedra, cellSetName);
448 
449 
450  // Scan for all possible lagrangian clouds
451  HashSet<fileName> allCloudDirs;
452  forAll(timeDirs, timeI)
453  {
454  runTime.setTime(timeDirs[timeI], timeI);
455  fileNameList cloudDirs
456  (
457  readDir
458  (
459  runTime.timePath()/regionPrefix/cloud::prefix,
461  )
462  );
463  forAll(cloudDirs, i)
464  {
465  IOobjectList sprayObjs
466  (
467  mesh,
468  runTime.name(),
469  cloud::prefix/cloudDirs[i]
470  );
471 
472  IOobject* positionsPtr = sprayObjs.lookup(word("positions"));
473 
474  if (positionsPtr)
475  {
476  if (allCloudDirs.insert(cloudDirs[i]))
477  {
478  Info<< "At time: " << runTime.name()
479  << " detected cloud directory : " << cloudDirs[i]
480  << endl;
481  }
482  }
483  }
484  }
485 
486 
487  forAll(timeDirs, timeI)
488  {
489  runTime.setTime(timeDirs[timeI], timeI);
490 
491  Info<< "Time: " << runTime.name() << endl;
492 
493  word timeDesc =
494  useTimeName ? runTime.name() : Foam::name(runTime.timeIndex());
495 
496  // Check for new polyMesh/ and update mesh, fvMeshSubset and cell
497  // decomposition.
498  fvMesh::readUpdateState meshState = vMesh.readUpdate();
499 
500  const fvMesh& mesh = vMesh.mesh();
501 
502  if (meshState >= fvMesh::TOPO_CHANGE)
503  {
504  Info<< " Read new mesh" << nl << endl;
505  }
506 
507  // If faceSet: write faceSet only (as polydata)
508  if (faceSetName.size())
509  {
510  // Load the faceSet
511  faceSet set(mesh, faceSetName);
512 
513  // Filename as if patch with same name.
514  mkDir(fvPath/set.name());
515 
516  fileName patchFileName
517  (
518  fvPath/set.name()/set.name()
519  + "_"
520  + timeDesc
521  + ".vtk"
522  );
523 
524  Info<< " FaceSet : " << patchFileName << endl;
525 
526  writeFaceSet(binary, vMesh, set, patchFileName);
527 
528  continue;
529  }
530 
531  // If pointSet: write pointSet only (as polydata)
532  if (pointSetName.size())
533  {
534  // Load the pointSet
535  pointSet set(mesh, pointSetName);
536 
537  // Filename as if patch with same name.
538  mkDir(fvPath/set.name());
539 
540  fileName patchFileName
541  (
542  fvPath/set.name()/set.name()
543  + "_"
544  + timeDesc
545  + ".vtk"
546  );
547 
548  Info<< " pointSet : " << patchFileName << endl;
549 
550  writePointSet(binary, vMesh, set, patchFileName);
551 
552  continue;
553  }
554 
555 
556  // Search for list of objects for this time
557  IOobjectList objects(mesh, runTime.name());
558 
559  HashSet<word> selectedFields;
560  bool specifiedFields = args.optionReadIfPresent
561  (
562  "fields",
563  selectedFields
564  );
565 
566 
567  // Construct the vol internal fields (on the original mesh if subsetted)
568 
574 
575  if (!specifiedFields || selectedFields.size())
576  {
577  readFields(vMesh, vMesh.baseMesh(), objects, selectedFields, visf);
578  print(" volScalarField::Internal :", Info, visf);
579 
580  readFields(vMesh, vMesh.baseMesh(), objects, selectedFields, vivf);
581  print(" volVectorField::Internal :", Info, vivf);
582 
583  readFields
584  (
585  vMesh,
586  vMesh.baseMesh(),
587  objects,
588  selectedFields,
589  visptf
590  );
591  print(" volSphericalTensorFields::Internal :", Info, visptf);
592 
593  readFields
594  (
595  vMesh,
596  vMesh.baseMesh(),
597  objects,
598  selectedFields,
599  visytf
600  );
601  print(" volSymmTensorFields::Internal :", Info, visytf);
602 
603  readFields(vMesh, vMesh.baseMesh(), objects, selectedFields, vitf);
604  print(" volTensorFields::Internal :", Info, vitf);
605  }
606 
607  label nVolInternalFields =
608  visf.size()
609  + vivf.size()
610  + visptf.size()
611  + visytf.size()
612  + vitf.size();
613 
614 
615  // Construct the vol fields (on the original mesh if subsetted)
616 
622 
623  if (!specifiedFields || selectedFields.size())
624  {
625  readFields(vMesh, vMesh.baseMesh(), objects, selectedFields, vsf);
626  print(" volScalarFields :", Info, vsf);
627 
628  readFields(vMesh, vMesh.baseMesh(), objects, selectedFields, vvf);
629  print(" volVectorFields :", Info, vvf);
630 
631  readFields
632  (
633  vMesh,
634  vMesh.baseMesh(),
635  objects,
636  selectedFields,
637  vsptf
638  );
639  print(" volSphericalTensorFields :", Info, vsptf);
640 
641  readFields
642  (
643  vMesh,
644  vMesh.baseMesh(),
645  objects,
646  selectedFields,
647  vsytf
648  );
649  print(" volSymmTensorFields :", Info, vsytf);
650 
651  readFields(vMesh, vMesh.baseMesh(), objects, selectedFields, vtf);
652  print(" volTensorFields :", Info, vtf);
653  }
654 
655  label nVolFields =
656  vsf.size()
657  + vvf.size()
658  + vsptf.size()
659  + vsytf.size()
660  + vtf.size();
661 
662 
663  // Construct pointMesh only if requested
664  if (noPointValues)
665  {
666  Info<< " pointScalarFields : switched off"
667  << " (\"-noPointValues\" (at your option)\n";
668  Info<< " pointVectorFields : switched off"
669  << " (\"-noPointValues\" (at your option)\n";
670  }
671 
677 
678  if (!noPointValues && !(specifiedFields && selectedFields.empty()))
679  {
680  readFields
681  (
682  vMesh,
683  pointMesh::New(vMesh.baseMesh()),
684  objects,
685  selectedFields,
686  psf
687  );
688  print(" pointScalarFields :", Info, psf);
689 
690  readFields
691  (
692  vMesh,
693  pointMesh::New(vMesh.baseMesh()),
694  objects,
695  selectedFields,
696  pvf
697  );
698  print(" pointVectorFields :", Info, pvf);
699 
700  readFields
701  (
702  vMesh,
703  pointMesh::New(vMesh.baseMesh()),
704  objects,
705  selectedFields,
706  psptf
707  );
708  print(" pointSphericalTensorFields :", Info, psptf);
709 
710  readFields
711  (
712  vMesh,
713  pointMesh::New(vMesh.baseMesh()),
714  objects,
715  selectedFields,
716  psytf
717  );
718  print(" pointSymmTensorFields :", Info, psytf);
719 
720  readFields
721  (
722  vMesh,
723  pointMesh::New(vMesh.baseMesh()),
724  objects,
725  selectedFields,
726  ptf
727  );
728  print(" pointTensorFields :", Info, ptf);
729  }
730  Info<< endl;
731 
732  label nPointFields =
733  psf.size()
734  + pvf.size()
735  + psptf.size()
736  + psytf.size()
737  + ptf.size();
738 
739  if (doWriteInternal)
740  {
741  // Create file and write header
742  fileName vtkFileName
743  (
744  fvPath/vtkName
745  + "_"
746  + timeDesc
747  + ".vtk"
748  );
749 
750  Info<< " Internal : " << vtkFileName << endl;
751 
752  // Write mesh
753  internalWriter writer(vMesh, binary, vtkFileName);
754 
755  // cellID + volFields::Internal + VolFields
757  (
758  writer.os(),
759  vMesh.nFieldCells(),
760  1 + nVolInternalFields + nVolFields
761  );
762 
763  // Write cellID field
764  writer.writeCellIndices();
765 
766  // Write volFields::Internal
767  writer.write(visf);
768  writer.write(vivf);
769  writer.write(visptf);
770  writer.write(visytf);
771  writer.write(vitf);
772 
773  // Write volFields
774  writer.write(vsf);
775  writer.write(vvf);
776  writer.write(vsptf);
777  writer.write(vsytf);
778  writer.write(vtf);
779 
780  if (!noPointValues)
781  {
783  (
784  writer.os(),
785  vMesh.nFieldPoints(),
786  nVolFields + nPointFields
787  );
788 
789  // pointFields
790  writer.write(psf);
791  writer.write(pvf);
792  writer.write(psptf);
793  writer.write(psytf);
794  writer.write(ptf);
795 
796  // Interpolated volFields
797  const volPointInterpolation& pInterp
798  (
800  );
801 
802  writer.write(pInterp, vsf);
803  writer.write(pInterp, vvf);
804  writer.write(pInterp, vsptf);
805  writer.write(pInterp, vsytf);
806  writer.write(pInterp, vtf);
807  }
808  }
809 
810  //---------------------------------------------------------------------
811  //
812  // Write surface fields
813  //
814  //---------------------------------------------------------------------
815 
816  if (args.optionFound("surfaceFields"))
817  {
819  readFields
820  (
821  vMesh,
822  vMesh.baseMesh(),
823  objects,
824  selectedFields,
825  ssf
826  );
827  print(" surfScalarFields :", Info, ssf);
828 
830  readFields
831  (
832  vMesh,
833  vMesh.baseMesh(),
834  objects,
835  selectedFields,
836  svf
837  );
838  print(" surfVectorFields :", Info, svf);
839 
840  if (ssf.size() + svf.size() > 0)
841  {
842  // Rework the scalar fields into vectorfields.
843  label sz = svf.size();
844 
845  svf.setSize(sz + ssf.size());
846 
847  surfaceVectorField n(mesh.Sf()/mesh.magSf());
848 
849  forAll(ssf, i)
850  {
851  surfaceVectorField* ssfiPtr = (ssf[i]*n).ptr();
852  ssfiPtr->rename(ssf[i].name());
853  svf.set(sz+i, ssfiPtr);
854  }
855  ssf.clear();
856 
857  mkDir(fvPath / "surfaceFields");
858 
859  fileName surfFileName
860  (
861  fvPath
862  /"surfaceFields"
863  /"surfaceFields"
864  + "_"
865  + timeDesc
866  + ".vtk"
867  );
868 
870  (
871  binary,
872  vMesh,
873  surfFileName,
874  svf
875  );
876  }
877  }
878 
879 
880  //---------------------------------------------------------------------
881  //
882  // Write patches (POLYDATA file, one for each patch)
883  //
884  //---------------------------------------------------------------------
885 
886  const polyBoundaryMesh& patches = mesh.boundaryMesh();
887 
888  const labelList patchIDs(getSelectedPatches(patches, excludePatches));
889 
890  if (allPatches)
891  {
892  mkDir(fvPath/"allPatches");
893 
894  fileName patchFileName;
895 
896  if (vMesh.useSubMesh())
897  {
898  patchFileName =
899  fvPath/"allPatches"/cellSetName
900  + "_"
901  + timeDesc
902  + ".vtk";
903  }
904  else
905  {
906  patchFileName =
907  fvPath/"allPatches"/"allPatches"
908  + "_"
909  + timeDesc
910  + ".vtk";
911  }
912 
913  Info<< " Combined patches : " << patchFileName << endl;
914 
915  patchWriter writer
916  (
917  vMesh,
918  binary,
919  nearCellValue,
920  patchFileName,
921  getSelectedPatches(patches, excludePatches)
922  );
923 
924  // VolFields + patchID
926  (
927  writer.os(),
928  writer.nFaces(),
929  1 + nVolFields
930  );
931 
932  // Write patchID field
933  writer.writePatchIndices();
934 
935  // Write volFields
936  writer.write(vsf);
937  writer.write(vvf);
938  writer.write(vsptf);
939  writer.write(vsytf);
940  writer.write(vtf);
941 
942  if (!noPointValues)
943  {
945  (
946  writer.os(),
947  writer.nPoints(),
948  nPointFields
949  );
950 
951  // Write pointFields
952  writer.write(psf);
953  writer.write(pvf);
954  writer.write(psptf);
955  writer.write(psytf);
956  writer.write(ptf);
957  }
958  }
959  else
960  {
961  forAll(patchIDs, i)
962  {
963  const polyPatch& pp = patches[patchIDs[i]];
964 
965  mkDir(fvPath/pp.name());
966 
967  fileName patchFileName;
968 
969  if (vMesh.useSubMesh())
970  {
971  patchFileName =
972  fvPath/pp.name()/cellSetName
973  + "_"
974  + timeDesc
975  + ".vtk";
976  }
977  else
978  {
979  patchFileName =
980  fvPath/pp.name()/pp.name()
981  + "_"
982  + timeDesc
983  + ".vtk";
984  }
985 
986  Info<< " Patch : " << patchFileName << endl;
987 
988  patchWriter writer
989  (
990  vMesh,
991  binary,
992  nearCellValue,
993  patchFileName,
994  labelList(1, patchIDs[i])
995  );
996 
997  // VolFields + patchID
999  (
1000  writer.os(),
1001  writer.nFaces(),
1002  1 + nVolFields
1003  );
1004 
1005  // Write patchID field
1006  writer.writePatchIndices();
1007 
1008  // Write volFields
1009  writer.write(vsf);
1010  writer.write(vvf);
1011  writer.write(vsptf);
1012  writer.write(vsytf);
1013  writer.write(vtf);
1014 
1015  if (!noPointValues)
1016  {
1018  (
1019  writer.os(),
1020  writer.nPoints(),
1021  nVolFields
1022  + nPointFields
1023  );
1024 
1025  // Write pointFields
1026  writer.write(psf);
1027  writer.write(pvf);
1028  writer.write(psptf);
1029  writer.write(psytf);
1030  writer.write(ptf);
1031 
1033  (
1034  pp
1035  );
1036 
1037  // Write interpolated volFields
1038  writer.write(pInter, vsf);
1039  writer.write(pInter, vvf);
1040  writer.write(pInter, vsptf);
1041  writer.write(pInter, vsytf);
1042  writer.write(pInter, vtf);
1043  }
1044  }
1045  }
1046 
1047  //---------------------------------------------------------------------
1048  //
1049  // Write faceZones (POLYDATA file, one for each zone)
1050  //
1051  //---------------------------------------------------------------------
1052 
1053  if (doFaceZones)
1054  {
1056  readFields
1057  (
1058  vMesh,
1059  vMesh.baseMesh(),
1060  objects,
1061  selectedFields,
1062  ssf
1063  );
1064  print(" surfScalarFields :", Info, ssf);
1065 
1067  readFields
1068  (
1069  vMesh,
1070  vMesh.baseMesh(),
1071  objects,
1072  selectedFields,
1073  svf
1074  );
1075  print(" surfVectorFields :", Info, svf);
1076 
1077  const faceZoneList& zones = mesh.faceZones();
1078 
1079  forAll(zones, zoneI)
1080  {
1081  const faceZone& fz = zones[zoneI];
1082 
1083  mkDir(fvPath/fz.name());
1084 
1085  fileName patchFileName;
1086 
1087  if (vMesh.useSubMesh())
1088  {
1089  patchFileName =
1090  fvPath/fz.name()/cellSetName
1091  + "_"
1092  + timeDesc
1093  + ".vtk";
1094  }
1095  else
1096  {
1097  patchFileName =
1098  fvPath/fz.name()/fz.name()
1099  + "_"
1100  + timeDesc
1101  + ".vtk";
1102  }
1103 
1104  Info<< " FaceZone : " << patchFileName << endl;
1105 
1107  (
1108  IndirectList<face>(mesh.faces(), fz),
1109  mesh.points()
1110  );
1111 
1112  surfaceMeshWriter writer
1113  (
1114  binary,
1115  pp,
1116  fz.name(),
1117  patchFileName
1118  );
1119 
1120  // Number of fields
1122  (
1123  writer.os(),
1124  pp.size(),
1125  ssf.size() + svf.size()
1126  );
1127 
1128  writer.write(ssf);
1129  writer.write(svf);
1130  }
1131  }
1132 
1133 
1134 
1135  //---------------------------------------------------------------------
1136  //
1137  // Write lagrangian data
1138  //
1139  //---------------------------------------------------------------------
1140 
1141  forAllConstIter(HashSet<fileName>, allCloudDirs, iter)
1142  {
1143  const fileName& cloudName = iter.key();
1144 
1145  // Always create the cloud directory.
1146  mkDir(fvPath/cloud::prefix/cloudName);
1147 
1148  fileName lagrFileName
1149  (
1151  + "_" + timeDesc + ".vtk"
1152  );
1153 
1154  Info<< " Lagrangian: " << lagrFileName << endl;
1155 
1156 
1157  IOobjectList sprayObjs
1158  (
1159  mesh,
1160  runTime.name(),
1162  );
1163 
1164  IOobject* positionsPtr = sprayObjs.lookup(word("positions"));
1165 
1166  if (positionsPtr)
1167  {
1168  wordList labelNames(sprayObjs.names(labelIOField::typeName));
1169  Info<< " labels :";
1170  print(Info, labelNames);
1171 
1172  wordList scalarNames(sprayObjs.names(scalarIOField::typeName));
1173  Info<< " scalars :";
1174  print(Info, scalarNames);
1175 
1176  wordList vectorNames(sprayObjs.names(vectorIOField::typeName));
1177  Info<< " vectors :";
1178  print(Info, vectorNames);
1179 
1180  wordList sphereNames
1181  (
1182  sprayObjs.names
1183  (
1185  )
1186  );
1187  Info<< " spherical tensors :";
1188  print(Info, sphereNames);
1189 
1190  wordList symmNames
1191  (
1192  sprayObjs.names
1193  (
1195  )
1196  );
1197  Info<< " symm tensors :";
1198  print(Info, symmNames);
1199 
1200  wordList tensorNames(sprayObjs.names(tensorIOField::typeName));
1201  Info<< " tensors :";
1202  print(Info, tensorNames);
1203 
1204  lagrangianWriter writer
1205  (
1206  vMesh,
1207  binary,
1208  lagrFileName,
1209  cloudName,
1210  false
1211  );
1212 
1213  // Write number of fields
1214  writer.writeFieldsHeader
1215  (
1216  labelNames.size()
1217  + scalarNames.size()
1218  + vectorNames.size()
1219  + sphereNames.size()
1220  + symmNames.size()
1221  + tensorNames.size()
1222  );
1223 
1224  // Fields
1225  writer.writeIOField<label>(labelNames);
1226  writer.writeIOField<scalar>(scalarNames);
1227  writer.writeIOField<vector>(vectorNames);
1228  writer.writeIOField<sphericalTensor>(sphereNames);
1229  writer.writeIOField<symmTensor>(symmNames);
1230  writer.writeIOField<tensor>(tensorNames);
1231  }
1232  else
1233  {
1234  lagrangianWriter writer
1235  (
1236  vMesh,
1237  binary,
1238  lagrFileName,
1239  cloudName,
1240  true
1241  );
1242 
1243  // Write number of fields
1244  writer.writeFieldsHeader(0);
1245  }
1246  }
1247  }
1248 
1249 
1250  //---------------------------------------------------------------------
1251  //
1252  // Link parallel outputs back to undecomposed case for ease of loading
1253  //
1254  //---------------------------------------------------------------------
1255 
1256  if (Pstream::parRun() && doLinks)
1257  {
1258  mkDir(runTime.globalPath()/"VTK");
1259  chDir(runTime.globalPath()/"VTK");
1260 
1261  Info<< "Linking all processor files to " << runTime.globalPath()/"VTK"
1262  << endl;
1263 
1264  // Get list of vtk files
1265  fileName procVTK
1266  (
1267  fileName("..")
1268  /"processor" + name(Pstream::myProcNo())
1269  /"VTK"
1270  );
1271 
1272  fileNameList dirs(readDir(procVTK, fileType::directory));
1273  label sz = dirs.size();
1274  dirs.setSize(sz + 1);
1275  dirs[sz] = ".";
1276 
1277  forAll(dirs, i)
1278  {
1279  fileNameList subFiles(readDir(procVTK/dirs[i], fileType::file));
1280 
1281  forAll(subFiles, j)
1282  {
1283  fileName procFile(procVTK/dirs[i]/subFiles[j]);
1284 
1285  if (exists(procFile))
1286  {
1287  ln
1288  (
1289  procFile,
1290  "processor"
1291  + name(Pstream::myProcNo())
1292  + "_"
1293  + procFile.name()
1294  );
1295  }
1296  }
1297  }
1298  }
1299 
1300  Info<< "End\n" << endl;
1301 
1302  return 0;
1303 }
1304 
1305 
1306 // ************************************************************************* //
graph_traits< Graph >::vertices_size_type size_type
Definition: SloanRenumber.C:73
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
static pointMesh & New(const word &name, const polyMesh &mesh)
Construct and return the named DemandDrivenMeshObject.
static const char *const typeName
Definition: Field.H:106
Generic GeometricField class.
A HashTable with keys but without contents.
Definition: HashSet.H:62
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:109
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
bool empty() const
Return true if the hash table is empty.
Definition: HashTableI.H:72
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
A List with indirect addressing.
Definition: IndirectList.H:105
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
Interpolation class within a primitive patch. Allows interpolation from points to faces and vice vers...
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
bool set(const label) const
Is element set.
Definition: PtrListI.H:62
void clear()
Clear the PtrList, i.e. set size to zero deleting all the.
Definition: PtrList.C:198
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
Definition: PtrList.C:131
Templated 3D SphericalTensor derived from VectorSpace adding construction from 1 component,...
label size() const
Return the number of elements in the UList.
Definition: UListI.H:311
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
const word & name() const
Return name.
Definition: Zone.H:184
static void addOption(const word &opt, const string &param="", const string &usage="")
Add to an option to validOptions with usage information.
Definition: argList.C:128
static void addNote(const string &)
Add extra notes for the usage information.
Definition: argList.C:159
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
Definition: argList.C:118
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:114
bool optionReadIfPresent(const word &opt, T &) const
Read a value from the named option if present.
Definition: argListI.H:204
IStringStream optionLookup(const word &opt) const
Return an IStringStream from the named option.
Definition: argListI.H:120
T optionLookupOrDefault(const word &opt, const T &deflt) const
Read a value from the named option if present.
Definition: argListI.H:243
static const word prefix
The prefix to local: lagrangian.
Definition: cloud.H:63
A list of face labels.
Definition: faceSet.H:51
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:59
A class for handling file names.
Definition: fileName.H:82
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: fvMesh.H:108
@ TOPO_CHANGE
Definition: fvMesh.H:112
const surfaceVectorField & Sf() const
Return cell face area vectors.
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
const polyMesh & mesh() const
Return reference to polyMesh.
Definition: fvMesh.H:441
Write fields (internal).
Writes Lagrangian points and fields.
const word & name() const
Return name.
Write patch fields.
Definition: patchWriter.H:60
A set of point labels.
Definition: pointSet.H:51
Foam::polyBoundaryMesh.
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:267
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1326
const faceZoneList & faceZones() const
Return face zones.
Definition: polyMesh.H:446
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:404
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1313
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:70
virtual void rename(const word &newName)
Rename.
Definition: regIOobject.C:417
Write faces with fields.
static void addOptions(const bool constant=true, const bool withZero=false)
Add the options handled by timeSelector to argList::validOptions.
Definition: timeSelector.C:114
static instantList select0(Time &runTime, const argList &args)
Return the set of times selected based on the argList options.
Definition: timeSelector.C:252
Interpolate from cell centres to points (vertices) using inverse distance weighting.
Encapsulation of VTK mesh data. Holds mesh or meshsubset and polyhedral-cell decomposition on it.
Definition: vtkMesh.H:54
static const NamedEnum< vtkPolyhedra, 3 > vtkPolyhedraNames_
Names of groups of cell types retain as polyhedra.
Definition: vtkTopo.H:66
vtkPolyhedra
Groups of cell types retain as polyhedra.
Definition: vtkTopo.H:59
A class for handling words, derived from string.
Definition: word.H:62
Foam::faceZoneList.
int main(int argc, char *argv[])
Definition: financialFoam.C:44
label patchi
static instantList timeDirs
Definition: globalFoam.H:44
const fvPatchList & patches
#define WarningInFunction
Report a warning using Foam::Warning.
void writePointDataHeader(std::ostream &, const label nPoints, const label nFields)
Write point data header.
void writeCellDataHeader(std::ostream &, const label nCells, const label nFields)
Write cell data header.
Namespace for OpenFOAM.
void writeSurfFields(const bool binary, const vtkMesh &vMesh, const fileName &fileName, const UPtrList< const surfaceVectorField > &surfVectorFields)
void writeFaceSet(const bool binary, const vtkMesh &vMesh, const faceSet &set, const fileName &fileName)
List< label > labelList
A List of labels.
Definition: labelList.H:56
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:59
bool findStrings(const wordReListMatcher &matcher, const std::string &str)
Return true if string matches one of the regular expressions.
Definition: stringListOps.H:52
bool mkDir(const fileName &, mode_t=0777)
Make a directory and return an error if it could not be created.
Definition: POSIX.C:290
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const HashSet< word > &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the specified type.
Definition: ReadFields.C:244
const word & regionName(const solver &region)
Definition: solver.H:209
bool exists(const fileName &, const bool checkVariants=true, const bool followLink=true)
Does the name exist (as directory or file) in the file system?
Definition: POSIX.C:520
messageStream Info
bool rmDir(const fileName &)
Remove a directory and its contents.
Definition: POSIX.C:1047
float floatScalar
Float precision floating point scalar type.
Definition: floatScalar.H:52
bool isDir(const fileName &, const bool followLink=true)
Does the name exist as a directory in the file system?
Definition: POSIX.C:539
static const char nl
Definition: Ostream.H:266
fileNameList readDir(const fileName &, const fileType=fileType::file, const bool filterVariants=true, const bool followLink=true)
Read a directory and return the entries as a string list.
Definition: POSIX.C:662
bool ln(const fileName &src, const fileName &dst)
Create a softlink. dst should not exist. Returns true if successful.
Definition: POSIX.C:908
void writePointSet(const bool binary, const vtkMesh &vMesh, const pointSet &set, const fileName &fileName)
bool chDir(const fileName &dir)
Change the current directory to the one given and return true,.
Definition: POSIX.C:284
objects
Foam::argList args(argc, argv)
Operations on lists of strings.
const word cloudName(propsDict.lookup("cloudName"))
Write faceSet to vtk polydata file. Only one data which is original faceID.
Write pointSet to vtk polydata file. Only one data which is original pointID.
Write a patch with its data.