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-2022 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 poly
82  write polyhedral cells without tet/pyramid decomposition
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 "fvCFD.H"
137 #include "pointMesh.H"
138 #include "volPointInterpolation.H"
139 #include "emptyPolyPatch.H"
140 #include "nonConformalPolyPatch.H"
141 #include "labelIOField.H"
142 #include "scalarIOField.H"
143 #include "sphericalTensorIOField.H"
144 #include "symmTensorIOField.H"
145 #include "tensorIOField.H"
146 #include "meshFaceZones.H"
147 #include "Cloud.H"
148 #include "passiveParticle.H"
149 #include "stringListOps.H"
150 
151 #include "vtkMesh.H"
152 #include "readFields.H"
153 #include "vtkWriteOps.H"
154 
155 #include "internalWriter.H"
156 #include "patchWriter.H"
157 #include "lagrangianWriter.H"
158 
159 #include "writeFaceSet.H"
160 #include "writePointSet.H"
161 #include "surfaceMeshWriter.H"
162 #include "writeSurfFields.H"
163 
164 
165 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
166 
167 template<class GeoField>
168 void print(const char* msg, Ostream& os, const PtrList<const GeoField>& flds)
169 {
170  if (flds.size())
171  {
172  os << msg;
173  forAll(flds, i)
174  {
175  os << ' ' << flds[i].name();
176  }
177  os << endl;
178  }
179 }
180 
181 
182 void print(Ostream& os, const wordList& flds)
183 {
184  forAll(flds, i)
185  {
186  os << ' ' << flds[i];
187  }
188  os << endl;
189 }
190 
191 
192 labelList getSelectedPatches
193 (
194  const polyBoundaryMesh& patches,
195  const List<wordRe>& excludePatches
196 )
197 {
198  DynamicList<label> patchIDs(patches.size());
199 
200  Info<< "Combining patches:" << endl;
201 
202  forAll(patches, patchi)
203  {
204  const polyPatch& pp = patches[patchi];
205 
206  if
207  (
208  isA<emptyPolyPatch>(pp)
209  || isA<nonConformalPolyPatch>(pp)
210  || (Pstream::parRun() && isType<processorPolyPatch>(pp))
211  )
212  {
213  Info<< " discarding empty/nonConformal/processor patch "
214  << patchi << " " << pp.name() << endl;
215  }
216  else if (findStrings(excludePatches, pp.name()))
217  {
218  Info<< " excluding patch " << patchi
219  << " " << pp.name() << endl;
220  }
221  else
222  {
223  patchIDs.append(patchi);
224  Info<< " patch " << patchi << " " << pp.name() << endl;
225  }
226  }
227 
228  return patchIDs.shrink();
229 }
230 
231 
232 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
233 
234 int main(int argc, char *argv[])
235 {
236  argList::addNote
237  (
238  "legacy VTK file format writer"
239  );
240  timeSelector::addOptions();
241 
242  #include "addRegionOption.H"
243  argList::addOption
244  (
245  "fields",
246  "wordList",
247  "only convert the specified fields - eg '(p T U)'"
248  );
249  argList::addOption
250  (
251  "cellSet",
252  "name",
253  "convert a mesh subset corresponding to the specified cellSet"
254  );
255  argList::addOption
256  (
257  "faceSet",
258  "name",
259  "restrict conversion to the specified faceSet"
260  );
261  argList::addOption
262  (
263  "pointSet",
264  "name",
265  "restrict conversion to the specified pointSet"
266  );
267  argList::addBoolOption
268  (
269  "ascii",
270  "write in ASCII format instead of binary"
271  );
272  argList::addBoolOption
273  (
274  "poly",
275  "write polyhedral cells without tet/pyramid decomposition"
276  );
277  argList::addBoolOption
278  (
279  "surfaceFields",
280  "write surfaceScalarFields (e.g., phi)"
281  );
282  argList::addBoolOption
283  (
284  "nearCellValue",
285  "use cell value on patches instead of patch value itself"
286  );
287  argList::addBoolOption
288  (
289  "noInternal",
290  "do not generate file for mesh, only for patches"
291  );
292  argList::addBoolOption
293  (
294  "noPointValues",
295  "no pointFields"
296  );
297  argList::addBoolOption
298  (
299  "allPatches",
300  "combine all patches into a single file"
301  );
302  argList::addOption
303  (
304  "excludePatches",
305  "wordReList",
306  "a list of patches to exclude - eg '( inlet \".*Wall\" )' "
307  );
308  argList::addBoolOption
309  (
310  "noFaceZones",
311  "no faceZones"
312  );
313  argList::addBoolOption
314  (
315  "noLinks",
316  "don't link processor VTK files - parallel only"
317  );
318  argList::addBoolOption
319  (
320  "useTimeName",
321  "use the time name instead of the time index when naming the files"
322  );
323 
324  #include "setRootCase.H"
325  #include "createTime.H"
326 
327  const bool doWriteInternal = !args.optionFound("noInternal");
328  const bool doFaceZones = !args.optionFound("noFaceZones");
329  const bool doLinks = !args.optionFound("noLinks");
330  bool binary = !args.optionFound("ascii");
331  const bool useTimeName = args.optionFound("useTimeName");
332 
333  // Decomposition of polyhedral cells into tets/pyramids cells
334  vtkTopo::decomposePoly = !args.optionFound("poly");
335 
336  if (binary && (sizeof(floatScalar) != 4 || sizeof(label) != 4))
337  {
339  << "Using ASCII rather than binary VTK format because "
340  "floatScalar and/or label are not 4 bytes in size."
341  << nl << endl;
342  binary = false;
343  }
344 
345  const bool nearCellValue = args.optionFound("nearCellValue");
346 
347  if (nearCellValue)
348  {
350  << "Using neighbouring cell value instead of patch value"
351  << nl << endl;
352  }
353 
354  const bool noPointValues = args.optionFound("noPointValues");
355 
356  if (noPointValues)
357  {
359  << "Outputting cell values only" << nl << endl;
360  }
361 
362  const bool allPatches = args.optionFound("allPatches");
363 
364  List<wordRe> excludePatches;
365  if (args.optionFound("excludePatches"))
366  {
367  args.optionLookup("excludePatches")() >> excludePatches;
368 
369  Info<< "Not including patches " << excludePatches << nl << endl;
370  }
371 
372  word cellSetName;
373  word faceSetName;
374  word pointSetName;
375  string vtkName = runTime.caseName();
376 
377  if (args.optionReadIfPresent("cellSet", cellSetName))
378  {
379  vtkName = cellSetName;
380  }
381  else if (Pstream::parRun())
382  {
383  // Strip off leading casename, leaving just processor_DDD ending.
384  string::size_type i = vtkName.rfind("processor");
385 
386  if (i != string::npos)
387  {
388  vtkName = vtkName.substr(i);
389  }
390  }
391  args.optionReadIfPresent("faceSet", faceSetName);
392  args.optionReadIfPresent("pointSet", pointSetName);
393 
394 
395 
396  instantList timeDirs = timeSelector::select0(runTime, args);
397 
398  #include "createNamedMesh.H"
399 
400  // VTK/ directory in the case
401  fileName fvPath(runTime.path()/"VTK");
402 
403  // Directory of mesh (region0 gets filtered out)
404  fileName regionPrefix = "";
405 
406  if (regionName != polyMesh::defaultRegion)
407  {
408  fvPath = fvPath/regionName;
409  regionPrefix = regionName;
410  }
411 
412  if (isDir(fvPath))
413  {
414  if
415  (
416  args.optionFound("time")
417  || args.optionFound("latestTime")
418  || cellSetName.size()
419  || faceSetName.size()
420  || pointSetName.size()
421  || regionName != polyMesh::defaultRegion
422  )
423  {
424  Info<< "Keeping old VTK files in " << fvPath << nl << endl;
425  }
426  else
427  {
428  Info<< "Deleting old VTK files in " << fvPath << nl << endl;
429 
430  rmDir(fvPath);
431  }
432  }
433 
434  mkDir(fvPath);
435 
436 
437  // Mesh wrapper; does subsetting and decomposition
438  vtkMesh vMesh(mesh, cellSetName);
439 
440 
441  // Scan for all possible lagrangian clouds
442  HashSet<fileName> allCloudDirs;
443  forAll(timeDirs, timeI)
444  {
445  runTime.setTime(timeDirs[timeI], timeI);
446  fileNameList cloudDirs
447  (
448  readDir
449  (
450  runTime.timePath()/regionPrefix/cloud::prefix,
451  fileType::directory
452  )
453  );
454  forAll(cloudDirs, i)
455  {
456  IOobjectList sprayObjs
457  (
458  mesh,
459  runTime.timeName(),
460  cloud::prefix/cloudDirs[i]
461  );
462 
463  IOobject* positionsPtr = sprayObjs.lookup(word("positions"));
464 
465  if (positionsPtr)
466  {
467  if (allCloudDirs.insert(cloudDirs[i]))
468  {
469  Info<< "At time: " << runTime.timeName()
470  << " detected cloud directory : " << cloudDirs[i]
471  << endl;
472  }
473  }
474  }
475  }
476 
477 
478  forAll(timeDirs, timeI)
479  {
480  runTime.setTime(timeDirs[timeI], timeI);
481 
482  Info<< "Time: " << runTime.timeName() << endl;
483 
484  word timeDesc =
485  useTimeName ? runTime.timeName() : Foam::name(runTime.timeIndex());
486 
487  // Check for new polyMesh/ and update mesh, fvMeshSubset and cell
488  // decomposition.
489  polyMesh::readUpdateState meshState = vMesh.readUpdate();
490 
491  const fvMesh& mesh = vMesh.mesh();
492 
493  if
494  (
495  meshState == polyMesh::TOPO_CHANGE
496  || meshState == polyMesh::TOPO_PATCH_CHANGE
497  )
498  {
499  Info<< " Read new mesh" << nl << endl;
500  }
501 
502  // If faceSet: write faceSet only (as polydata)
503  if (faceSetName.size())
504  {
505  // Load the faceSet
506  faceSet set(mesh, faceSetName);
507 
508  // Filename as if patch with same name.
509  mkDir(fvPath/set.name());
510 
511  fileName patchFileName
512  (
513  fvPath/set.name()/set.name()
514  + "_"
515  + timeDesc
516  + ".vtk"
517  );
518 
519  Info<< " FaceSet : " << patchFileName << endl;
520 
521  writeFaceSet(binary, vMesh, set, patchFileName);
522 
523  continue;
524  }
525 
526  // If pointSet: write pointSet only (as polydata)
527  if (pointSetName.size())
528  {
529  // Load the pointSet
530  pointSet set(mesh, pointSetName);
531 
532  // Filename as if patch with same name.
533  mkDir(fvPath/set.name());
534 
535  fileName patchFileName
536  (
537  fvPath/set.name()/set.name()
538  + "_"
539  + timeDesc
540  + ".vtk"
541  );
542 
543  Info<< " pointSet : " << patchFileName << endl;
544 
545  writePointSet(binary, vMesh, set, patchFileName);
546 
547  continue;
548  }
549 
550 
551  // Search for list of objects for this time
552  IOobjectList objects(mesh, runTime.timeName());
553 
554  HashSet<word> selectedFields;
555  bool specifiedFields = args.optionReadIfPresent
556  (
557  "fields",
558  selectedFields
559  );
560 
561 
562  // Construct the vol internal fields (on the original mesh if subsetted)
563 
564  PtrList<const volScalarField::Internal> visf;
565  PtrList<const volVectorField::Internal> vivf;
566  PtrList<const volSphericalTensorField::Internal> visptf;
567  PtrList<const volSymmTensorField::Internal> visytf;
568  PtrList<const volTensorField::Internal> vitf;
569 
570  if (!specifiedFields || selectedFields.size())
571  {
572  readFields(vMesh, vMesh.baseMesh(), objects, selectedFields, visf);
573  print(" volScalarField::Internal :", Info, visf);
574 
575  readFields(vMesh, vMesh.baseMesh(), objects, selectedFields, vivf);
576  print(" volVectorField::Internal :", Info, vivf);
577 
578  readFields
579  (
580  vMesh,
581  vMesh.baseMesh(),
582  objects,
583  selectedFields,
584  visptf
585  );
586  print(" volSphericalTensorFields::Internal :", Info, visptf);
587 
588  readFields
589  (
590  vMesh,
591  vMesh.baseMesh(),
592  objects,
593  selectedFields,
594  visytf
595  );
596  print(" volSymmTensorFields::Internal :", Info, visytf);
597 
598  readFields(vMesh, vMesh.baseMesh(), objects, selectedFields, vitf);
599  print(" volTensorFields::Internal :", Info, vitf);
600  }
601 
602  label nVolInternalFields =
603  visf.size()
604  + vivf.size()
605  + visptf.size()
606  + visytf.size()
607  + vitf.size();
608 
609 
610  // Construct the vol fields (on the original mesh if subsetted)
611 
612  PtrList<const volScalarField> vsf;
613  PtrList<const volVectorField> vvf;
614  PtrList<const volSphericalTensorField> vsptf;
615  PtrList<const volSymmTensorField> vsytf;
616  PtrList<const volTensorField> vtf;
617 
618  if (!specifiedFields || selectedFields.size())
619  {
620  readFields(vMesh, vMesh.baseMesh(), objects, selectedFields, vsf);
621  print(" volScalarFields :", Info, vsf);
622 
623  readFields(vMesh, vMesh.baseMesh(), objects, selectedFields, vvf);
624  print(" volVectorFields :", Info, vvf);
625 
626  readFields
627  (
628  vMesh,
629  vMesh.baseMesh(),
630  objects,
631  selectedFields,
632  vsptf
633  );
634  print(" volSphericalTensorFields :", Info, vsptf);
635 
636  readFields
637  (
638  vMesh,
639  vMesh.baseMesh(),
640  objects,
641  selectedFields,
642  vsytf
643  );
644  print(" volSymmTensorFields :", Info, vsytf);
645 
646  readFields(vMesh, vMesh.baseMesh(), objects, selectedFields, vtf);
647  print(" volTensorFields :", Info, vtf);
648  }
649 
650  label nVolFields =
651  vsf.size()
652  + vvf.size()
653  + vsptf.size()
654  + vsytf.size()
655  + vtf.size();
656 
657 
658  // Construct pointMesh only if requested
659  if (noPointValues)
660  {
661  Info<< " pointScalarFields : switched off"
662  << " (\"-noPointValues\" (at your option)\n";
663  Info<< " pointVectorFields : switched off"
664  << " (\"-noPointValues\" (at your option)\n";
665  }
666 
667  PtrList<const pointScalarField> psf;
668  PtrList<const pointVectorField> pvf;
669  PtrList<const pointSphericalTensorField> psptf;
670  PtrList<const pointSymmTensorField> psytf;
671  PtrList<const pointTensorField> ptf;
672 
673  if (!noPointValues && !(specifiedFields && selectedFields.empty()))
674  {
675  readFields
676  (
677  vMesh,
678  pointMesh::New(vMesh.baseMesh()),
679  objects,
680  selectedFields,
681  psf
682  );
683  print(" pointScalarFields :", Info, psf);
684 
685  readFields
686  (
687  vMesh,
688  pointMesh::New(vMesh.baseMesh()),
689  objects,
690  selectedFields,
691  pvf
692  );
693  print(" pointVectorFields :", Info, pvf);
694 
695  readFields
696  (
697  vMesh,
698  pointMesh::New(vMesh.baseMesh()),
699  objects,
700  selectedFields,
701  psptf
702  );
703  print(" pointSphericalTensorFields :", Info, psptf);
704 
705  readFields
706  (
707  vMesh,
708  pointMesh::New(vMesh.baseMesh()),
709  objects,
710  selectedFields,
711  psytf
712  );
713  print(" pointSymmTensorFields :", Info, psytf);
714 
715  readFields
716  (
717  vMesh,
718  pointMesh::New(vMesh.baseMesh()),
719  objects,
720  selectedFields,
721  ptf
722  );
723  print(" pointTensorFields :", Info, ptf);
724  }
725  Info<< endl;
726 
727  label nPointFields =
728  psf.size()
729  + pvf.size()
730  + psptf.size()
731  + psytf.size()
732  + ptf.size();
733 
734  if (doWriteInternal)
735  {
736  // Create file and write header
737  fileName vtkFileName
738  (
739  fvPath/vtkName
740  + "_"
741  + timeDesc
742  + ".vtk"
743  );
744 
745  Info<< " Internal : " << vtkFileName << endl;
746 
747  // Write mesh
748  internalWriter writer(vMesh, binary, vtkFileName);
749 
750  // cellID + volFields::Internal + VolFields
752  (
753  writer.os(),
754  vMesh.nFieldCells(),
755  1 + nVolInternalFields + nVolFields
756  );
757 
758  // Write cellID field
759  writer.writeCellIDs();
760 
761  // Write volFields::Internal
762  writer.write(visf);
763  writer.write(vivf);
764  writer.write(visptf);
765  writer.write(visytf);
766  writer.write(vitf);
767 
768  // Write volFields
769  writer.write(vsf);
770  writer.write(vvf);
771  writer.write(vsptf);
772  writer.write(vsytf);
773  writer.write(vtf);
774 
775  if (!noPointValues)
776  {
778  (
779  writer.os(),
780  vMesh.nFieldPoints(),
781  nVolFields + nPointFields
782  );
783 
784  // pointFields
785  writer.write(psf);
786  writer.write(pvf);
787  writer.write(psptf);
788  writer.write(psytf);
789  writer.write(ptf);
790 
791  // Interpolated volFields
792  volPointInterpolation pInterp(mesh);
793  writer.write(pInterp, vsf);
794  writer.write(pInterp, vvf);
795  writer.write(pInterp, vsptf);
796  writer.write(pInterp, vsytf);
797  writer.write(pInterp, vtf);
798  }
799  }
800 
801  //---------------------------------------------------------------------
802  //
803  // Write surface fields
804  //
805  //---------------------------------------------------------------------
806 
807  if (args.optionFound("surfaceFields"))
808  {
809  PtrList<const surfaceScalarField> ssf;
810  readFields
811  (
812  vMesh,
813  vMesh.baseMesh(),
814  objects,
815  selectedFields,
816  ssf
817  );
818  print(" surfScalarFields :", Info, ssf);
819 
820  PtrList<const surfaceVectorField> svf;
821  readFields
822  (
823  vMesh,
824  vMesh.baseMesh(),
825  objects,
826  selectedFields,
827  svf
828  );
829  print(" surfVectorFields :", Info, svf);
830 
831  if (ssf.size() + svf.size() > 0)
832  {
833  // Rework the scalar fields into vectorfields.
834  label sz = svf.size();
835 
836  svf.setSize(sz + ssf.size());
837 
838  surfaceVectorField n(mesh.Sf()/mesh.magSf());
839 
840  forAll(ssf, i)
841  {
842  surfaceVectorField* ssfiPtr = (ssf[i]*n).ptr();
843  ssfiPtr->rename(ssf[i].name());
844  svf.set(sz+i, ssfiPtr);
845  }
846  ssf.clear();
847 
848  mkDir(fvPath / "surfaceFields");
849 
850  fileName surfFileName
851  (
852  fvPath
853  /"surfaceFields"
854  /"surfaceFields"
855  + "_"
856  + timeDesc
857  + ".vtk"
858  );
859 
861  (
862  binary,
863  vMesh,
864  surfFileName,
865  svf
866  );
867  }
868  }
869 
870 
871  //---------------------------------------------------------------------
872  //
873  // Write patches (POLYDATA file, one for each patch)
874  //
875  //---------------------------------------------------------------------
876 
877  const polyBoundaryMesh& patches = mesh.boundaryMesh();
878 
879  const labelList patchIDs(getSelectedPatches(patches, excludePatches));
880 
881  if (allPatches)
882  {
883  mkDir(fvPath/"allPatches");
884 
885  fileName patchFileName;
886 
887  if (vMesh.useSubMesh())
888  {
889  patchFileName =
890  fvPath/"allPatches"/cellSetName
891  + "_"
892  + timeDesc
893  + ".vtk";
894  }
895  else
896  {
897  patchFileName =
898  fvPath/"allPatches"/"allPatches"
899  + "_"
900  + timeDesc
901  + ".vtk";
902  }
903 
904  Info<< " Combined patches : " << patchFileName << endl;
905 
906  patchWriter writer
907  (
908  vMesh,
909  binary,
910  nearCellValue,
911  patchFileName,
912  getSelectedPatches(patches, excludePatches)
913  );
914 
915  // VolFields + patchID
917  (
918  writer.os(),
919  writer.nFaces(),
920  1 + nVolFields
921  );
922 
923  // Write patchID field
924  writer.writePatchIDs();
925 
926  // Write volFields
927  writer.write(vsf);
928  writer.write(vvf);
929  writer.write(vsptf);
930  writer.write(vsytf);
931  writer.write(vtf);
932 
933  if (!noPointValues)
934  {
936  (
937  writer.os(),
938  writer.nPoints(),
939  nPointFields
940  );
941 
942  // Write pointFields
943  writer.write(psf);
944  writer.write(pvf);
945  writer.write(psptf);
946  writer.write(psytf);
947  writer.write(ptf);
948  }
949  }
950  else
951  {
952  forAll(patchIDs, i)
953  {
954  const polyPatch& pp = patches[patchIDs[i]];
955 
956  mkDir(fvPath/pp.name());
957 
958  fileName patchFileName;
959 
960  if (vMesh.useSubMesh())
961  {
962  patchFileName =
963  fvPath/pp.name()/cellSetName
964  + "_"
965  + timeDesc
966  + ".vtk";
967  }
968  else
969  {
970  patchFileName =
971  fvPath/pp.name()/pp.name()
972  + "_"
973  + timeDesc
974  + ".vtk";
975  }
976 
977  Info<< " Patch : " << patchFileName << endl;
978 
979  patchWriter writer
980  (
981  vMesh,
982  binary,
983  nearCellValue,
984  patchFileName,
985  labelList(1, patchIDs[i])
986  );
987 
988  // VolFields + patchID
990  (
991  writer.os(),
992  writer.nFaces(),
993  1 + nVolFields
994  );
995 
996  // Write patchID field
997  writer.writePatchIDs();
998 
999  // Write volFields
1000  writer.write(vsf);
1001  writer.write(vvf);
1002  writer.write(vsptf);
1003  writer.write(vsytf);
1004  writer.write(vtf);
1005 
1006  if (!noPointValues)
1007  {
1009  (
1010  writer.os(),
1011  writer.nPoints(),
1012  nVolFields
1013  + nPointFields
1014  );
1015 
1016  // Write pointFields
1017  writer.write(psf);
1018  writer.write(pvf);
1019  writer.write(psptf);
1020  writer.write(psytf);
1021  writer.write(ptf);
1022 
1023  PrimitivePatchInterpolation<primitivePatch> pInter
1024  (
1025  pp
1026  );
1027 
1028  // Write interpolated volFields
1029  writer.write(pInter, vsf);
1030  writer.write(pInter, vvf);
1031  writer.write(pInter, vsptf);
1032  writer.write(pInter, vsytf);
1033  writer.write(pInter, vtf);
1034  }
1035  }
1036  }
1037 
1038  //---------------------------------------------------------------------
1039  //
1040  // Write faceZones (POLYDATA file, one for each zone)
1041  //
1042  //---------------------------------------------------------------------
1043 
1044  if (doFaceZones)
1045  {
1046  PtrList<const surfaceScalarField> ssf;
1047  readFields
1048  (
1049  vMesh,
1050  vMesh.baseMesh(),
1051  objects,
1052  selectedFields,
1053  ssf
1054  );
1055  print(" surfScalarFields :", Info, ssf);
1056 
1057  PtrList<const surfaceVectorField> svf;
1058  readFields
1059  (
1060  vMesh,
1061  vMesh.baseMesh(),
1062  objects,
1063  selectedFields,
1064  svf
1065  );
1066  print(" surfVectorFields :", Info, svf);
1067 
1068  const meshFaceZones& zones = mesh.faceZones();
1069 
1070  forAll(zones, zoneI)
1071  {
1072  const faceZone& fz = zones[zoneI];
1073 
1074  mkDir(fvPath/fz.name());
1075 
1076  fileName patchFileName;
1077 
1078  if (vMesh.useSubMesh())
1079  {
1080  patchFileName =
1081  fvPath/fz.name()/cellSetName
1082  + "_"
1083  + timeDesc
1084  + ".vtk";
1085  }
1086  else
1087  {
1088  patchFileName =
1089  fvPath/fz.name()/fz.name()
1090  + "_"
1091  + timeDesc
1092  + ".vtk";
1093  }
1094 
1095  Info<< " FaceZone : " << patchFileName << endl;
1096 
1098  (
1099  IndirectList<face>(mesh.faces(), fz),
1100  mesh.points()
1101  );
1102 
1103  surfaceMeshWriter writer
1104  (
1105  binary,
1106  pp,
1107  fz.name(),
1108  patchFileName
1109  );
1110 
1111  // Number of fields
1113  (
1114  writer.os(),
1115  pp.size(),
1116  ssf.size() + svf.size()
1117  );
1118 
1119  writer.write(ssf);
1120  writer.write(svf);
1121  }
1122  }
1123 
1124 
1125 
1126  //---------------------------------------------------------------------
1127  //
1128  // Write lagrangian data
1129  //
1130  //---------------------------------------------------------------------
1131 
1132  forAllConstIter(HashSet<fileName>, allCloudDirs, iter)
1133  {
1134  const fileName& cloudName = iter.key();
1135 
1136  // Always create the cloud directory.
1137  mkDir(fvPath/cloud::prefix/cloudName);
1138 
1139  fileName lagrFileName
1140  (
1141  fvPath/cloud::prefix/cloudName/cloudName
1142  + "_" + timeDesc + ".vtk"
1143  );
1144 
1145  Info<< " Lagrangian: " << lagrFileName << endl;
1146 
1147 
1148  IOobjectList sprayObjs
1149  (
1150  mesh,
1151  runTime.timeName(),
1152  cloud::prefix/cloudName
1153  );
1154 
1155  IOobject* positionsPtr = sprayObjs.lookup(word("positions"));
1156 
1157  if (positionsPtr)
1158  {
1159  wordList labelNames(sprayObjs.names(labelIOField::typeName));
1160  Info<< " labels :";
1161  print(Info, labelNames);
1162 
1163  wordList scalarNames(sprayObjs.names(scalarIOField::typeName));
1164  Info<< " scalars :";
1165  print(Info, scalarNames);
1166 
1167  wordList vectorNames(sprayObjs.names(vectorIOField::typeName));
1168  Info<< " vectors :";
1169  print(Info, vectorNames);
1170 
1171  wordList sphereNames
1172  (
1173  sprayObjs.names
1174  (
1175  sphericalTensorIOField::typeName
1176  )
1177  );
1178  Info<< " spherical tensors :";
1179  print(Info, sphereNames);
1180 
1181  wordList symmNames
1182  (
1183  sprayObjs.names
1184  (
1185  symmTensorIOField::typeName
1186  )
1187  );
1188  Info<< " symm tensors :";
1189  print(Info, symmNames);
1190 
1191  wordList tensorNames(sprayObjs.names(tensorIOField::typeName));
1192  Info<< " tensors :";
1193  print(Info, tensorNames);
1194 
1195  lagrangianWriter writer
1196  (
1197  vMesh,
1198  binary,
1199  lagrFileName,
1200  cloudName,
1201  false
1202  );
1203 
1204  // Write number of fields
1205  writer.writeFieldsHeader
1206  (
1207  labelNames.size()
1208  + scalarNames.size()
1209  + vectorNames.size()
1210  + sphereNames.size()
1211  + symmNames.size()
1212  + tensorNames.size()
1213  );
1214 
1215  // Fields
1216  writer.writeIOField<label>(labelNames);
1217  writer.writeIOField<scalar>(scalarNames);
1218  writer.writeIOField<vector>(vectorNames);
1219  writer.writeIOField<sphericalTensor>(sphereNames);
1220  writer.writeIOField<symmTensor>(symmNames);
1221  writer.writeIOField<tensor>(tensorNames);
1222  }
1223  else
1224  {
1225  lagrangianWriter writer
1226  (
1227  vMesh,
1228  binary,
1229  lagrFileName,
1230  cloudName,
1231  true
1232  );
1233 
1234  // Write number of fields
1235  writer.writeFieldsHeader(0);
1236  }
1237  }
1238  }
1239 
1240 
1241  //---------------------------------------------------------------------
1242  //
1243  // Link parallel outputs back to undecomposed case for ease of loading
1244  //
1245  //---------------------------------------------------------------------
1246 
1247  if (Pstream::parRun() && doLinks)
1248  {
1249  mkDir(runTime.globalPath()/"VTK");
1250  chDir(runTime.globalPath()/"VTK");
1251 
1252  Info<< "Linking all processor files to " << runTime.globalPath()/"VTK"
1253  << endl;
1254 
1255  // Get list of vtk files
1256  fileName procVTK
1257  (
1258  fileName("..")
1259  /"processor" + name(Pstream::myProcNo())
1260  /"VTK"
1261  );
1262 
1263  fileNameList dirs(readDir(procVTK, fileType::directory));
1264  label sz = dirs.size();
1265  dirs.setSize(sz + 1);
1266  dirs[sz] = ".";
1267 
1268  forAll(dirs, i)
1269  {
1270  fileNameList subFiles(readDir(procVTK/dirs[i], fileType::file));
1271 
1272  forAll(subFiles, j)
1273  {
1274  fileName procFile(procVTK/dirs[i]/subFiles[j]);
1275 
1276  if (exists(procFile))
1277  {
1278  ln
1279  (
1280  procFile,
1281  "processor"
1282  + name(Pstream::myProcNo())
1283  + "_"
1284  + procFile.name()
1285  );
1286  }
1287  }
1288  }
1289  }
1290 
1291  Info<< "End\n" << endl;
1292 
1293  return 0;
1294 }
1295 
1296 
1297 // ************************************************************************* //
List< instant > instantList
List of instants.
Definition: instantList.H:42
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
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Foam::word regionName
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
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
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
Foam::meshFaceZones.
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:114
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
void writeCellDataHeader(std::ostream &, const label nCells, const label nFields)
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
Operations on lists of strings.
PrimitivePatch< IndirectList< face >, const pointField & > indirectPrimitivePatch
Foam::indirectPrimitivePatch.
fvMesh & mesh
bool optionReadIfPresent(const word &opt, T &) const
Read a value from the named option if present.
Definition: argListI.H:204
bool findStrings(const wordReListMatcher &matcher, const std::string &str)
Return true if string matches one of the regular expressions.
Definition: stringListOps.H:52
void writeSurfFields(const bool binary, const vtkMesh &vMesh, const fileName &fileName, const UPtrList< const surfaceVectorField > &surfVectorFields)
SymmTensor< scalar > symmTensor
SymmTensor of scalars.
Definition: symmTensor.H:48
bool isDir(const fileName &, const bool followLink=true)
Does the name exist as a directory in the file system?
Definition: POSIX.C:539
mkDir(pdfPath)
float floatScalar
Float precision floating point scalar type.
Definition: floatScalar.H:52
Write faceSet to vtk polydata file. Only one data which is original faceID.
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
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
graph_traits< Graph >::vertices_size_type size_type
Definition: SloanRenumber.C:73
bool ln(const fileName &src, const fileName &dst)
Create a softlink. dst should not exist. Returns true if successful.
Definition: POSIX.C:908
bool rmDir(const fileName &)
Remove a directory and its contents.
Definition: POSIX.C:1047
Write a patch with its data.
static const char nl
Definition: Ostream.H:260
objects
virtual void rename(const word &newName)
Rename.
Definition: regIOobject.C:417
Write pointSet to vtk polydata file. Only one data which is original pointID.
MeshZones< faceZone, polyMesh > meshFaceZones
A MeshZones with the type faceZone.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
List< word > wordList
A List of words.
Definition: fileName.H:54
label patchi
#define WarningInFunction
Report a warning using Foam::Warning.
void writePointDataHeader(std::ostream &, const label nPoints, const label nFields)
label n
const word cloudName(propsDict.lookup("cloudName"))
SphericalTensor< scalar > sphericalTensor
SphericalTensor of scalars.
Foam::argList args(argc, argv)
List< fileName > fileNameList
A List of fileNames.
Definition: fileNameList.H:50
Tensor< scalar > tensor
Tensor of scalars.
Definition: tensor.H:51
IStringStream optionLookup(const word &opt) const
Return an IStringStream from the named option.
Definition: argListI.H:120