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