BUG: UListIO: byteSize overflowing on really big faceLists
[OpenFOAM-2.0.x.git] / applications / utilities / postProcessing / dataConversion / foamToVTK / foamToVTK.C
blob652591d38e70a1e24abc9afb3ab3e428b709eb95
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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.
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.
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/>.
24 Application
25     foamToVTK
27 Description
28     Legacy VTK file format writer.
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.
39 Usage
41     - foamToVTK [OPTION]
43     \param -ascii \n
44     Write VTK data in ASCII format instead of binary.
46     \param -mesh \<name\>\n
47     Use a different mesh name (instead of -region)
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.
57     \param -surfaceFields \n
58     Write surfaceScalarFields (e.g., phi)
60     \param -cellSet \<name\>\n
61     \param -faceSet \<name\>\n
62     \param -pointSet \<name\>\n
63     Restrict conversion to the cellSet, faceSet or pointSet.
65     \param -nearCellValue \n
66     Output cell value on patches instead of patch value itself
68     \param -noInternal \n
69     Do not generate file for mesh, only for patches
71     \param -noPointValues \n
72     No pointFields
74     \param -noFaceZones \n
75     No faceZones
77     \param -noLinks \n
78     (in parallel) do not link processor files to master
80     \param poly \n
81     write polyhedral cells without tet/pyramid decomposition
83     \param -allPatches \n
84     Combine all patches into a single file
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.
95     \param -useTimeName \n
96     use the time index in the VTK file name instead of the time index
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.
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:
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
134 \*---------------------------------------------------------------------------*/
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"
150 #include "vtkMesh.H"
151 #include "readFields.H"
152 #include "writeFuns.H"
154 #include "internalWriter.H"
155 #include "patchWriter.H"
156 #include "lagrangianWriter.H"
158 #include "writeFaceSet.H"
159 #include "writePointSet.H"
160 #include "surfaceMeshWriter.H"
161 #include "writeSurfFields.H"
164 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
166 template<class GeoField>
167 void print(const char* msg, Ostream& os, const PtrList<GeoField>& flds)
169     if (flds.size())
170     {
171         os  << msg;
172         forAll(flds, i)
173         {
174             os  << ' ' << flds[i].name();
175         }
176         os  << endl;
177     }
181 void print(Ostream& os, const wordList& flds)
183     forAll(flds, i)
184     {
185         os  << ' ' << flds[i];
186     }
187     os  << endl;
191 labelList getSelectedPatches
193     const polyBoundaryMesh& patches,
194     const List<wordRe>& excludePatches  //HashSet<word>& excludePatches
197     DynamicList<label> patchIDs(patches.size());
199     Info<< "Combining patches:" << endl;
201     forAll(patches, patchI)
202     {
203         const polyPatch& pp = patches[patchI];
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();
229 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
230 // Main program:
232 int main(int argc, char *argv[])
234     argList::addNote
235     (
236         "legacy VTK file format writer"
237     );
238     timeSelector::addOptions();
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     );
322     #include "setRootCase.H"
323     #include "createTime.H"
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");
331     // decomposition of polyhedral cells into tets/pyramids cells
332     vtkTopo::decomposePoly     = !args.optionFound("poly");
334     if (binary && (sizeof(floatScalar) != 4 || sizeof(label) != 4))
335     {
336         FatalErrorIn(args.executable())
337             << "floatScalar and/or label are not 4 bytes in size" << nl
338             << "Hence cannot use binary VTK format. Please use -ascii"
339             << exit(FatalError);
340     }
342     const bool nearCellValue = args.optionFound("nearCellValue");
344     if (nearCellValue)
345     {
346         WarningIn(args.executable())
347             << "Using neighbouring cell value instead of patch value"
348             << nl << endl;
349     }
351     const bool noPointValues = args.optionFound("noPointValues");
353     if (noPointValues)
354     {
355         WarningIn(args.executable())
356             << "Outputting cell values only" << nl << endl;
357     }
359     const bool allPatches = args.optionFound("allPatches");
361     List<wordRe> excludePatches;
362     if (args.optionFound("excludePatches"))
363     {
364         args.optionLookup("excludePatches")() >> excludePatches;
366         Info<< "Not including patches " << excludePatches << nl << endl;
367     }
369     word cellSetName;
370     string vtkName = runTime.caseName();
372     if (args.optionReadIfPresent("cellSet", cellSetName))
373     {
374         vtkName = cellSetName;
375     }
376     else if (Pstream::parRun())
377     {
378         // Strip off leading casename, leaving just processor_DDD ending.
379         string::size_type i = vtkName.rfind("processor");
381         if (i != string::npos)
382         {
383             vtkName = vtkName.substr(i);
384         }
385     }
388     instantList timeDirs = timeSelector::select0(runTime, args);
390 #   include "createNamedMesh.H"
392     // VTK/ directory in the case
393     fileName fvPath(runTime.path()/"VTK");
394     // Directory of mesh (region0 gets filtered out)
395     fileName regionPrefix = "";
397     if (regionName != polyMesh::defaultRegion)
398     {
399         fvPath = fvPath/regionName;
400         regionPrefix = regionName;
401     }
403     if (isDir(fvPath))
404     {
405         if
406         (
407             args.optionFound("time")
408          || args.optionFound("latestTime")
409          || cellSetName.size()
410          || regionName != polyMesh::defaultRegion
411         )
412         {
413             Info<< "Keeping old VTK files in " << fvPath << nl << endl;
414         }
415         else
416         {
417             Info<< "Deleting old VTK files in " << fvPath << nl << endl;
419             rmDir(fvPath);
420         }
421     }
423     mkDir(fvPath);
426     // mesh wrapper; does subsetting and decomposition
427     vtkMesh vMesh(mesh, cellSetName);
430     // Scan for all possible lagrangian clouds
431     HashSet<fileName> allCloudDirs;
432     forAll(timeDirs, timeI)
433     {
434         runTime.setTime(timeDirs[timeI], timeI);
435         fileNameList cloudDirs
436         (
437             readDir
438             (
439                 runTime.timePath()/regionPrefix/cloud::prefix,
440                 fileName::DIRECTORY
441             )
442         );
443         forAll(cloudDirs, i)
444         {
445             IOobjectList sprayObjs
446             (
447                 mesh,
448                 runTime.timeName(),
449                 cloud::prefix/cloudDirs[i]
450             );
452             IOobject* positionsPtr = sprayObjs.lookup("positions");
454             if (positionsPtr)
455             {
456                 if (allCloudDirs.insert(cloudDirs[i]))
457                 {
458                     Info<< "At time: " << runTime.timeName()
459                         << " detected cloud directory : " << cloudDirs[i]
460                         << endl;
461                 }
462             }
463         }
464     }
467     forAll(timeDirs, timeI)
468     {
469         runTime.setTime(timeDirs[timeI], timeI);
471         Info<< "Time: " << runTime.timeName() << endl;
473         word timeDesc =
474             useTimeName ? runTime.timeName() : Foam::name(runTime.timeIndex());
476         // Check for new polyMesh/ and update mesh, fvMeshSubset and cell
477         // decomposition.
478         polyMesh::readUpdateState meshState = vMesh.readUpdate();
480         const fvMesh& mesh = vMesh.mesh();
482         if
483         (
484             meshState == polyMesh::TOPO_CHANGE
485          || meshState == polyMesh::TOPO_PATCH_CHANGE
486         )
487         {
488             Info<< "    Read new mesh" << nl << endl;
489         }
492         // If faceSet: write faceSet only (as polydata)
493         if (args.optionFound("faceSet"))
494         {
495             // Load the faceSet
496             faceSet set(mesh, args["faceSet"]);
498             // Filename as if patch with same name.
499             mkDir(fvPath/set.name());
501             fileName patchFileName
502             (
503                 fvPath/set.name()/set.name()
504               + "_"
505               + timeDesc
506               + ".vtk"
507             );
509             Info<< "    FaceSet   : " << patchFileName << endl;
511             writeFaceSet(binary, vMesh, set, patchFileName);
513             continue;
514         }
515         // If pointSet: write pointSet only (as polydata)
516         if (args.optionFound("pointSet"))
517         {
518             // Load the pointSet
519             pointSet set(mesh, args["pointSet"]);
521             // Filename as if patch with same name.
522             mkDir(fvPath/set.name());
524             fileName patchFileName
525             (
526                 fvPath/set.name()/set.name()
527               + "_"
528               + timeDesc
529               + ".vtk"
530             );
532             Info<< "    pointSet   : " << patchFileName << endl;
534             writePointSet(binary, vMesh, set, patchFileName);
536             continue;
537         }
540         // Search for list of objects for this time
541         IOobjectList objects(mesh, runTime.timeName());
543         HashSet<word> selectedFields;
544         args.optionReadIfPresent("fields", selectedFields);
546         // Construct the vol fields (on the original mesh if subsetted)
548         PtrList<volScalarField> vsf;
549         readFields(vMesh, vMesh.baseMesh(), objects, selectedFields, vsf);
550         print("    volScalarFields            :", Info, vsf);
552         PtrList<volVectorField> vvf;
553         readFields(vMesh, vMesh.baseMesh(), objects, selectedFields, vvf);
554         print("    volVectorFields            :", Info, vvf);
556         PtrList<volSphericalTensorField> vSpheretf;
557         readFields(vMesh, vMesh.baseMesh(), objects, selectedFields, vSpheretf);
558         print("    volSphericalTensorFields   :", Info, vSpheretf);
560         PtrList<volSymmTensorField> vSymmtf;
561         readFields(vMesh, vMesh.baseMesh(), objects, selectedFields, vSymmtf);
562         print("    volSymmTensorFields        :", Info, vSymmtf);
564         PtrList<volTensorField> vtf;
565         readFields(vMesh, vMesh.baseMesh(), objects, selectedFields, vtf);
566         print("    volTensorFields            :", Info, vtf);
568         label nVolFields =
569                 vsf.size()
570               + vvf.size()
571               + vSpheretf.size()
572               + vSymmtf.size()
573               + vtf.size();
576         // Construct pointMesh only if nessecary since constructs edge
577         // addressing (expensive on polyhedral meshes)
578         if (noPointValues)
579         {
580             Info<< "    pointScalarFields : switched off"
581                 << " (\"-noPointValues\" (at your option)\n";
582             Info<< "    pointVectorFields : switched off"
583                 << " (\"-noPointValues\" (at your option)\n";
584         }
586         PtrList<pointScalarField> psf;
587         PtrList<pointVectorField> pvf;
588         PtrList<pointSphericalTensorField> pSpheretf;
589         PtrList<pointSymmTensorField> pSymmtf;
590         PtrList<pointTensorField> ptf;
592         if (!noPointValues)
593         {
594             readFields
595             (
596                 vMesh,
597                 pointMesh::New(vMesh.baseMesh()),
598                 objects,
599                 selectedFields,
600                 psf
601             );
602             print("    pointScalarFields          :", Info, psf);
604             readFields
605             (
606                 vMesh,
607                 pointMesh::New(vMesh.baseMesh()),
608                 objects,
609                 selectedFields,
610                 pvf
611             );
612             print("    pointVectorFields          :", Info, pvf);
614             readFields
615             (
616                 vMesh,
617                 pointMesh::New(vMesh.baseMesh()),
618                 objects,
619                 selectedFields,
620                 pSpheretf
621             );
622             print("    pointSphericalTensorFields :", Info, pSpheretf);
624             readFields
625             (
626                 vMesh,
627                 pointMesh::New(vMesh.baseMesh()),
628                 objects,
629                 selectedFields,
630                 pSymmtf
631             );
632             print("    pointSymmTensorFields      :", Info, pSymmtf);
634             readFields
635             (
636                 vMesh,
637                 pointMesh::New(vMesh.baseMesh()),
638                 objects,
639                 selectedFields,
640                 ptf
641             );
642             print("    pointTensorFields          :", Info, ptf);
643         }
644         Info<< endl;
646         label nPointFields =
647             psf.size()
648           + pvf.size()
649           + pSpheretf.size()
650           + pSymmtf.size()
651           + ptf.size();
653         if (doWriteInternal)
654         {
655             //
656             // Create file and write header
657             //
658             fileName vtkFileName
659             (
660                 fvPath/vtkName
661               + "_"
662               + timeDesc
663               + ".vtk"
664             );
666             Info<< "    Internal  : " << vtkFileName << endl;
668             // Write mesh
669             internalWriter writer(vMesh, binary, vtkFileName);
671             // VolFields + cellID
672             writeFuns::writeCellDataHeader
673             (
674                 writer.os(),
675                 vMesh.nFieldCells(),
676                 1+nVolFields
677             );
679             // Write cellID field
680             writer.writeCellIDs();
682             // Write volFields
683             writer.write(vsf);
684             writer.write(vvf);
685             writer.write(vSpheretf);
686             writer.write(vSymmtf);
687             writer.write(vtf);
689             if (!noPointValues)
690             {
691                 writeFuns::writePointDataHeader
692                 (
693                     writer.os(),
694                     vMesh.nFieldPoints(),
695                     nVolFields+nPointFields
696                 );
698                 // pointFields
699                 writer.write(psf);
700                 writer.write(pvf);
701                 writer.write(pSpheretf);
702                 writer.write(pSymmtf);
703                 writer.write(ptf);
705                 // Interpolated volFields
706                 volPointInterpolation pInterp(mesh);
707                 writer.write(pInterp, vsf);
708                 writer.write(pInterp, vvf);
709                 writer.write(pInterp, vSpheretf);
710                 writer.write(pInterp, vSymmtf);
711                 writer.write(pInterp, vtf);
712             }
713         }
715         //---------------------------------------------------------------------
716         //
717         // Write surface fields
718         //
719         //---------------------------------------------------------------------
721         if (args.optionFound("surfaceFields"))
722         {
723             PtrList<surfaceScalarField> ssf;
724             readFields
725             (
726                 vMesh,
727                 vMesh.baseMesh(),
728                 objects,
729                 selectedFields,
730                 ssf
731             );
732             print("    surfScalarFields  :", Info, ssf);
734             PtrList<surfaceVectorField> svf;
735             readFields
736             (
737                 vMesh,
738                 vMesh.baseMesh(),
739                 objects,
740                 selectedFields,
741                 svf
742             );
743             print("    surfVectorFields  :", Info, svf);
745             if (ssf.size() + svf.size() > 0)
746             {
747                 // Rework the scalar fields into vectorfields.
748                 label sz = svf.size();
750                 svf.setSize(sz+ssf.size());
752                 surfaceVectorField n(mesh.Sf()/mesh.magSf());
754                 forAll(ssf, i)
755                 {
756                     svf.set(sz+i, ssf[i]*n);
757                     svf[sz+i].rename(ssf[i].name());
758                 }
759                 ssf.clear();
761                 mkDir(fvPath / "surfaceFields");
763                 fileName surfFileName
764                 (
765                     fvPath
766                    /"surfaceFields"
767                    /"surfaceFields"
768                    + "_"
769                    + timeDesc
770                    + ".vtk"
771                 );
773                 writeSurfFields
774                 (
775                     binary,
776                     vMesh,
777                     surfFileName,
778                     svf
779                 );
780             }
781         }
784         //---------------------------------------------------------------------
785         //
786         // Write patches (POLYDATA file, one for each patch)
787         //
788         //---------------------------------------------------------------------
790         const polyBoundaryMesh& patches = mesh.boundaryMesh();
792         if (allPatches)
793         {
794             mkDir(fvPath/"allPatches");
796             fileName patchFileName;
798             if (vMesh.useSubMesh())
799             {
800                 patchFileName =
801                     fvPath/"allPatches"/cellSetName
802                   + "_"
803                   + timeDesc
804                   + ".vtk";
805             }
806             else
807             {
808                 patchFileName =
809                     fvPath/"allPatches"/"allPatches"
810                   + "_"
811                   + timeDesc
812                   + ".vtk";
813             }
815             Info<< "    Combined patches     : " << patchFileName << endl;
817             patchWriter writer
818             (
819                 vMesh,
820                 binary,
821                 nearCellValue,
822                 patchFileName,
823                 getSelectedPatches(patches, excludePatches)
824             );
826             // VolFields + patchID
827             writeFuns::writeCellDataHeader
828             (
829                 writer.os(),
830                 writer.nFaces(),
831                 1+nVolFields
832             );
834             // Write patchID field
835             writer.writePatchIDs();
837             // Write volFields
838             writer.write(vsf);
839             writer.write(vvf);
840             writer.write(vSpheretf);
841             writer.write(vSymmtf);
842             writer.write(vtf);
844             if (!noPointValues)
845             {
846                 writeFuns::writePointDataHeader
847                 (
848                     writer.os(),
849                     writer.nPoints(),
850                     nPointFields
851                 );
853                 // Write pointFields
854                 writer.write(psf);
855                 writer.write(pvf);
856                 writer.write(pSpheretf);
857                 writer.write(pSymmtf);
858                 writer.write(ptf);
860                 // no interpolated volFields since I cannot be bothered to
861                 // create the patchInterpolation for all subpatches.
862             }
863         }
864         else
865         {
866             forAll(patches, patchI)
867             {
868                 const polyPatch& pp = patches[patchI];
870                 if (!findStrings(excludePatches, pp.name()))
871                 {
872                     mkDir(fvPath/pp.name());
874                     fileName patchFileName;
876                     if (vMesh.useSubMesh())
877                     {
878                         patchFileName =
879                             fvPath/pp.name()/cellSetName
880                           + "_"
881                           + timeDesc
882                           + ".vtk";
883                     }
884                     else
885                     {
886                         patchFileName =
887                             fvPath/pp.name()/pp.name()
888                           + "_"
889                           + timeDesc
890                           + ".vtk";
891                     }
893                     Info<< "    Patch     : " << patchFileName << endl;
895                     patchWriter writer
896                     (
897                         vMesh,
898                         binary,
899                         nearCellValue,
900                         patchFileName,
901                         labelList(1, patchI)
902                     );
904                     if (!isA<emptyPolyPatch>(pp))
905                     {
906                         // VolFields + patchID
907                         writeFuns::writeCellDataHeader
908                         (
909                             writer.os(),
910                             writer.nFaces(),
911                             1+nVolFields
912                         );
914                         // Write patchID field
915                         writer.writePatchIDs();
917                         // Write volFields
918                         writer.write(vsf);
919                         writer.write(vvf);
920                         writer.write(vSpheretf);
921                         writer.write(vSymmtf);
922                         writer.write(vtf);
924                         if (!noPointValues)
925                         {
926                             writeFuns::writePointDataHeader
927                             (
928                                 writer.os(),
929                                 writer.nPoints(),
930                                 nVolFields
931                               + nPointFields
932                             );
934                             // Write pointFields
935                             writer.write(psf);
936                             writer.write(pvf);
937                             writer.write(pSpheretf);
938                             writer.write(pSymmtf);
939                             writer.write(ptf);
941                             PrimitivePatchInterpolation<primitivePatch> pInter
942                             (
943                                 pp
944                             );
946                             // Write interpolated volFields
947                             writer.write(pInter, vsf);
948                             writer.write(pInter, vvf);
949                             writer.write(pInter, vSpheretf);
950                             writer.write(pInter, vSymmtf);
951                             writer.write(pInter, vtf);
952                         }
953                     }
954                 }
955             }
956         }
958         //---------------------------------------------------------------------
959         //
960         // Write faceZones (POLYDATA file, one for each zone)
961         //
962         //---------------------------------------------------------------------
964         if (doFaceZones)
965         {
966             PtrList<surfaceScalarField> ssf;
967             readFields
968             (
969                 vMesh,
970                 vMesh.baseMesh(),
971                 objects,
972                 selectedFields,
973                 ssf
974             );
975             print("    surfScalarFields  :", Info, ssf);
977             PtrList<surfaceVectorField> svf;
978             readFields
979             (
980                 vMesh,
981                 vMesh.baseMesh(),
982                 objects,
983                 selectedFields,
984                 svf
985             );
986             print("    surfVectorFields  :", Info, svf);
988             const faceZoneMesh& zones = mesh.faceZones();
990             forAll(zones, zoneI)
991             {
992                 const faceZone& fz = zones[zoneI];
994                 mkDir(fvPath/fz.name());
996                 fileName patchFileName;
998                 if (vMesh.useSubMesh())
999                 {
1000                     patchFileName =
1001                         fvPath/fz.name()/cellSetName
1002                       + "_"
1003                       + timeDesc
1004                       + ".vtk";
1005                 }
1006                 else
1007                 {
1008                     patchFileName =
1009                         fvPath/fz.name()/fz.name()
1010                       + "_"
1011                       + timeDesc
1012                       + ".vtk";
1013                 }
1015                 Info<< "    FaceZone  : " << patchFileName << endl;
1017                 indirectPrimitivePatch pp
1018                 (
1019                     IndirectList<face>(mesh.faces(), fz),
1020                     mesh.points()
1021                 );
1023                 surfaceMeshWriter writer
1024                 (
1025                     vMesh,
1026                     binary,
1027                     pp,
1028                     fz.name(),
1029                     patchFileName
1030                 );
1032                 // Number of fields
1033                 writeFuns::writeCellDataHeader
1034                 (
1035                     writer.os(),
1036                     pp.size(),
1037                     ssf.size()+svf.size()
1038                 );
1040                 writer.write(ssf);
1041                 writer.write(svf);
1042             }
1043         }
1047         //---------------------------------------------------------------------
1048         //
1049         // Write lagrangian data
1050         //
1051         //---------------------------------------------------------------------
1053         forAllConstIter(HashSet<fileName>, allCloudDirs, iter)
1054         {
1055             const fileName& cloudName = iter.key();
1057             // Always create the cloud directory.
1058             mkDir(fvPath/cloud::prefix/cloudName);
1060             fileName lagrFileName
1061             (
1062                 fvPath/cloud::prefix/cloudName/cloudName
1063               + "_" + timeDesc + ".vtk"
1064             );
1066             Info<< "    Lagrangian: " << lagrFileName << endl;
1069             IOobjectList sprayObjs
1070             (
1071                 mesh,
1072                 runTime.timeName(),
1073                 cloud::prefix/cloudName
1074             );
1076             IOobject* positionsPtr = sprayObjs.lookup("positions");
1078             if (positionsPtr)
1079             {
1080                 wordList labelNames(sprayObjs.names(labelIOField::typeName));
1081                 Info<< "        labels            :";
1082                 print(Info, labelNames);
1084                 wordList scalarNames(sprayObjs.names(scalarIOField::typeName));
1085                 Info<< "        scalars           :";
1086                 print(Info, scalarNames);
1088                 wordList vectorNames(sprayObjs.names(vectorIOField::typeName));
1089                 Info<< "        vectors           :";
1090                 print(Info, vectorNames);
1092                 wordList sphereNames
1093                 (
1094                     sprayObjs.names
1095                     (
1096                         sphericalTensorIOField::typeName
1097                     )
1098                 );
1099                 Info<< "        spherical tensors :";
1100                 print(Info, sphereNames);
1102                 wordList symmNames
1103                 (
1104                     sprayObjs.names
1105                     (
1106                         symmTensorIOField::typeName
1107                     )
1108                 );
1109                 Info<< "        symm tensors      :";
1110                 print(Info, symmNames);
1112                 wordList tensorNames(sprayObjs.names(tensorIOField::typeName));
1113                 Info<< "        tensors           :";
1114                 print(Info, tensorNames);
1116                 lagrangianWriter writer
1117                 (
1118                     vMesh,
1119                     binary,
1120                     lagrFileName,
1121                     cloudName,
1122                     false
1123                 );
1125                 // Write number of fields
1126                 writer.writeParcelHeader
1127                 (
1128                     labelNames.size()
1129                   + scalarNames.size()
1130                   + vectorNames.size()
1131                   + sphereNames.size()
1132                   + symmNames.size()
1133                   + tensorNames.size()
1134                 );
1136                 // Fields
1137                 writer.writeIOField<label>(labelNames);
1138                 writer.writeIOField<scalar>(scalarNames);
1139                 writer.writeIOField<vector>(vectorNames);
1140                 writer.writeIOField<sphericalTensor>(sphereNames);
1141                 writer.writeIOField<symmTensor>(symmNames);
1142                 writer.writeIOField<tensor>(tensorNames);
1143             }
1144             else
1145             {
1146                 lagrangianWriter writer
1147                 (
1148                     vMesh,
1149                     binary,
1150                     lagrFileName,
1151                     cloudName,
1152                     true
1153                 );
1155                 // Write number of fields
1156                 writer.writeParcelHeader(0);
1157             }
1158         }
1159     }
1162     //---------------------------------------------------------------------
1163     //
1164     // Link parallel outputs back to undecomposed case for ease of loading
1165     //
1166     //---------------------------------------------------------------------
1168     if (Pstream::parRun() && doLinks)
1169     {
1170         mkDir(runTime.path()/".."/"VTK");
1171         chDir(runTime.path()/".."/"VTK");
1173         Info<< "Linking all processor files to " << runTime.path()/".."/"VTK"
1174             << endl;
1176         // Get list of vtk files
1177         fileName procVTK
1178         (
1179             fileName("..")
1180            /"processor" + name(Pstream::myProcNo())
1181            /"VTK"
1182         );
1184         fileNameList dirs(readDir(procVTK, fileName::DIRECTORY));
1185         label sz = dirs.size();
1186         dirs.setSize(sz+1);
1187         dirs[sz] = ".";
1189         forAll(dirs, i)
1190         {
1191             fileNameList subFiles(readDir(procVTK/dirs[i], fileName::FILE));
1193             forAll(subFiles, j)
1194             {
1195                 fileName procFile(procVTK/dirs[i]/subFiles[j]);
1197                 if (exists(procFile))
1198                 {
1199                     string cmd
1200                     (
1201                         "ln -s "
1202                       + procFile
1203                       + " "
1204                       + "processor"
1205                       + name(Pstream::myProcNo())
1206                       + "_"
1207                       + procFile.name()
1208                     );
1209                     if (system(cmd.c_str()) == -1)
1210                     {
1211                         WarningIn(args.executable())
1212                             << "Could not execute command " << cmd << endl;
1213                     }
1214                 }
1215             }
1216         }
1217     }
1219     Info<< "End\n" << endl;
1221     return 0;
1225 // ************************************************************************* //