1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 1991-2009 OpenCFD Ltd.
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 the
13 Free Software Foundation; either version 2 of the License, or (at your
14 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, write to the Free Software Foundation,
23 Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
29 Tecplot binary file format writer.
33 - foamToTecplot360 [OPTION]
35 @param -fields \<fields\>\n
36 Convert selected fields only. For example,
40 The quoting is required to avoid shell expansions and to pass the
41 information as a single argument.
43 @param -cellSet \<name\>\n
44 @param -faceSet \<name\>\n
45 Restrict conversion to the cellSet, faceSet.
47 @param -nearCellValue \n
48 Output cell value on patches instead of patch value itself
51 Do not generate file for mesh, only for patches
53 @param -noPointValues \n
56 @param -noFaceZones \n
59 @param -excludePatches \<patchNames\>\n
60 Specify patches (wildcards) to exclude. For example,
62 -excludePatches '( inlet_1 inlet_2 "proc.*")'
64 The quoting is required to avoid shell expansions and to pass the
65 information as a single argument. The double quotes denote a regular
68 @param -useTimeName \n
69 use the time index in the VTK file name instead of the time index
71 \*---------------------------------------------------------------------------*/
73 #include "pointMesh.H"
74 #include "volPointInterpolation.H"
75 #include "emptyPolyPatch.H"
76 #include "labelIOField.H"
77 #include "scalarIOField.H"
78 #include "sphericalTensorIOField.H"
79 #include "symmTensorIOField.H"
80 #include "tensorIOField.H"
81 #include "passiveParticleCloud.H"
83 #include "stringListOps.H"
87 #include "readFields.H"
88 #include "tecplotWriter.H"
92 // Note: needs to be after TECIO to prevent Foam::Time conflicting with
96 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
98 template<class GeoField>
99 void print(const char* msg, Ostream& os, const PtrList<GeoField>& flds)
106 os<< ' ' << flds[i].name();
113 void print(Ostream& os, const wordList& flds)
123 labelList getSelectedPatches
125 const polyBoundaryMesh& patches,
126 const List<wordRe>& excludePatches //HashSet<word>& excludePatches
129 DynamicList<label> patchIDs(patches.size());
131 Info<< "Combining patches:" << endl;
133 forAll(patches, patchI)
135 const polyPatch& pp = patches[patchI];
139 isType<emptyPolyPatch>(pp)
140 || (Pstream::parRun() && isType<processorPolyPatch>(pp))
143 Info<< " discarding empty/processor patch " << patchI
144 << " " << pp.name() << endl;
146 else if (findStrings(excludePatches, pp.name()))
148 Info<< " excluding patch " << patchI
149 << " " << pp.name() << endl;
153 patchIDs.append(patchI);
154 Info<< " patch " << patchI << " " << pp.name() << endl;
157 return patchIDs.shrink();
165 int main(int argc, char *argv[])
167 timeSelector::addOptions();
169 # include "addRegionOption.H"
171 argList::validOptions.insert("fields", "fields");
172 argList::validOptions.insert("cellSet", "cellSet name");
173 argList::validOptions.insert("faceSet", "faceSet name");
174 argList::validOptions.insert("nearCellValue","");
175 argList::validOptions.insert("noInternal","");
176 argList::validOptions.insert("noPointValues","");
177 argList::validOptions.insert
180 "patches (wildcards) to exclude"
182 argList::validOptions.insert("noFaceZones","");
184 # include "setRootCase.H"
185 # include "createTime.H"
187 bool doWriteInternal = !args.optionFound("noInternal");
188 bool doFaceZones = !args.optionFound("noFaceZones");
190 bool nearCellValue = args.optionFound("nearCellValue");
194 WarningIn(args.executable())
195 << "Using neighbouring cell value instead of patch value"
199 bool noPointValues = args.optionFound("noPointValues");
203 WarningIn(args.executable())
204 << "Outputting cell values only" << nl << endl;
207 List<wordRe> excludePatches;
208 if (args.optionFound("excludePatches"))
210 args.optionLookup("excludePatches")() >> excludePatches;
212 Info<< "Not including patches " << excludePatches << nl << endl;
218 if (args.optionFound("cellSet"))
220 cellSetName = args.option("cellSet");
221 vtkName = cellSetName;
223 else if (Pstream::parRun())
225 // Strip off leading casename, leaving just processor_DDD ending.
226 vtkName = runTime.caseName();
228 string::size_type i = vtkName.rfind("processor");
230 if (i != string::npos)
232 vtkName = vtkName.substr(i);
237 vtkName = runTime.caseName();
241 instantList timeDirs = timeSelector::select0(runTime, args);
243 # include "createNamedMesh.H"
245 // TecplotData/ directory in the case
246 fileName fvPath(runTime.path()/"Tecplot360");
247 // Directory of mesh (region0 gets filtered out)
248 fileName regionPrefix = "";
250 if (regionName != polyMesh::defaultRegion)
252 fvPath = fvPath/regionName;
253 regionPrefix = regionName;
260 args.optionFound("time")
261 || args.optionFound("latestTime")
262 || cellSetName.size()
263 || regionName != polyMesh::defaultRegion
266 Info<< "Keeping old files in " << fvPath << nl << endl;
270 Info<< "Deleting old VTK files in " << fvPath << nl << endl;
279 // mesh wrapper; does subsetting and decomposition
280 vtkMesh vMesh(mesh, cellSetName);
282 forAll(timeDirs, timeI)
284 runTime.setTime(timeDirs[timeI], timeI);
286 Info<< "Time: " << runTime.timeName() << endl;
288 const word timeDesc = name(timeI); //name(runTime.timeIndex());
290 // Check for new polyMesh/ and update mesh, fvMeshSubset and cell
292 polyMesh::readUpdateState meshState = vMesh.readUpdate();
294 const fvMesh& mesh = vMesh.mesh();
296 INTEGER4 nFaceNodes = 0;
297 forAll(mesh.faces(), faceI)
299 nFaceNodes += mesh.faces()[faceI].size();
303 // Read all fields on the new mesh
304 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
306 // Search for list of objects for this time
307 IOobjectList objects(mesh, runTime.timeName());
309 HashSet<word> selectedFields;
310 if (args.optionFound("fields"))
312 args.optionLookup("fields")() >> selectedFields;
315 // Construct the vol fields (on the original mesh if subsetted)
317 PtrList<volScalarField> vsf;
318 readFields(vMesh, vMesh.baseMesh(), objects, selectedFields, vsf);
319 print(" volScalarFields :", Info, vsf);
321 PtrList<volVectorField> vvf;
322 readFields(vMesh, vMesh.baseMesh(), objects, selectedFields, vvf);
323 print(" volVectorFields :", Info, vvf);
325 PtrList<volSphericalTensorField> vSpheretf;
326 readFields(vMesh, vMesh.baseMesh(), objects, selectedFields, vSpheretf);
327 print(" volSphericalTensorFields :", Info, vSpheretf);
329 PtrList<volSymmTensorField> vSymmtf;
330 readFields(vMesh, vMesh.baseMesh(), objects, selectedFields, vSymmtf);
331 print(" volSymmTensorFields :", Info, vSymmtf);
333 PtrList<volTensorField> vtf;
334 readFields(vMesh, vMesh.baseMesh(), objects, selectedFields, vtf);
335 print(" volTensorFields :", Info, vtf);
338 // Construct pointMesh only if nessecary since constructs edge
339 // addressing (expensive on polyhedral meshes)
342 Info<< " pointScalarFields : switched off"
343 << " (\"-noPointValues\" option)\n";
344 Info<< " pointVectorFields : switched off"
345 << " (\"-noPointValues\" option)\n";
348 PtrList<pointScalarField> psf;
349 PtrList<pointVectorField> pvf;
350 //PtrList<pointSphericalTensorField> pSpheretf;
351 //PtrList<pointSymmTensorField> pSymmtf;
352 //PtrList<pointTensorField> ptf;
357 //// Add interpolated volFields
358 //const volPointInterpolation& pInterp = volPointInterpolation::New
363 //label nPsf = psf.size();
364 //psf.setSize(nPsf+vsf.size());
367 // Info<< "Interpolating " << vsf[i].name() << endl;
368 // tmp<pointScalarField> tvsf(pInterp.interpolate(vsf[i]));
369 // tvsf().rename(vsf[i].name() + "_point");
370 // psf.set(nPsf, tvsf);
374 //label nPvf = pvf.size();
375 //pvf.setSize(vvf.size());
378 // Info<< "Interpolating " << vvf[i].name() << endl;
379 // tmp<pointVectorField> tvvf(pInterp.interpolate(vvf[i]));
380 // tvvf().rename(vvf[i].name() + "_point");
381 // pvf.set(nPvf, tvvf);
388 pointMesh::New(vMesh.baseMesh()),
393 print(" pointScalarFields :", Info, psf);
398 pointMesh::New(vMesh.baseMesh()),
403 print(" pointVectorFields :", Info, pvf);
408 // pointMesh::New(vMesh.baseMesh()),
413 //print(" pointSphericalTensorFields :", Info, pSpheretf);
418 // pointMesh::New(vMesh.baseMesh()),
423 //print(" pointSymmTensorFields :", Info, pSymmtf);
428 // pointMesh::New(vMesh.baseMesh()),
433 //print(" pointTensorFields :", Info, ptf);
442 DynamicList<INTEGER4> varLocation;
445 DynamicList<INTEGER4> cellVarLocation;
448 tecplotWriter::getTecplotNames
451 ValueLocation_CellCentered,
455 tecplotWriter::getTecplotNames
458 ValueLocation_CellCentered,
463 tecplotWriter::getTecplotNames
466 ValueLocation_CellCentered,
470 tecplotWriter::getTecplotNames
473 ValueLocation_CellCentered,
478 tecplotWriter::getTecplotNames
481 ValueLocation_CellCentered,
485 tecplotWriter::getTecplotNames
488 ValueLocation_CellCentered,
493 tecplotWriter::getTecplotNames
496 ValueLocation_CellCentered,
500 tecplotWriter::getTecplotNames
503 ValueLocation_CellCentered,
508 tecplotWriter::getTecplotNames
511 ValueLocation_CellCentered,
515 tecplotWriter::getTecplotNames
518 ValueLocation_CellCentered,
526 tecplotWriter::getTecplotNames
534 tecplotWriter::getTecplotNames
542 // strandID (= piece id. Gets incremented for every piece of geometry
544 INTEGER4 strandID = 1;
547 if (meshState != polyMesh::UNCHANGED)
551 // Output mesh and fields
560 tecplotWriter writer(runTime);
562 string allVarNames = string("X Y Z ") + varNames;
563 DynamicList<INTEGER4> allVarLocation;
564 allVarLocation.append(ValueLocation_Nodal);
565 allVarLocation.append(ValueLocation_Nodal);
566 allVarLocation.append(ValueLocation_Nodal);
567 allVarLocation.append(varLocation);
577 writer.writePolyhedralZone
579 mesh.name(), // regionName
580 strandID++, // strandID
587 writer.writeField(mesh.points().component(0)());
588 writer.writeField(mesh.points().component(1)());
589 writer.writeField(mesh.points().component(2)());
594 writer.writeField(vsf[i]);
598 writer.writeField(vvf[i]);
602 writer.writeField(vSpheretf[i]);
606 writer.writeField(vSymmtf[i]);
610 writer.writeField(vtf[i]);
615 writer.writeField(psf[i]);
619 writer.writeField(pvf[i]);
622 writer.writeConnectivity(mesh);
632 // Output static mesh only
641 tecplotWriter writer(runTime);
651 writer.writePolyhedralZone
653 mesh.name(), // regionName
654 strandID, // strandID
656 List<INTEGER4>(3, ValueLocation_Nodal),
661 writer.writeField(mesh.points().component(0)());
662 writer.writeField(mesh.points().component(1)());
663 writer.writeField(mesh.points().component(2)());
665 writer.writeConnectivity(mesh);
669 // Output solution file
678 tecplotWriter writer(runTime);
685 DataFileType_Solution
688 writer.writePolyhedralZone
690 mesh.name(), // regionName
691 strandID++, // strandID
700 writer.writeField(vsf[i]);
704 writer.writeField(vvf[i]);
708 writer.writeField(vSpheretf[i]);
712 writer.writeField(vSymmtf[i]);
716 writer.writeField(vtf[i]);
721 writer.writeField(psf[i]);
725 writer.writeField(pvf[i]);
732 //---------------------------------------------------------------------
736 //---------------------------------------------------------------------
738 if (args.optionFound("faceSet"))
741 word setName(args.option("faceSet"));
742 labelList faceLabels(faceSet(mesh, setName).toc());
744 // Filename as if patch with same name.
745 mkDir(fvPath/setName);
747 fileName patchFileName
749 fvPath/setName/setName
755 Info<< " FaceSet : " << patchFileName << endl;
757 tecplotWriter writer(runTime);
759 string allVarNames = string("X Y Z ") + cellVarNames;
760 DynamicList<INTEGER4> allVarLocation;
761 allVarLocation.append(ValueLocation_Nodal);
762 allVarLocation.append(ValueLocation_Nodal);
763 allVarLocation.append(ValueLocation_Nodal);
764 allVarLocation.append(cellVarLocation);
774 const indirectPrimitivePatch ipp
776 IndirectList<face>(mesh.faces(), faceLabels),
780 writer.writePolygonalZone
789 writer.writeField(ipp.localPoints().component(0)());
790 writer.writeField(ipp.localPoints().component(1)());
791 writer.writeField(ipp.localPoints().component(2)());
793 // Write all volfields
800 fvc::interpolate(vsf[i])(),
811 fvc::interpolate(vvf[i])(),
822 fvc::interpolate(vSpheretf[i])(),
833 fvc::interpolate(vSymmtf[i])(),
844 fvc::interpolate(vtf[i])(),
849 writer.writeConnectivity(ipp);
856 //---------------------------------------------------------------------
858 // Write patches as multi-zone file
860 //---------------------------------------------------------------------
862 const polyBoundaryMesh& patches = mesh.boundaryMesh();
864 labelList patchIDs(getSelectedPatches(patches, excludePatches));
866 mkDir(fvPath/"boundaryMesh");
868 fileName patchFileName;
870 if (vMesh.useSubMesh())
873 fvPath/"boundaryMesh"/cellSetName
881 fvPath/"boundaryMesh"/"boundaryMesh"
887 Info<< " Combined patches : " << patchFileName << endl;
889 tecplotWriter writer(runTime);
891 string allVarNames = string("X Y Z ") + varNames;
892 DynamicList<INTEGER4> allVarLocation;
893 allVarLocation.append(ValueLocation_Nodal);
894 allVarLocation.append(ValueLocation_Nodal);
895 allVarLocation.append(ValueLocation_Nodal);
896 allVarLocation.append(varLocation);
908 label patchID = patchIDs[i];
909 const polyPatch& pp = patches[patchID];
910 //INTEGER4 strandID = 1 + i;
912 Info<< " Writing patch " << patchID << "\t" << pp.name()
913 << "\tstrand:" << strandID << nl << endl;
915 const indirectPrimitivePatch ipp
917 IndirectList<face>(pp, identity(pp.size())),
921 writer.writePolygonalZone
924 strandID++, //strandID,
930 writer.writeField(ipp.localPoints().component(0)());
931 writer.writeField(ipp.localPoints().component(1)());
932 writer.writeField(ipp.localPoints().component(2)());
1000 psf[i].boundaryField()[patchID].patchInternalField()()
1007 pvf[i].boundaryField()[patchID].patchInternalField()()
1011 writer.writeConnectivity(ipp);
1017 //---------------------------------------------------------------------
1019 // Write faceZones as multi-zone file
1021 //---------------------------------------------------------------------
1023 const faceZoneMesh& zones = mesh.faceZones();
1025 if (doFaceZones && zones.size() > 0)
1027 mkDir(fvPath/"faceZoneMesh");
1029 fileName patchFileName;
1031 if (vMesh.useSubMesh())
1034 fvPath/"faceZoneMesh"/cellSetName
1042 fvPath/"faceZoneMesh"/"faceZoneMesh"
1048 Info<< " FaceZone : " << patchFileName << endl;
1050 tecplotWriter writer(runTime);
1052 string allVarNames = string("X Y Z ") + cellVarNames;
1053 DynamicList<INTEGER4> allVarLocation;
1054 allVarLocation.append(ValueLocation_Nodal);
1055 allVarLocation.append(ValueLocation_Nodal);
1056 allVarLocation.append(ValueLocation_Nodal);
1057 allVarLocation.append(cellVarLocation);
1067 forAll(zones, zoneI)
1069 const faceZone& pp = zones[zoneI];
1071 const indirectPrimitivePatch ipp
1073 IndirectList<face>(mesh.faces(), pp),
1077 writer.writePolygonalZone
1080 strandID++, //1+patchIDs.size()+zoneI, //strandID,
1085 // Write coordinates
1086 writer.writeField(ipp.localPoints().component(0)());
1087 writer.writeField(ipp.localPoints().component(1)());
1088 writer.writeField(ipp.localPoints().component(2)());
1090 // Write all volfields
1097 fvc::interpolate(vsf[i])(),
1108 fvc::interpolate(vvf[i])(),
1113 forAll(vSpheretf, i)
1119 fvc::interpolate(vSpheretf[i])(),
1130 fvc::interpolate(vSymmtf[i])(),
1141 fvc::interpolate(vtf[i])(),
1147 writer.writeConnectivity(ipp);
1153 //---------------------------------------------------------------------
1155 // Write lagrangian data
1157 //---------------------------------------------------------------------
1159 fileNameList cloudDirs
1163 runTime.timePath()/regionPrefix/cloud::prefix,
1168 forAll(cloudDirs, cloudI)
1170 IOobjectList sprayObjs
1174 cloud::prefix/cloudDirs[cloudI]
1177 IOobject* positionsPtr = sprayObjs.lookup("positions");
1181 mkDir(fvPath/cloud::prefix/cloudDirs[cloudI]);
1183 fileName lagrFileName
1185 fvPath/cloud::prefix/cloudDirs[cloudI]/cloudDirs[cloudI]
1186 + "_" + timeDesc + ".plt"
1189 Info<< " Lagrangian: " << lagrFileName << endl;
1191 wordList labelNames(sprayObjs.names(labelIOField::typeName));
1193 print(Info, labelNames);
1195 wordList scalarNames(sprayObjs.names(scalarIOField::typeName));
1196 Info<< " scalars :";
1197 print(Info, scalarNames);
1199 wordList vectorNames(sprayObjs.names(vectorIOField::typeName));
1200 Info<< " vectors :";
1201 print(Info, vectorNames);
1203 //wordList sphereNames
1207 // sphericalTensorIOField::typeName
1210 //Info<< " spherical tensors :";
1211 //print(Info, sphereNames);
1213 //wordList symmNames
1217 // symmTensorIOField::typeName
1220 //Info<< " symm tensors :";
1221 //print(Info, symmNames);
1223 //wordList tensorNames(sprayObjs.names(tensorIOField::typeName));
1224 //Info<< " tensors :";
1225 //print(Info, tensorNames);
1228 // Load cloud positions
1229 passiveParticleCloud parcels(mesh, cloudDirs[cloudI]);
1231 // Get positions as pointField
1232 pointField positions(parcels.size());
1234 forAllConstIter(passiveParticleCloud, parcels, elmnt)
1236 positions[n++] = elmnt().position();
1240 string allVarNames = string("X Y Z");
1241 DynamicList<INTEGER4> allVarLocation;
1242 allVarLocation.append(ValueLocation_Nodal);
1243 allVarLocation.append(ValueLocation_Nodal);
1244 allVarLocation.append(ValueLocation_Nodal);
1246 tecplotWriter::getTecplotNames<label>
1249 ValueLocation_Nodal,
1254 tecplotWriter::getTecplotNames<scalar>
1257 ValueLocation_Nodal,
1262 tecplotWriter::getTecplotNames<vector>
1265 ValueLocation_Nodal,
1271 tecplotWriter writer(runTime);
1281 writer.writeOrderedZone
1284 strandID++, //strandID,
1289 // Write coordinates
1290 writer.writeField(positions.component(0)());
1291 writer.writeField(positions.component(1)());
1292 writer.writeField(positions.component(2)());
1295 forAll(labelNames, i)
1302 mesh.time().timeName(),
1303 cloud::prefix/cloudDirs[cloudI],
1305 IOobject::MUST_READ,
1311 scalarField sfld(fld.size());
1314 sfld[j] = scalar(fld[j]);
1316 writer.writeField(sfld);
1319 forAll(scalarNames, i)
1326 mesh.time().timeName(),
1327 cloud::prefix/cloudDirs[cloudI],
1329 IOobject::MUST_READ,
1334 writer.writeField(fld);
1337 forAll(vectorNames, i)
1344 mesh.time().timeName(),
1345 cloud::prefix/cloudDirs[cloudI],
1347 IOobject::MUST_READ,
1352 writer.writeField(fld);
1360 Info<< "End\n" << endl;
1366 // ************************************************************************* //