1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
7 -------------------------------------------------------------------------------
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
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/>.
28 Legacy VTK file format writer.
30 - handles volScalar, volVector, pointScalar, pointVector, surfaceScalar
33 - both ascii and binary.
34 - single time step writing.
36 - automatic decomposition of cells; polygons on boundary undecomposed since
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,
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
69 Do not generate file for mesh, only for patches
71 \param -noPointValues \n
74 \param -noFaceZones \n
78 (in parallel) do not link processor files to master
81 write polyhedral cells without tet/pyramid decomposition
84 Combine all patches into a single file
86 \param -excludePatches \<patchNames\>\n
87 Specify patches (wildcards) to exclude. For example,
89 -excludePatches '( inlet_1 inlet_2 "proc.*")'
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
95 \param -useTimeName \n
96 use the time index in the VTK file name instead of the time index
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.
109 no automatic timestep recognition.
110 However can have .pvd file format which refers to time simulation
111 if XML *.vtu files are available:
114 <?xml version="1.0"?>
115 <VTKFile type="Collection"
117 byte_order="LittleEndian"
118 compressor="vtkZLibDataCompressor">
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"/>
134 \*---------------------------------------------------------------------------*/
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"
147 #include "passiveParticle.H"
148 #include "stringListOps.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)
174 os << ' ' << flds[i].name();
181 void print(Ostream& os, const wordList& flds)
185 os << ' ' << flds[i];
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)
203 const polyPatch& pp = patches[patchI];
207 isType<emptyPolyPatch>(pp)
208 || (Pstream::parRun() && isType<processorPolyPatch>(pp))
211 Info<< " discarding empty/processor patch " << patchI
212 << " " << pp.name() << endl;
214 else if (findStrings(excludePatches, pp.name()))
216 Info<< " excluding patch " << patchI
217 << " " << pp.name() << endl;
221 patchIDs.append(patchI);
222 Info<< " patch " << patchI << " " << pp.name() << endl;
225 return patchIDs.shrink();
229 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
232 int main(int argc, char *argv[])
236 "legacy VTK file format writer"
238 timeSelector::addOptions();
240 #include "addRegionOption.H"
245 "only convert the specified fields - eg '(p T U)'"
251 "convert a mesh subset corresponding to the specified cellSet"
257 "restrict conversion to the specified faceSet"
263 "restrict conversion to the specified pointSet"
265 argList::addBoolOption
268 "write in ASCII format instead of binary"
270 argList::addBoolOption
273 "write polyhedral cells without tet/pyramid decomposition"
275 argList::addBoolOption
278 "write surfaceScalarFields (e.g., phi)"
280 argList::addBoolOption
283 "use cell value on patches instead of patch value itself"
285 argList::addBoolOption
288 "do not generate file for mesh, only for patches"
290 argList::addBoolOption
295 argList::addBoolOption
298 "combine all patches into a single file"
304 "a list of patches to exclude - eg '( inlet \".*Wall\" )' "
306 argList::addBoolOption
311 argList::addBoolOption
314 "don't link processor VTK files - parallel only"
316 argList::addBoolOption
319 "use the time name instead of the time index when naming the files"
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))
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"
342 const bool nearCellValue = args.optionFound("nearCellValue");
346 WarningIn(args.executable())
347 << "Using neighbouring cell value instead of patch value"
351 const bool noPointValues = args.optionFound("noPointValues");
355 WarningIn(args.executable())
356 << "Outputting cell values only" << nl << endl;
359 const bool allPatches = args.optionFound("allPatches");
361 List<wordRe> excludePatches;
362 if (args.optionFound("excludePatches"))
364 args.optionLookup("excludePatches")() >> excludePatches;
366 Info<< "Not including patches " << excludePatches << nl << endl;
370 string vtkName = runTime.caseName();
372 if (args.optionReadIfPresent("cellSet", cellSetName))
374 vtkName = cellSetName;
376 else if (Pstream::parRun())
378 // Strip off leading casename, leaving just processor_DDD ending.
379 string::size_type i = vtkName.rfind("processor");
381 if (i != string::npos)
383 vtkName = vtkName.substr(i);
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)
399 fvPath = fvPath/regionName;
400 regionPrefix = regionName;
407 args.optionFound("time")
408 || args.optionFound("latestTime")
409 || cellSetName.size()
410 || regionName != polyMesh::defaultRegion
413 Info<< "Keeping old VTK files in " << fvPath << nl << endl;
417 Info<< "Deleting old VTK files in " << fvPath << nl << endl;
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)
434 runTime.setTime(timeDirs[timeI], timeI);
435 fileNameList cloudDirs
439 runTime.timePath()/regionPrefix/cloud::prefix,
445 IOobjectList sprayObjs
449 cloud::prefix/cloudDirs[i]
452 IOobject* positionsPtr = sprayObjs.lookup("positions");
456 if (allCloudDirs.insert(cloudDirs[i]))
458 Info<< "At time: " << runTime.timeName()
459 << " detected cloud directory : " << cloudDirs[i]
467 forAll(timeDirs, timeI)
469 runTime.setTime(timeDirs[timeI], timeI);
471 Info<< "Time: " << runTime.timeName() << endl;
474 useTimeName ? runTime.timeName() : Foam::name(runTime.timeIndex());
476 // Check for new polyMesh/ and update mesh, fvMeshSubset and cell
478 polyMesh::readUpdateState meshState = vMesh.readUpdate();
480 const fvMesh& mesh = vMesh.mesh();
484 meshState == polyMesh::TOPO_CHANGE
485 || meshState == polyMesh::TOPO_PATCH_CHANGE
488 Info<< " Read new mesh" << nl << endl;
492 // If faceSet: write faceSet only (as polydata)
493 if (args.optionFound("faceSet"))
496 faceSet set(mesh, args["faceSet"]);
498 // Filename as if patch with same name.
499 mkDir(fvPath/set.name());
501 fileName patchFileName
503 fvPath/set.name()/set.name()
509 Info<< " FaceSet : " << patchFileName << endl;
511 writeFaceSet(binary, vMesh, set, patchFileName);
515 // If pointSet: write pointSet only (as polydata)
516 if (args.optionFound("pointSet"))
519 pointSet set(mesh, args["pointSet"]);
521 // Filename as if patch with same name.
522 mkDir(fvPath/set.name());
524 fileName patchFileName
526 fvPath/set.name()/set.name()
532 Info<< " pointSet : " << patchFileName << endl;
534 writePointSet(binary, vMesh, set, patchFileName);
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);
576 // Construct pointMesh only if nessecary since constructs edge
577 // addressing (expensive on polyhedral meshes)
580 Info<< " pointScalarFields : switched off"
581 << " (\"-noPointValues\" (at your option)\n";
582 Info<< " pointVectorFields : switched off"
583 << " (\"-noPointValues\" (at your option)\n";
586 PtrList<pointScalarField> psf;
587 PtrList<pointVectorField> pvf;
588 PtrList<pointSphericalTensorField> pSpheretf;
589 PtrList<pointSymmTensorField> pSymmtf;
590 PtrList<pointTensorField> ptf;
597 pointMesh::New(vMesh.baseMesh()),
602 print(" pointScalarFields :", Info, psf);
607 pointMesh::New(vMesh.baseMesh()),
612 print(" pointVectorFields :", Info, pvf);
617 pointMesh::New(vMesh.baseMesh()),
622 print(" pointSphericalTensorFields :", Info, pSpheretf);
627 pointMesh::New(vMesh.baseMesh()),
632 print(" pointSymmTensorFields :", Info, pSymmtf);
637 pointMesh::New(vMesh.baseMesh()),
642 print(" pointTensorFields :", Info, ptf);
656 // Create file and write header
666 Info<< " Internal : " << vtkFileName << endl;
669 internalWriter writer(vMesh, binary, vtkFileName);
671 // VolFields + cellID
672 writeFuns::writeCellDataHeader
679 // Write cellID field
680 writer.writeCellIDs();
685 writer.write(vSpheretf);
686 writer.write(vSymmtf);
691 writeFuns::writePointDataHeader
694 vMesh.nFieldPoints(),
695 nVolFields+nPointFields
701 writer.write(pSpheretf);
702 writer.write(pSymmtf);
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);
715 //---------------------------------------------------------------------
717 // Write surface fields
719 //---------------------------------------------------------------------
721 if (args.optionFound("surfaceFields"))
723 PtrList<surfaceScalarField> ssf;
732 print(" surfScalarFields :", Info, ssf);
734 PtrList<surfaceVectorField> svf;
743 print(" surfVectorFields :", Info, svf);
745 if (ssf.size() + svf.size() > 0)
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());
756 svf.set(sz+i, ssf[i]*n);
757 svf[sz+i].rename(ssf[i].name());
761 mkDir(fvPath / "surfaceFields");
763 fileName surfFileName
784 //---------------------------------------------------------------------
786 // Write patches (POLYDATA file, one for each patch)
788 //---------------------------------------------------------------------
790 const polyBoundaryMesh& patches = mesh.boundaryMesh();
794 mkDir(fvPath/"allPatches");
796 fileName patchFileName;
798 if (vMesh.useSubMesh())
801 fvPath/"allPatches"/cellSetName
809 fvPath/"allPatches"/"allPatches"
815 Info<< " Combined patches : " << patchFileName << endl;
823 getSelectedPatches(patches, excludePatches)
826 // VolFields + patchID
827 writeFuns::writeCellDataHeader
834 // Write patchID field
835 writer.writePatchIDs();
840 writer.write(vSpheretf);
841 writer.write(vSymmtf);
846 writeFuns::writePointDataHeader
856 writer.write(pSpheretf);
857 writer.write(pSymmtf);
860 // no interpolated volFields since I cannot be bothered to
861 // create the patchInterpolation for all subpatches.
866 forAll(patches, patchI)
868 const polyPatch& pp = patches[patchI];
870 if (!findStrings(excludePatches, pp.name()))
872 mkDir(fvPath/pp.name());
874 fileName patchFileName;
876 if (vMesh.useSubMesh())
879 fvPath/pp.name()/cellSetName
887 fvPath/pp.name()/pp.name()
893 Info<< " Patch : " << patchFileName << endl;
904 if (!isA<emptyPolyPatch>(pp))
906 // VolFields + patchID
907 writeFuns::writeCellDataHeader
914 // Write patchID field
915 writer.writePatchIDs();
920 writer.write(vSpheretf);
921 writer.write(vSymmtf);
926 writeFuns::writePointDataHeader
937 writer.write(pSpheretf);
938 writer.write(pSymmtf);
941 PrimitivePatchInterpolation<primitivePatch> pInter
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);
958 //---------------------------------------------------------------------
960 // Write faceZones (POLYDATA file, one for each zone)
962 //---------------------------------------------------------------------
966 PtrList<surfaceScalarField> ssf;
975 print(" surfScalarFields :", Info, ssf);
977 PtrList<surfaceVectorField> svf;
986 print(" surfVectorFields :", Info, svf);
988 const faceZoneMesh& zones = mesh.faceZones();
992 const faceZone& fz = zones[zoneI];
994 mkDir(fvPath/fz.name());
996 fileName patchFileName;
998 if (vMesh.useSubMesh())
1001 fvPath/fz.name()/cellSetName
1009 fvPath/fz.name()/fz.name()
1015 Info<< " FaceZone : " << patchFileName << endl;
1017 indirectPrimitivePatch pp
1019 IndirectList<face>(mesh.faces(), fz),
1023 surfaceMeshWriter writer
1033 writeFuns::writeCellDataHeader
1037 ssf.size()+svf.size()
1047 //---------------------------------------------------------------------
1049 // Write lagrangian data
1051 //---------------------------------------------------------------------
1053 forAllConstIter(HashSet<fileName>, allCloudDirs, iter)
1055 const fileName& cloudName = iter.key();
1057 // Always create the cloud directory.
1058 mkDir(fvPath/cloud::prefix/cloudName);
1060 fileName lagrFileName
1062 fvPath/cloud::prefix/cloudName/cloudName
1063 + "_" + timeDesc + ".vtk"
1066 Info<< " Lagrangian: " << lagrFileName << endl;
1069 IOobjectList sprayObjs
1073 cloud::prefix/cloudName
1076 IOobject* positionsPtr = sprayObjs.lookup("positions");
1080 wordList labelNames(sprayObjs.names(labelIOField::typeName));
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
1096 sphericalTensorIOField::typeName
1099 Info<< " spherical tensors :";
1100 print(Info, sphereNames);
1106 symmTensorIOField::typeName
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
1125 // Write number of fields
1126 writer.writeParcelHeader
1129 + scalarNames.size()
1130 + vectorNames.size()
1131 + sphereNames.size()
1133 + tensorNames.size()
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);
1146 lagrangianWriter writer
1155 // Write number of fields
1156 writer.writeParcelHeader(0);
1162 //---------------------------------------------------------------------
1164 // Link parallel outputs back to undecomposed case for ease of loading
1166 //---------------------------------------------------------------------
1168 if (Pstream::parRun() && doLinks)
1170 mkDir(runTime.path()/".."/"VTK");
1171 chDir(runTime.path()/".."/"VTK");
1173 Info<< "Linking all processor files to " << runTime.path()/".."/"VTK"
1176 // Get list of vtk files
1180 /"processor" + name(Pstream::myProcNo())
1184 fileNameList dirs(readDir(procVTK, fileName::DIRECTORY));
1185 label sz = dirs.size();
1191 fileNameList subFiles(readDir(procVTK/dirs[i], fileName::FILE));
1195 fileName procFile(procVTK/dirs[i]/subFiles[j]);
1197 if (exists(procFile))
1205 + name(Pstream::myProcNo())
1209 if (system(cmd.c_str()) == -1)
1211 WarningIn(args.executable())
1212 << "Could not execute command " << cmd << endl;
1219 Info<< "End\n" << endl;
1225 // ************************************************************************* //