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-2023 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 "meshFaceZones.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 
407  #include "createNamedMesh.H"
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  polyMesh::readUpdateState meshState = vMesh.readUpdate();
499 
500  const fvMesh& mesh = vMesh.mesh();
501 
502  if
503  (
504  meshState == polyMesh::TOPO_CHANGE
505  || meshState == polyMesh::TOPO_PATCH_CHANGE
506  )
507  {
508  Info<< " Read new mesh" << nl << endl;
509  }
510 
511  // If faceSet: write faceSet only (as polydata)
512  if (faceSetName.size())
513  {
514  // Load the faceSet
515  faceSet set(mesh, faceSetName);
516 
517  // Filename as if patch with same name.
518  mkDir(fvPath/set.name());
519 
520  fileName patchFileName
521  (
522  fvPath/set.name()/set.name()
523  + "_"
524  + timeDesc
525  + ".vtk"
526  );
527 
528  Info<< " FaceSet : " << patchFileName << endl;
529 
530  writeFaceSet(binary, vMesh, set, patchFileName);
531 
532  continue;
533  }
534 
535  // If pointSet: write pointSet only (as polydata)
536  if (pointSetName.size())
537  {
538  // Load the pointSet
539  pointSet set(mesh, pointSetName);
540 
541  // Filename as if patch with same name.
542  mkDir(fvPath/set.name());
543 
544  fileName patchFileName
545  (
546  fvPath/set.name()/set.name()
547  + "_"
548  + timeDesc
549  + ".vtk"
550  );
551 
552  Info<< " pointSet : " << patchFileName << endl;
553 
554  writePointSet(binary, vMesh, set, patchFileName);
555 
556  continue;
557  }
558 
559 
560  // Search for list of objects for this time
561  IOobjectList objects(mesh, runTime.name());
562 
563  HashSet<word> selectedFields;
564  bool specifiedFields = args.optionReadIfPresent
565  (
566  "fields",
567  selectedFields
568  );
569 
570 
571  // Construct the vol internal fields (on the original mesh if subsetted)
572 
578 
579  if (!specifiedFields || selectedFields.size())
580  {
581  readFields(vMesh, vMesh.baseMesh(), objects, selectedFields, visf);
582  print(" volScalarField::Internal :", Info, visf);
583 
584  readFields(vMesh, vMesh.baseMesh(), objects, selectedFields, vivf);
585  print(" volVectorField::Internal :", Info, vivf);
586 
587  readFields
588  (
589  vMesh,
590  vMesh.baseMesh(),
591  objects,
592  selectedFields,
593  visptf
594  );
595  print(" volSphericalTensorFields::Internal :", Info, visptf);
596 
597  readFields
598  (
599  vMesh,
600  vMesh.baseMesh(),
601  objects,
602  selectedFields,
603  visytf
604  );
605  print(" volSymmTensorFields::Internal :", Info, visytf);
606 
607  readFields(vMesh, vMesh.baseMesh(), objects, selectedFields, vitf);
608  print(" volTensorFields::Internal :", Info, vitf);
609  }
610 
611  label nVolInternalFields =
612  visf.size()
613  + vivf.size()
614  + visptf.size()
615  + visytf.size()
616  + vitf.size();
617 
618 
619  // Construct the vol fields (on the original mesh if subsetted)
620 
626 
627  if (!specifiedFields || selectedFields.size())
628  {
629  readFields(vMesh, vMesh.baseMesh(), objects, selectedFields, vsf);
630  print(" volScalarFields :", Info, vsf);
631 
632  readFields(vMesh, vMesh.baseMesh(), objects, selectedFields, vvf);
633  print(" volVectorFields :", Info, vvf);
634 
635  readFields
636  (
637  vMesh,
638  vMesh.baseMesh(),
639  objects,
640  selectedFields,
641  vsptf
642  );
643  print(" volSphericalTensorFields :", Info, vsptf);
644 
645  readFields
646  (
647  vMesh,
648  vMesh.baseMesh(),
649  objects,
650  selectedFields,
651  vsytf
652  );
653  print(" volSymmTensorFields :", Info, vsytf);
654 
655  readFields(vMesh, vMesh.baseMesh(), objects, selectedFields, vtf);
656  print(" volTensorFields :", Info, vtf);
657  }
658 
659  label nVolFields =
660  vsf.size()
661  + vvf.size()
662  + vsptf.size()
663  + vsytf.size()
664  + vtf.size();
665 
666 
667  // Construct pointMesh only if requested
668  if (noPointValues)
669  {
670  Info<< " pointScalarFields : switched off"
671  << " (\"-noPointValues\" (at your option)\n";
672  Info<< " pointVectorFields : switched off"
673  << " (\"-noPointValues\" (at your option)\n";
674  }
675 
681 
682  if (!noPointValues && !(specifiedFields && selectedFields.empty()))
683  {
684  readFields
685  (
686  vMesh,
687  pointMesh::New(vMesh.baseMesh()),
688  objects,
689  selectedFields,
690  psf
691  );
692  print(" pointScalarFields :", Info, psf);
693 
694  readFields
695  (
696  vMesh,
697  pointMesh::New(vMesh.baseMesh()),
698  objects,
699  selectedFields,
700  pvf
701  );
702  print(" pointVectorFields :", Info, pvf);
703 
704  readFields
705  (
706  vMesh,
707  pointMesh::New(vMesh.baseMesh()),
708  objects,
709  selectedFields,
710  psptf
711  );
712  print(" pointSphericalTensorFields :", Info, psptf);
713 
714  readFields
715  (
716  vMesh,
717  pointMesh::New(vMesh.baseMesh()),
718  objects,
719  selectedFields,
720  psytf
721  );
722  print(" pointSymmTensorFields :", Info, psytf);
723 
724  readFields
725  (
726  vMesh,
727  pointMesh::New(vMesh.baseMesh()),
728  objects,
729  selectedFields,
730  ptf
731  );
732  print(" pointTensorFields :", Info, ptf);
733  }
734  Info<< endl;
735 
736  label nPointFields =
737  psf.size()
738  + pvf.size()
739  + psptf.size()
740  + psytf.size()
741  + ptf.size();
742 
743  if (doWriteInternal)
744  {
745  // Create file and write header
746  fileName vtkFileName
747  (
748  fvPath/vtkName
749  + "_"
750  + timeDesc
751  + ".vtk"
752  );
753 
754  Info<< " Internal : " << vtkFileName << endl;
755 
756  // Write mesh
757  internalWriter writer(vMesh, binary, vtkFileName);
758 
759  // cellID + volFields::Internal + VolFields
761  (
762  writer.os(),
763  vMesh.nFieldCells(),
764  1 + nVolInternalFields + nVolFields
765  );
766 
767  // Write cellID field
768  writer.writeCellIDs();
769 
770  // Write volFields::Internal
771  writer.write(visf);
772  writer.write(vivf);
773  writer.write(visptf);
774  writer.write(visytf);
775  writer.write(vitf);
776 
777  // Write volFields
778  writer.write(vsf);
779  writer.write(vvf);
780  writer.write(vsptf);
781  writer.write(vsytf);
782  writer.write(vtf);
783 
784  if (!noPointValues)
785  {
787  (
788  writer.os(),
789  vMesh.nFieldPoints(),
790  nVolFields + nPointFields
791  );
792 
793  // pointFields
794  writer.write(psf);
795  writer.write(pvf);
796  writer.write(psptf);
797  writer.write(psytf);
798  writer.write(ptf);
799 
800  // Interpolated volFields
801  const volPointInterpolation& pInterp
802  (
804  );
805 
806  writer.write(pInterp, vsf);
807  writer.write(pInterp, vvf);
808  writer.write(pInterp, vsptf);
809  writer.write(pInterp, vsytf);
810  writer.write(pInterp, vtf);
811  }
812  }
813 
814  //---------------------------------------------------------------------
815  //
816  // Write surface fields
817  //
818  //---------------------------------------------------------------------
819 
820  if (args.optionFound("surfaceFields"))
821  {
823  readFields
824  (
825  vMesh,
826  vMesh.baseMesh(),
827  objects,
828  selectedFields,
829  ssf
830  );
831  print(" surfScalarFields :", Info, ssf);
832 
834  readFields
835  (
836  vMesh,
837  vMesh.baseMesh(),
838  objects,
839  selectedFields,
840  svf
841  );
842  print(" surfVectorFields :", Info, svf);
843 
844  if (ssf.size() + svf.size() > 0)
845  {
846  // Rework the scalar fields into vectorfields.
847  label sz = svf.size();
848 
849  svf.setSize(sz + ssf.size());
850 
851  surfaceVectorField n(mesh.Sf()/mesh.magSf());
852 
853  forAll(ssf, i)
854  {
855  surfaceVectorField* ssfiPtr = (ssf[i]*n).ptr();
856  ssfiPtr->rename(ssf[i].name());
857  svf.set(sz+i, ssfiPtr);
858  }
859  ssf.clear();
860 
861  mkDir(fvPath / "surfaceFields");
862 
863  fileName surfFileName
864  (
865  fvPath
866  /"surfaceFields"
867  /"surfaceFields"
868  + "_"
869  + timeDesc
870  + ".vtk"
871  );
872 
874  (
875  binary,
876  vMesh,
877  surfFileName,
878  svf
879  );
880  }
881  }
882 
883 
884  //---------------------------------------------------------------------
885  //
886  // Write patches (POLYDATA file, one for each patch)
887  //
888  //---------------------------------------------------------------------
889 
890  const polyBoundaryMesh& patches = mesh.boundaryMesh();
891 
892  const labelList patchIDs(getSelectedPatches(patches, excludePatches));
893 
894  if (allPatches)
895  {
896  mkDir(fvPath/"allPatches");
897 
898  fileName patchFileName;
899 
900  if (vMesh.useSubMesh())
901  {
902  patchFileName =
903  fvPath/"allPatches"/cellSetName
904  + "_"
905  + timeDesc
906  + ".vtk";
907  }
908  else
909  {
910  patchFileName =
911  fvPath/"allPatches"/"allPatches"
912  + "_"
913  + timeDesc
914  + ".vtk";
915  }
916 
917  Info<< " Combined patches : " << patchFileName << endl;
918 
919  patchWriter writer
920  (
921  vMesh,
922  binary,
923  nearCellValue,
924  patchFileName,
925  getSelectedPatches(patches, excludePatches)
926  );
927 
928  // VolFields + patchID
930  (
931  writer.os(),
932  writer.nFaces(),
933  1 + nVolFields
934  );
935 
936  // Write patchID field
937  writer.writePatchIDs();
938 
939  // Write volFields
940  writer.write(vsf);
941  writer.write(vvf);
942  writer.write(vsptf);
943  writer.write(vsytf);
944  writer.write(vtf);
945 
946  if (!noPointValues)
947  {
949  (
950  writer.os(),
951  writer.nPoints(),
952  nPointFields
953  );
954 
955  // Write pointFields
956  writer.write(psf);
957  writer.write(pvf);
958  writer.write(psptf);
959  writer.write(psytf);
960  writer.write(ptf);
961  }
962  }
963  else
964  {
965  forAll(patchIDs, i)
966  {
967  const polyPatch& pp = patches[patchIDs[i]];
968 
969  mkDir(fvPath/pp.name());
970 
971  fileName patchFileName;
972 
973  if (vMesh.useSubMesh())
974  {
975  patchFileName =
976  fvPath/pp.name()/cellSetName
977  + "_"
978  + timeDesc
979  + ".vtk";
980  }
981  else
982  {
983  patchFileName =
984  fvPath/pp.name()/pp.name()
985  + "_"
986  + timeDesc
987  + ".vtk";
988  }
989 
990  Info<< " Patch : " << patchFileName << endl;
991 
992  patchWriter writer
993  (
994  vMesh,
995  binary,
996  nearCellValue,
997  patchFileName,
998  labelList(1, patchIDs[i])
999  );
1000 
1001  // VolFields + patchID
1003  (
1004  writer.os(),
1005  writer.nFaces(),
1006  1 + nVolFields
1007  );
1008 
1009  // Write patchID field
1010  writer.writePatchIDs();
1011 
1012  // Write volFields
1013  writer.write(vsf);
1014  writer.write(vvf);
1015  writer.write(vsptf);
1016  writer.write(vsytf);
1017  writer.write(vtf);
1018 
1019  if (!noPointValues)
1020  {
1022  (
1023  writer.os(),
1024  writer.nPoints(),
1025  nVolFields
1026  + nPointFields
1027  );
1028 
1029  // Write pointFields
1030  writer.write(psf);
1031  writer.write(pvf);
1032  writer.write(psptf);
1033  writer.write(psytf);
1034  writer.write(ptf);
1035 
1037  (
1038  pp
1039  );
1040 
1041  // Write interpolated volFields
1042  writer.write(pInter, vsf);
1043  writer.write(pInter, vvf);
1044  writer.write(pInter, vsptf);
1045  writer.write(pInter, vsytf);
1046  writer.write(pInter, vtf);
1047  }
1048  }
1049  }
1050 
1051  //---------------------------------------------------------------------
1052  //
1053  // Write faceZones (POLYDATA file, one for each zone)
1054  //
1055  //---------------------------------------------------------------------
1056 
1057  if (doFaceZones)
1058  {
1060  readFields
1061  (
1062  vMesh,
1063  vMesh.baseMesh(),
1064  objects,
1065  selectedFields,
1066  ssf
1067  );
1068  print(" surfScalarFields :", Info, ssf);
1069 
1071  readFields
1072  (
1073  vMesh,
1074  vMesh.baseMesh(),
1075  objects,
1076  selectedFields,
1077  svf
1078  );
1079  print(" surfVectorFields :", Info, svf);
1080 
1081  const meshFaceZones& zones = mesh.faceZones();
1082 
1083  forAll(zones, zoneI)
1084  {
1085  const faceZone& fz = zones[zoneI];
1086 
1087  mkDir(fvPath/fz.name());
1088 
1089  fileName patchFileName;
1090 
1091  if (vMesh.useSubMesh())
1092  {
1093  patchFileName =
1094  fvPath/fz.name()/cellSetName
1095  + "_"
1096  + timeDesc
1097  + ".vtk";
1098  }
1099  else
1100  {
1101  patchFileName =
1102  fvPath/fz.name()/fz.name()
1103  + "_"
1104  + timeDesc
1105  + ".vtk";
1106  }
1107 
1108  Info<< " FaceZone : " << patchFileName << endl;
1109 
1111  (
1112  IndirectList<face>(mesh.faces(), fz),
1113  mesh.points()
1114  );
1115 
1116  surfaceMeshWriter writer
1117  (
1118  binary,
1119  pp,
1120  fz.name(),
1121  patchFileName
1122  );
1123 
1124  // Number of fields
1126  (
1127  writer.os(),
1128  pp.size(),
1129  ssf.size() + svf.size()
1130  );
1131 
1132  writer.write(ssf);
1133  writer.write(svf);
1134  }
1135  }
1136 
1137 
1138 
1139  //---------------------------------------------------------------------
1140  //
1141  // Write lagrangian data
1142  //
1143  //---------------------------------------------------------------------
1144 
1145  forAllConstIter(HashSet<fileName>, allCloudDirs, iter)
1146  {
1147  const fileName& cloudName = iter.key();
1148 
1149  // Always create the cloud directory.
1150  mkDir(fvPath/cloud::prefix/cloudName);
1151 
1152  fileName lagrFileName
1153  (
1155  + "_" + timeDesc + ".vtk"
1156  );
1157 
1158  Info<< " Lagrangian: " << lagrFileName << endl;
1159 
1160 
1161  IOobjectList sprayObjs
1162  (
1163  mesh,
1164  runTime.name(),
1166  );
1167 
1168  IOobject* positionsPtr = sprayObjs.lookup(word("positions"));
1169 
1170  if (positionsPtr)
1171  {
1172  wordList labelNames(sprayObjs.names(labelIOField::typeName));
1173  Info<< " labels :";
1174  print(Info, labelNames);
1175 
1176  wordList scalarNames(sprayObjs.names(scalarIOField::typeName));
1177  Info<< " scalars :";
1178  print(Info, scalarNames);
1179 
1180  wordList vectorNames(sprayObjs.names(vectorIOField::typeName));
1181  Info<< " vectors :";
1182  print(Info, vectorNames);
1183 
1184  wordList sphereNames
1185  (
1186  sprayObjs.names
1187  (
1189  )
1190  );
1191  Info<< " spherical tensors :";
1192  print(Info, sphereNames);
1193 
1194  wordList symmNames
1195  (
1196  sprayObjs.names
1197  (
1199  )
1200  );
1201  Info<< " symm tensors :";
1202  print(Info, symmNames);
1203 
1204  wordList tensorNames(sprayObjs.names(tensorIOField::typeName));
1205  Info<< " tensors :";
1206  print(Info, tensorNames);
1207 
1208  lagrangianWriter writer
1209  (
1210  vMesh,
1211  binary,
1212  lagrFileName,
1213  cloudName,
1214  false
1215  );
1216 
1217  // Write number of fields
1218  writer.writeFieldsHeader
1219  (
1220  labelNames.size()
1221  + scalarNames.size()
1222  + vectorNames.size()
1223  + sphereNames.size()
1224  + symmNames.size()
1225  + tensorNames.size()
1226  );
1227 
1228  // Fields
1229  writer.writeIOField<label>(labelNames);
1230  writer.writeIOField<scalar>(scalarNames);
1231  writer.writeIOField<vector>(vectorNames);
1232  writer.writeIOField<sphericalTensor>(sphereNames);
1233  writer.writeIOField<symmTensor>(symmNames);
1234  writer.writeIOField<tensor>(tensorNames);
1235  }
1236  else
1237  {
1238  lagrangianWriter writer
1239  (
1240  vMesh,
1241  binary,
1242  lagrFileName,
1243  cloudName,
1244  true
1245  );
1246 
1247  // Write number of fields
1248  writer.writeFieldsHeader(0);
1249  }
1250  }
1251  }
1252 
1253 
1254  //---------------------------------------------------------------------
1255  //
1256  // Link parallel outputs back to undecomposed case for ease of loading
1257  //
1258  //---------------------------------------------------------------------
1259 
1260  if (Pstream::parRun() && doLinks)
1261  {
1262  mkDir(runTime.globalPath()/"VTK");
1263  chDir(runTime.globalPath()/"VTK");
1264 
1265  Info<< "Linking all processor files to " << runTime.globalPath()/"VTK"
1266  << endl;
1267 
1268  // Get list of vtk files
1269  fileName procVTK
1270  (
1271  fileName("..")
1272  /"processor" + name(Pstream::myProcNo())
1273  /"VTK"
1274  );
1275 
1276  fileNameList dirs(readDir(procVTK, fileType::directory));
1277  label sz = dirs.size();
1278  dirs.setSize(sz + 1);
1279  dirs[sz] = ".";
1280 
1281  forAll(dirs, i)
1282  {
1283  fileNameList subFiles(readDir(procVTK/dirs[i], fileType::file));
1284 
1285  forAll(subFiles, j)
1286  {
1287  fileName procFile(procVTK/dirs[i]/subFiles[j]);
1288 
1289  if (exists(procFile))
1290  {
1291  ln
1292  (
1293  procFile,
1294  "processor"
1295  + name(Pstream::myProcNo())
1296  + "_"
1297  + procFile.name()
1298  );
1299  }
1300  }
1301  }
1302  }
1303 
1304  Info<< "End\n" << endl;
1305 
1306  return 0;
1307 }
1308 
1309 
1310 // ************************************************************************* //
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 const char *const typeName
Definition: Field.H:105
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:111
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:65
void clear()
Clear the PtrList, i.e. set size to zero deleting all the.
Definition: PtrList.C:174
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,...
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
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:68
A class for handling file names.
Definition: fileName.H:82
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
const surfaceVectorField & Sf() const
Return cell face area vectors.
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
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:268
const meshFaceZones & faceZones() const
Return face zones.
Definition: polyMesh.H:447
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1373
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: polyMesh.H:91
@ TOPO_PATCH_CHANGE
Definition: polyMesh.H:95
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:405
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1360
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:53
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
const word & name() const
Return name.
Definition: zone.H:147
Foam::word regionName
int main(int argc, char *argv[])
Definition: financialFoam.C:44
label patchi
static instantList timeDirs
Definition: globalFoam.H:44
const fvPatchList & patches
Foam::meshFaceZones.
#define WarningInFunction
Report a warning using Foam::Warning.
void writePointDataHeader(std::ostream &, const label nPoints, const label nFields)
void writeCellDataHeader(std::ostream &, const label nCells, const label nFields)
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:251
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
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:260
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.