foamToVTK.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2015 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 volScalar, volVector, pointScalar, pointVector, surfaceScalar
31  fields.
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 
41  - foamToVTK [OPTION]
42 
43  \param -ascii \n
44  Write VTK data in ASCII format instead of binary.
45 
46  \param -mesh <name>\n
47  Use a different mesh name (instead of -region)
48 
49  \param -fields <fields>\n
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  \param -surfaceFields \n
58  Write surfaceScalarFields (e.g., phi)
59 
60  \param -cellSet <name>\n
61  \param -faceSet <name>\n
62  \param -pointSet <name>\n
63  Restrict conversion to the cellSet, faceSet or pointSet.
64 
65  \param -nearCellValue \n
66  Output cell value on patches instead of patch value itself
67 
68  \param -noInternal \n
69  Do not generate file for mesh, only for patches
70 
71  \param -noPointValues \n
72  No pointFields
73 
74  \param -noFaceZones \n
75  No faceZones
76 
77  \param -noLinks \n
78  (in parallel) do not link processor files to master
79 
80  \param poly \n
81  write polyhedral cells without tet/pyramid decomposition
82 
83  \param -allPatches \n
84  Combine all patches into a single file
85 
86  \param -excludePatches <patchNames>\n
87  Specify patches (wildcards) to exclude. For example,
88  \verbatim
89  -excludePatches '( inlet_1 inlet_2 "proc.*")'
90  \endverbatim
91  The quoting is required to avoid shell expansions and to pass the
92  information as a single argument. The double quotes denote a regular
93  expression.
94 
95  \param -useTimeName \n
96  use the time index in the VTK file name instead of the time index
97 
98 Note
99  mesh subset is handled by vtkMesh. Slight inconsistency in
100  interpolation: on the internal field it interpolates the whole volField
101  to the whole-mesh pointField and then selects only those values it
102  needs for the subMesh (using the fvMeshSubset cellMap(), pointMap()
103  functions). For the patches however it uses the
104  fvMeshSubset.interpolate function to directly interpolate the
105  whole-mesh values onto the subset patch.
106 
107 Note
108  new file format: \n
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 "labelIOField.H"
141 #include "scalarIOField.H"
142 #include "sphericalTensorIOField.H"
143 #include "symmTensorIOField.H"
144 #include "tensorIOField.H"
145 #include "faceZoneMesh.H"
146 #include "Cloud.H"
147 #include "passiveParticle.H"
148 #include "stringListOps.H"
149 
150 #include "vtkMesh.H"
151 #include "readFields.H"
152 #include "writeFuns.H"
153 
154 #include "internalWriter.H"
155 #include "patchWriter.H"
156 #include "lagrangianWriter.H"
157 
158 #include "writeFaceSet.H"
159 #include "writePointSet.H"
160 #include "surfaceMeshWriter.H"
161 #include "writeSurfFields.H"
162 
163 
164 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
165 
166 template<class GeoField>
167 void print(const char* msg, Ostream& os, const PtrList<GeoField>& flds)
168 {
169  if (flds.size())
170  {
171  os << msg;
172  forAll(flds, i)
173  {
174  os << ' ' << flds[i].name();
175  }
176  os << endl;
177  }
178 }
179 
180 
181 void print(Ostream& os, const wordList& flds)
182 {
183  forAll(flds, i)
184  {
185  os << ' ' << flds[i];
186  }
187  os << endl;
188 }
189 
190 
192 (
193  const polyBoundaryMesh& patches,
194  const List<wordRe>& excludePatches //HashSet<word>& excludePatches
195 )
196 {
197  DynamicList<label> patchIDs(patches.size());
198 
199  Info<< "Combining patches:" << endl;
200 
201  forAll(patches, patchI)
202  {
203  const polyPatch& pp = patches[patchI];
204 
205  if
206  (
207  isType<emptyPolyPatch>(pp)
208  || (Pstream::parRun() && isType<processorPolyPatch>(pp))
209  )
210  {
211  Info<< " discarding empty/processor patch " << patchI
212  << " " << pp.name() << endl;
213  }
214  else if (findStrings(excludePatches, pp.name()))
215  {
216  Info<< " excluding patch " << patchI
217  << " " << pp.name() << endl;
218  }
219  else
220  {
221  patchIDs.append(patchI);
222  Info<< " patch " << patchI << " " << pp.name() << endl;
223  }
224  }
225  return patchIDs.shrink();
226 }
227 
228 
229 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
230 
231 int main(int argc, char *argv[])
232 {
233  cout.setf(ios_base::fixed, ios_base::floatfield);//十进制计数法,不是科学计数法
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  const 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(doubleScalar) != 4 || sizeof(label) != 4))
335  {
337  << "doubleScalar and/or label are not 4 bytes in size" << nl
338  << "Hence cannot use binary VTK format. Please use -ascii"
339  << exit(FatalError);
340  }
341 
342  const bool nearCellValue = args.optionFound("nearCellValue");
343 
344  if (nearCellValue)
345  {
347  << "Using neighbouring cell value instead of patch value"
348  << nl << endl;
349  }
350 
351  const bool noPointValues = args.optionFound("noPointValues");
352 
353  if (noPointValues)
354  {
356  << "Outputting cell values only" << nl << endl;
357  }
358 
359  const bool allPatches = args.optionFound("allPatches");
360 
361  List<wordRe> excludePatches;
362  if (args.optionFound("excludePatches"))
363  {
364  args.optionLookup("excludePatches")() >> excludePatches;
365 
366  Info<< "Not including patches " << excludePatches << nl << endl;
367  }
368 
369  word cellSetName;
370  word faceSetName;
371  word pointSetName;
372  string vtkName = runTime.caseName();
373 
374  if (args.optionReadIfPresent("cellSet", cellSetName))
375  {
376  vtkName = cellSetName;
377  }
378  else if (Pstream::parRun())
379  {
380  // Strip off leading casename, leaving just processor_DDD ending.
381  string::size_type i = vtkName.rfind("processor");
382 
383  if (i != string::npos)
384  {
385  vtkName = vtkName.substr(i);
386  }
387  }
388  args.optionReadIfPresent("faceSet", faceSetName);
389  args.optionReadIfPresent("pointSet", pointSetName);
390 
391 
392 
393  instantList timeDirs = timeSelector::select0(runTime, args);
394 
395  #include "createNamedMesh.H"
396 
397  // VTK/ directory in the case
398  fileName fvPath(runTime.path()/"VTK");
399  // Directory of mesh (region0 gets filtered out)
400  fileName regionPrefix = "";
401 
402  if (regionName != polyMesh::defaultRegion)
403  {
404  fvPath = fvPath/regionName;
405  regionPrefix = regionName;
406  }
407 
408  if (isDir(fvPath))
409  {
410  if
411  (
412  args.optionFound("time")
413  || args.optionFound("latestTime")
414  || cellSetName.size()
415  || faceSetName.size()
416  || pointSetName.size()
417  || regionName != polyMesh::defaultRegion
418  )
419  {
420  Info<< "Keeping old VTK files in " << fvPath << nl << endl;
421  }
422  else
423  {
424  Info<< "Deleting old VTK files in " << fvPath << nl << endl;
425 
426  rmDir(fvPath);
427  }
428  }
429 
430  mkDir(fvPath);
431 
432 
433  // Mesh wrapper; does subsetting and decomposition
434  vtkMesh vMesh(mesh, cellSetName);
435 
436 
437  // Scan for all possible lagrangian clouds
438  HashSet<fileName> allCloudDirs;
439  forAll(timeDirs, timeI)
440  {
441  runTime.setTime(timeDirs[timeI], timeI);
442  fileNameList cloudDirs
443  (
444  readDir
445  (
446  runTime.timePath()/regionPrefix/cloud::prefix,
447  fileName::DIRECTORY
448  )
449  );
450  forAll(cloudDirs, i)
451  {
452  IOobjectList sprayObjs
453  (
454  mesh,
455  runTime.timeName(),
456  cloud::prefix/cloudDirs[i]
457  );
458 
459  IOobject* positionsPtr = sprayObjs.lookup(word("positions"));
460 
461  if (positionsPtr)
462  {
463  if (allCloudDirs.insert(cloudDirs[i]))
464  {
465  Info<< "At time: " << runTime.timeName()
466  << " detected cloud directory : " << cloudDirs[i]
467  << endl;
468  }
469  }
470  }
471  }
472 
473 
474  forAll(timeDirs, timeI)
475  {
476  runTime.setTime(timeDirs[timeI], timeI);
477 
478  Info<< "Time: " << runTime.timeName() << endl;
479 
480  word timeDesc =
481  useTimeName ? runTime.timeName() : Foam::name(runTime.timeIndex());
482 
483  // Check for new polyMesh/ and update mesh, fvMeshSubset and cell
484  // decomposition.
485  polyMesh::readUpdateState meshState = vMesh.readUpdate();
486 
487  const fvMesh& mesh = vMesh.mesh();
488 
489  if
490  (
491  meshState == polyMesh::TOPO_CHANGE
492  || meshState == polyMesh::TOPO_PATCH_CHANGE
493  )
494  {
495  Info<< " Read new mesh" << nl << endl;
496  }
497 
498 
499  // If faceSet: write faceSet only (as polydata)
500  if (faceSetName.size())
501  {
502  // Load the faceSet
503  faceSet set(mesh, faceSetName);
504 
505  // Filename as if patch with same name.
506  mkDir(fvPath/set.name());
507 
508  fileName patchFileName
509  (
510  fvPath/set.name()/set.name()
511  + "_"
512  + timeDesc
513  + ".vtk"
514  );
515 
516  Info<< " FaceSet : " << patchFileName << endl;
517 
518  writeFaceSet(binary, vMesh, set, patchFileName);
519 
520  continue;
521  }
522  // If pointSet: write pointSet only (as polydata)
523  if (pointSetName.size())
524  {
525  // Load the pointSet
526  pointSet set(mesh, pointSetName);
527 
528  // Filename as if patch with same name.
529  mkDir(fvPath/set.name());
530 
531  fileName patchFileName
532  (
533  fvPath/set.name()/set.name()
534  + "_"
535  + timeDesc
536  + ".vtk"
537  );
538 
539  Info<< " pointSet : " << patchFileName << endl;
540 
541  writePointSet(binary, vMesh, set, patchFileName);
542 
543  continue;
544  }
545 
546 
547  // Search for list of objects for this time
548  IOobjectList objects(mesh, runTime.timeName());
549 
550  HashSet<word> selectedFields;
551  bool specifiedFields = args.optionReadIfPresent
552  (
553  "fields",
554  selectedFields
555  );
556 
557  // Construct the vol fields (on the original mesh if subsetted)
558 
559  PtrList<volScalarField> vsf;
560  PtrList<volVectorField> vvf;
561  PtrList<volSphericalTensorField> vSpheretf;
562  PtrList<volSymmTensorField> vSymmtf;
563  PtrList<volTensorField> vtf;
564 
565  if (!specifiedFields || selectedFields.size())
566  {
567  readFields(vMesh, vMesh.baseMesh(), objects, selectedFields, vsf);
568  print(" volScalarFields :", Info, vsf);
569 
570  readFields(vMesh, vMesh.baseMesh(), objects, selectedFields, vvf);
571  print(" volVectorFields :", Info, vvf);
572 
573  readFields
574  (
575  vMesh,
576  vMesh.baseMesh(),
577  objects,
578  selectedFields,
579  vSpheretf
580  );
581  print(" volSphericalTensorFields :", Info, vSpheretf);
582 
583  readFields
584  (
585  vMesh,
586  vMesh.baseMesh(),
587  objects,
588  selectedFields,
589  vSymmtf
590  );
591  print(" volSymmTensorFields :", Info, vSymmtf);
592 
593  readFields(vMesh, vMesh.baseMesh(), objects, selectedFields, vtf);
594  print(" volTensorFields :", Info, vtf);
595  }
596 
597  label nVolFields =
598  vsf.size()
599  + vvf.size()
600  + vSpheretf.size()
601  + vSymmtf.size()
602  + vtf.size();
603 
604 
605  // Construct dimensioned fields
606  PtrList<volScalarField::DimensionedInternalField> dsf;
607  PtrList<volVectorField::DimensionedInternalField> dvf;
608  PtrList<volSphericalTensorField::DimensionedInternalField> dSpheretf;
609  PtrList<volSymmTensorField::DimensionedInternalField> dSymmtf;
610  PtrList<volTensorField::DimensionedInternalField> dtf;
611 
612  if (!specifiedFields || selectedFields.size())
613  {
614  readFields(vMesh, vMesh.baseMesh(), objects, selectedFields, dsf);
615  print(" volScalarFields::Internal :", Info, dsf);
616 
617  readFields(vMesh, vMesh.baseMesh(), objects, selectedFields, dvf);
618  print(" volVectorFields::Internal :", Info, dvf);
619 
620  readFields
621  (
622  vMesh,
623  vMesh.baseMesh(),
624  objects,
625  selectedFields,
626  dSpheretf
627  );
628  print(" volSphericalTensorFields::Internal :", Info, dSpheretf);
629 
630  readFields
631  (
632  vMesh,
633  vMesh.baseMesh(),
634  objects,
635  selectedFields,
636  dSymmtf
637  );
638  print(" volSymmTensorFields::Internal :", Info, dSymmtf);
639 
640  readFields(vMesh, vMesh.baseMesh(), objects, selectedFields, dtf);
641  print(" volTensorFields::Internal :", Info, dtf);
642  }
643 
644  label nDimFields =
645  dsf.size()
646  + dvf.size()
647  + dSpheretf.size()
648  + dSymmtf.size()
649  + dtf.size();
650 
651 
652  // Construct pointMesh only if nessecary since constructs edge
653  // addressing (expensive on polyhedral meshes)
654  if (noPointValues)
655  {
656  Info<< " pointScalarFields : switched off"
657  << " (\"-noPointValues\" (at your option)\n";
658  Info<< " pointVectorFields : switched off"
659  << " (\"-noPointValues\" (at your option)\n";
660  }
661 
662  PtrList<pointScalarField> psf;
663  PtrList<pointVectorField> pvf;
664  PtrList<pointSphericalTensorField> pSpheretf;
665  PtrList<pointSymmTensorField> pSymmtf;
666  PtrList<pointTensorField> ptf;
667 
668  if (!noPointValues && !(specifiedFields && selectedFields.empty()))
669  {
670  readFields
671  (
672  vMesh,
673  pointMesh::New(vMesh.baseMesh()),
674  objects,
675  selectedFields,
676  psf
677  );
678  print(" pointScalarFields :", Info, psf);
679 
680  readFields
681  (
682  vMesh,
683  pointMesh::New(vMesh.baseMesh()),
684  objects,
685  selectedFields,
686  pvf
687  );
688  print(" pointVectorFields :", Info, pvf);
689 
690  readFields
691  (
692  vMesh,
693  pointMesh::New(vMesh.baseMesh()),
694  objects,
695  selectedFields,
696  pSpheretf
697  );
698  print(" pointSphericalTensorFields :", Info, pSpheretf);
699 
700  readFields
701  (
702  vMesh,
703  pointMesh::New(vMesh.baseMesh()),
704  objects,
705  selectedFields,
706  pSymmtf
707  );
708  print(" pointSymmTensorFields :", Info, pSymmtf);
709 
710  readFields
711  (
712  vMesh,
713  pointMesh::New(vMesh.baseMesh()),
714  objects,
715  selectedFields,
716  ptf
717  );
718  print(" pointTensorFields :", Info, ptf);
719  }
720  Info<< endl;
721 
722  label nPointFields =
723  psf.size()
724  + pvf.size()
725  + pSpheretf.size()
726  + pSymmtf.size()
727  + ptf.size();
728 
729  if (doWriteInternal)
730  {
731  //
732  // Create file and write header
733  //
734  fileName vtkFileName
735  (
736  fvPath/vtkName
737  + "_"
738  + timeDesc
739  + ".vtk"
740  );
741 
742  Info<< " Internal : " << vtkFileName << endl;
743 
744  // Write mesh
745  internalWriter writer(vMesh, binary, vtkFileName);
746 
747  // VolFields + cellID
748  writeFuns::writeCellDataHeader
749  (
750  writer.os(),
751  vMesh.nFieldCells(),
752  1 + nVolFields + nDimFields
753  );
754 
755  // Write cellID field
756  writer.writeCellIDs();
757 
758  // Write volFields
759  writer.write(vsf);
760  writer.write(vvf);
761  writer.write(vSpheretf);
762  writer.write(vSymmtf);
763  writer.write(vtf);
764 
765  // Write dimensionedFields
766  writer.write<scalar, volMesh>(dsf);
767  writer.write<vector, volMesh>(dvf);
768  writer.write<sphericalTensor, volMesh>(dSpheretf);
769  writer.write<symmTensor, volMesh>(dSymmtf);
770  writer.write<tensor, volMesh>(dtf);
771 
772  if (!noPointValues)
773  {
774  writeFuns::writePointDataHeader
775  (
776  writer.os(),
777  vMesh.nFieldPoints(),
778  nVolFields + nDimFields + nPointFields
779  );
780 
781  // pointFields
782  writer.write(psf);
783  writer.write(pvf);
784  writer.write(pSpheretf);
785  writer.write(pSymmtf);
786  writer.write(ptf);
787 
788  // Interpolated volFields
789  volPointInterpolation pInterp(mesh);
790  writer.write(pInterp, vsf);
791  writer.write(pInterp, vvf);
792  writer.write(pInterp, vSpheretf);
793  writer.write(pInterp, vSymmtf);
794  writer.write(pInterp, vtf);
795 
796  writer.write<scalar, volMesh>(pInterp, dsf);
797  writer.write<vector, volMesh>(pInterp, dvf);
798  writer.write<sphericalTensor, volMesh>(pInterp, dSpheretf);
799  writer.write<symmTensor, volMesh>(pInterp, dSymmtf);
800  writer.write<tensor, volMesh>(pInterp, dtf);
801  }
802  }
803 
804  //---------------------------------------------------------------------
805  //
806  // Write surface fields
807  //
808  //---------------------------------------------------------------------
809 
810  if (args.optionFound("surfaceFields"))
811  {
812  PtrList<surfaceScalarField> ssf;
813  readFields
814  (
815  vMesh,
816  vMesh.baseMesh(),
817  objects,
818  selectedFields,
819  ssf
820  );
821  print(" surfScalarFields :", Info, ssf);
822 
823  PtrList<surfaceVectorField> svf;
824  readFields
825  (
826  vMesh,
827  vMesh.baseMesh(),
828  objects,
829  selectedFields,
830  svf
831  );
832  print(" surfVectorFields :", Info, svf);
833 
834  if (ssf.size() + svf.size() > 0)
835  {
836  // Rework the scalar fields into vectorfields.
837  label sz = svf.size();
838 
839  svf.setSize(sz+ssf.size());
840 
841  surfaceVectorField n(mesh.Sf()/mesh.magSf());
842 
843  forAll(ssf, i)
844  {
845  svf.set(sz+i, ssf[i]*n);
846  svf[sz+i].rename(ssf[i].name());
847  }
848  ssf.clear();
849 
850  mkDir(fvPath / "surfaceFields");
851 
852  fileName surfFileName
853  (
854  fvPath
855  /"surfaceFields"
856  /"surfaceFields"
857  + "_"
858  + timeDesc
859  + ".vtk"
860  );
861 
863  (
864  binary,
865  vMesh,
866  surfFileName,
867  svf
868  );
869  }
870  }
871 
872 
873  //---------------------------------------------------------------------
874  //
875  // Write patches (POLYDATA file, one for each patch)
876  //
877  //---------------------------------------------------------------------
878 
879  const polyBoundaryMesh& patches = mesh.boundaryMesh();
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
916  writeFuns::writeCellDataHeader
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(vSpheretf);
930  writer.write(vSymmtf);
931  writer.write(vtf);
932 
933  if (!noPointValues)
934  {
935  writeFuns::writePointDataHeader
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(pSpheretf);
946  writer.write(pSymmtf);
947  writer.write(ptf);
948 
949  // no interpolated volFields since I cannot be bothered to
950  // create the patchInterpolation for all subpatches.
951  }
952  }
953  else
954  {
955  forAll(patches, patchI)
956  {
957  const polyPatch& pp = patches[patchI];
958 
959  if (!findStrings(excludePatches, pp.name()))
960  {
961  mkDir(fvPath/pp.name());
962 
963  fileName patchFileName;
964 
965  if (vMesh.useSubMesh())
966  {
967  patchFileName =
968  fvPath/pp.name()/cellSetName
969  + "_"
970  + timeDesc
971  + ".vtk";
972  }
973  else
974  {
975  patchFileName =
976  fvPath/pp.name()/pp.name()
977  + "_"
978  + timeDesc
979  + ".vtk";
980  }
981 
982  Info<< " Patch : " << patchFileName << endl;
983 
984  patchWriter writer
985  (
986  vMesh,
987  binary,
988  nearCellValue,
989  patchFileName,
990  labelList(1, patchI)
991  );
992 
993  if (!isA<emptyPolyPatch>(pp))
994  {
995  // VolFields + patchID
996  writeFuns::writeCellDataHeader
997  (
998  writer.os(),
999  writer.nFaces(),
1000  1+nVolFields
1001  );
1002 
1003  // Write patchID field
1004  writer.writePatchIDs();
1005 
1006  // Write volFields
1007  writer.write(vsf);
1008  writer.write(vvf);
1009  writer.write(vSpheretf);
1010  writer.write(vSymmtf);
1011  writer.write(vtf);
1012 
1013  if (!noPointValues)
1014  {
1015  writeFuns::writePointDataHeader
1016  (
1017  writer.os(),
1018  writer.nPoints(),
1019  nVolFields
1020  + nPointFields
1021  );
1022 
1023  // Write pointFields
1024  writer.write(psf);
1025  writer.write(pvf);
1026  writer.write(pSpheretf);
1027  writer.write(pSymmtf);
1028  writer.write(ptf);
1029 
1030  PrimitivePatchInterpolation<primitivePatch> pInter
1031  (
1032  pp
1033  );
1034 
1035  // Write interpolated volFields
1036  writer.write(pInter, vsf);
1037  writer.write(pInter, vvf);
1038  writer.write(pInter, vSpheretf);
1039  writer.write(pInter, vSymmtf);
1040  writer.write(pInter, vtf);
1041  }
1042  }
1043  }
1044  }
1045  }
1046 
1047  //---------------------------------------------------------------------
1048  //
1049  // Write faceZones (POLYDATA file, one for each zone)
1050  //
1051  //---------------------------------------------------------------------
1052 
1053  if (doFaceZones)
1054  {
1055  PtrList<surfaceScalarField> ssf;
1056  readFields
1057  (
1058  vMesh,
1059  vMesh.baseMesh(),
1060  objects,
1061  selectedFields,
1062  ssf
1063  );
1064  print(" surfScalarFields :", Info, ssf);
1065 
1066  PtrList<surfaceVectorField> svf;
1067  readFields
1068  (
1069  vMesh,
1070  vMesh.baseMesh(),
1071  objects,
1072  selectedFields,
1073  svf
1074  );
1075  print(" surfVectorFields :", Info, svf);
1076 
1077  const faceZoneMesh& zones = mesh.faceZones();
1078 
1079  forAll(zones, zoneI)
1080  {
1081  const faceZone& fz = zones[zoneI];
1082 
1083  mkDir(fvPath/fz.name());
1084 
1085  fileName patchFileName;
1086 
1087  if (vMesh.useSubMesh())
1088  {
1089  patchFileName =
1090  fvPath/fz.name()/cellSetName
1091  + "_"
1092  + timeDesc
1093  + ".vtk";
1094  }
1095  else
1096  {
1097  patchFileName =
1098  fvPath/fz.name()/fz.name()
1099  + "_"
1100  + timeDesc
1101  + ".vtk";
1102  }
1103 
1104  Info<< " FaceZone : " << patchFileName << endl;
1105 
1107  (
1108  IndirectList<face>(mesh.faces(), fz),
1109  mesh.points()
1110  );
1111 
1112  surfaceMeshWriter writer
1113  (
1114  binary,
1115  pp,
1116  fz.name(),
1117  patchFileName
1118  );
1119 
1120  // Number of fields
1121  writeFuns::writeCellDataHeader
1122  (
1123  writer.os(),
1124  pp.size(),
1125  ssf.size()+svf.size()
1126  );
1127 
1128  writer.write(ssf);
1129  writer.write(svf);
1130  }
1131  }
1132 
1133 
1134 
1135  //---------------------------------------------------------------------
1136  //
1137  // Write lagrangian data
1138  //
1139  //---------------------------------------------------------------------
1140 
1141  forAllConstIter(HashSet<fileName>, allCloudDirs, iter)
1142  {
1143  const fileName& cloudName = iter.key();
1144 
1145  // Always create the cloud directory.
1146  mkDir(fvPath/cloud::prefix/cloudName);
1147 
1148  fileName lagrFileName
1149  (
1150  fvPath/cloud::prefix/cloudName/cloudName
1151  + "_" + timeDesc + ".vtk"
1152  );
1153 
1154  Info<< " Lagrangian: " << lagrFileName << endl;
1155 
1156 
1157  IOobjectList sprayObjs
1158  (
1159  mesh,
1160  runTime.timeName(),
1161  cloud::prefix/cloudName
1162  );
1163 
1164  IOobject* positionsPtr = sprayObjs.lookup(word("positions"));
1165 
1166  if (positionsPtr)
1167  {
1168  wordList labelNames(sprayObjs.names(labelIOField::typeName));
1169  Info<< " labels :";
1170  print(Info, labelNames);
1171 
1172  wordList scalarNames(sprayObjs.names(scalarIOField::typeName));
1173  Info<< " scalars :";
1174  print(Info, scalarNames);
1175 
1176  wordList vectorNames(sprayObjs.names(vectorIOField::typeName));
1177  Info<< " vectors :";
1178  print(Info, vectorNames);
1179 
1180  wordList sphereNames
1181  (
1182  sprayObjs.names
1183  (
1184  sphericalTensorIOField::typeName
1185  )
1186  );
1187  Info<< " spherical tensors :";
1188  print(Info, sphereNames);
1189 
1190  wordList symmNames
1191  (
1192  sprayObjs.names
1193  (
1194  symmTensorIOField::typeName
1195  )
1196  );
1197  Info<< " symm tensors :";
1198  print(Info, symmNames);
1199 
1200  wordList tensorNames(sprayObjs.names(tensorIOField::typeName));
1201  Info<< " tensors :";
1202  print(Info, tensorNames);
1203 
1204  lagrangianWriter writer
1205  (
1206  vMesh,
1207  binary,
1208  lagrFileName,
1209  cloudName,
1210  false
1211  );
1212 
1213  // Write number of fields
1214  writer.writeParcelHeader
1215  (
1216  labelNames.size()
1217  + scalarNames.size()
1218  + vectorNames.size()
1219  + sphereNames.size()
1220  + symmNames.size()
1221  + tensorNames.size()
1222  );
1223 
1224  // Fields
1225  writer.writeIOField<label>(labelNames);
1226  writer.writeIOField<scalar>(scalarNames);
1227  writer.writeIOField<vector>(vectorNames);
1228  writer.writeIOField<sphericalTensor>(sphereNames);
1229  writer.writeIOField<symmTensor>(symmNames);
1230  writer.writeIOField<tensor>(tensorNames);
1231  }
1232  else
1233  {
1234  lagrangianWriter writer
1235  (
1236  vMesh,
1237  binary,
1238  lagrFileName,
1239  cloudName,
1240  true
1241  );
1242 
1243  // Write number of fields
1244  writer.writeParcelHeader(0);
1245  }
1246  }
1247  }
1248 
1249 
1250  //---------------------------------------------------------------------
1251  //
1252  // Link parallel outputs back to undecomposed case for ease of loading
1253  //
1254  //---------------------------------------------------------------------
1255 
1256  if (Pstream::parRun() && doLinks)
1257  {
1258  mkDir(runTime.path()/".."/"VTK");
1259  chDir(runTime.path()/".."/"VTK");
1260 
1261  Info<< "Linking all processor files to " << runTime.path()/".."/"VTK"
1262  << endl;
1263 
1264  // Get list of vtk files
1265  fileName procVTK
1266  (
1267  fileName("..")
1268  /"processor" + name(Pstream::myProcNo())
1269  /"VTK"
1270  );
1271 
1272  fileNameList dirs(readDir(procVTK, fileName::DIRECTORY));
1273  label sz = dirs.size();
1274  dirs.setSize(sz+1);
1275  dirs[sz] = ".";
1276 
1277  forAll(dirs, i)
1278  {
1279  fileNameList subFiles(readDir(procVTK/dirs[i], fileName::FILE));
1280 
1281  forAll(subFiles, j)
1282  {
1283  fileName procFile(procVTK/dirs[i]/subFiles[j]);
1284 
1285  if (exists(procFile))
1286  {
1287  string cmd
1288  (
1289  "ln -s "
1290  + procFile
1291  + " "
1292  + "processor"
1293  + name(Pstream::myProcNo())
1294  + "_"
1295  + procFile.name()
1296  );
1297  if (system(cmd.c_str()) == -1)
1298  {
1300  << "Could not execute command " << cmd << endl;
1301  }
1302  }
1303  }
1304  }
1305  }
1306 
1307  Info<< "End\n" << endl;
1308 
1309  return 0;
1310 }
1311 
1312 
1313 // ************************************************************************* //
Foam::doubleScalar
double doubleScalar
Double precision floating point scalar type.
Definition: doubleScalar.H:49
Foam::sphericalTensor
SphericalTensor< scalar > sphericalTensor
SphericalTensor of scalars.
Definition: sphericalTensor.H:48
Foam::compressible::New
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Definition: turbulentFluidThermoModel.C:36
print
void print(const char *msg, Ostream &os, const PtrList< GeoField > &flds)
Definition: foamToVTK.C:164
Foam::exists
bool exists(const fileName &, const bool checkGzip=true)
Does the name exist (as DIRECTORY or FILE) in the file system?
Definition: POSIX.C:608
Foam::labelList
List< label > labelList
A List of labels.
Definition: labelList.H:56
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tensor
Tensor< scalar > tensor
Tensor of scalars.
Definition: tensor.H:51
Foam::indirectPrimitivePatch
PrimitivePatch< face, IndirectList, const pointField & > indirectPrimitivePatch
Foam::indirectPrimitivePatch.
Definition: indirectPrimitivePatch.H:45
Cloud.H
tensorIOField.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
scalarIOField.H
Foam::instantList
List< instant > instantList
List of instants.
Definition: instantList.H:42
Foam::writeFaceSet
void writeFaceSet(const bool binary, const vtkMesh &vMesh, const faceSet &set, const fileName &fileName)
Definition: writeFaceSet.C:33
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
Foam::wordList
List< word > wordList
A List of words.
Definition: fileName.H:54
Foam::faceZoneMesh
ZoneMesh< faceZone, polyMesh > faceZoneMesh
A ZoneMesh with the type faceZone.
Definition: faceZoneMeshFwd.H:41
n
label n
Definition: TABSMDCalcMethod2.H:31
regionName
Foam::word regionName
Definition: createNamedDynamicFvMesh.H:1
Foam::label
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
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
passiveParticle.H
addRegionOption.H
Foam::symmTensor
SymmTensor< scalar > symmTensor
SymmTensor of scalars.
Definition: symmTensor.H:48
size_type
graph_traits< Graph >::vertices_size_type size_type
Definition: SloanRenumber.C:73
Foam::fixed
IOstream & fixed(IOstream &io)
Definition: IOstream.H:576
Foam::FatalError
error FatalError
createNamedMesh.H
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
faceZoneMesh.H
Foam::faceZoneMesh.
Foam::isDir
bool isDir(const fileName &)
Does the name exist as a DIRECTORY in the file system?
Definition: POSIX.C:615
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
emptyPolyPatch.H
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::writeSurfFields
void writeSurfFields(const bool binary, const vtkMesh &vMesh, const fileName &fileName, const PtrList< surfaceVectorField > &surfVectorFields)
Definition: writeSurfFields.C:41
Foam::fileNameList
List< fileName > fileNameList
A List of fileNames.
Definition: fileNameList.H:50
volPointInterpolation.H
Foam::chDir
bool chDir(const fileName &dir)
Change the current directory to the one given and return true,.
Definition: POSIX.C:264
setRootCase.H
timeDirs
static instantList timeDirs
Definition: globalFoam.H:44
Foam::findStrings
bool findStrings(const wordReListMatcher &matcher, const std::string &str)
Return true if string matches one of the regular expressions.
Definition: stringListOps.H:52
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
readFields
void readFields(const boolList &haveMesh, const fvMesh &mesh, const autoPtr< fvMeshSubset > &subsetterPtr, IOobjectList &allObjects, PtrList< GeoField > &fields)
Definition: redistributePar.C:589
createTime.H
cloudName
const word cloudName(propsDict.lookup("cloudName"))
Foam::argList::optionFound
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:108
symmTensorIOField.H
fvCFD.H
stringListOps.H
Operations on lists of strings.
patches
patches[0]
Definition: createSingleCellMesh.H:36
List
Definition: Test.C:19
main
int main(int argc, char *argv[])
Definition: foamToVTK.C:228
sphericalTensorIOField.H
Foam::readDir
fileNameList readDir(const fileName &, const fileName::Type=fileName::FILE, const bool filtergz=true)
Read a directory and return the entries as a string list.
Definition: POSIX.C:660
Foam::surfaceVectorField
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
Definition: surfaceFieldsFwd.H:55
args
Foam::argList args(argc, argv)
Foam::mkDir
bool mkDir(const fileName &, mode_t=0777)
Make a directory and return an error if it could not be created.
Definition: POSIX.C:419
Foam::argList::optionReadIfPresent
bool optionReadIfPresent(const word &opt, T &) const
Read a value from the named option if present.
Definition: argListI.H:198
getSelectedPatches
labelList getSelectedPatches(const polyBoundaryMesh &patches, const List< wordRe > &excludePatches)
Definition: foamToVTK.C:189
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:259
Foam::system
int system(const std::string &command)
Execute the specified command.
Definition: POSIX.C:1155
Foam::rmDir
bool rmDir(const fileName &)
Remove a dirctory and its contents.
Definition: POSIX.C:974
pointMesh.H
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Foam::writePointSet
void writePointSet(const bool binary, const primitiveMesh &mesh, const topoSet &set, const fileName &fileName)
Write pointSet to vtk polydata file.
labelIOField.H
Foam::argList::optionLookup
IStringStream optionLookup(const word &opt) const
Return an IStringStream from the named option.
Definition: argListI.H:114