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 Tecplot binary file format writer.
32 - foamToTecplot360 [OPTION]
34 \param -fields \<names\>\n
35 Convert selected fields only. For example,
39 The quoting is required to avoid shell expansions and to pass the
40 information as a single argument.
42 \param -cellSet \<name\>\n
43 \param -faceSet \<name\>\n
44 Restrict conversion to the cellSet, faceSet.
46 \param -nearCellValue \n
47 Output cell value on patches instead of patch value itself
50 Do not generate file for mesh, only for patches
52 \param -noPointValues \n
55 \param -noFaceZones \n
58 \param -excludePatches \<patchNames\>\n
59 Specify patches (wildcards) to exclude. For example,
61 -excludePatches '( inlet_1 inlet_2 "proc.*")'
63 The quoting is required to avoid shell expansions and to pass the
64 information as a single argument. The double quotes denote a regular
67 \param -useTimeName \n
68 use the time index in the VTK file name instead of the time index
70 \*---------------------------------------------------------------------------*/
72 #include "pointMesh.H"
73 #include "volPointInterpolation.H"
74 #include "emptyPolyPatch.H"
75 #include "labelIOField.H"
76 #include "scalarIOField.H"
77 #include "sphericalTensorIOField.H"
78 #include "symmTensorIOField.H"
79 #include "tensorIOField.H"
80 #include "passiveParticleCloud.H"
82 #include "stringListOps.H"
86 #include "readFields.H"
87 #include "tecplotWriter.H"
91 // Note: needs to be after TECIO to prevent Foam::Time conflicting with
95 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
97 template<class GeoField>
98 void print(const char* msg, Ostream& os, const PtrList<GeoField>& flds)
105 os << ' ' << flds[i].name();
112 void print(Ostream& os, const wordList& flds)
116 os << ' ' << flds[i];
122 labelList getSelectedPatches
124 const polyBoundaryMesh& patches,
125 const List<wordRe>& excludePatches //HashSet<word>& excludePatches
128 DynamicList<label> patchIDs(patches.size());
130 Info<< "Combining patches:" << endl;
132 forAll(patches, patchI)
134 const polyPatch& pp = patches[patchI];
138 isType<emptyPolyPatch>(pp)
139 || (Pstream::parRun() && isType<processorPolyPatch>(pp))
142 Info<< " discarding empty/processor patch " << patchI
143 << " " << pp.name() << endl;
145 else if (findStrings(excludePatches, pp.name()))
147 Info<< " excluding patch " << patchI
148 << " " << pp.name() << endl;
152 patchIDs.append(patchI);
153 Info<< " patch " << patchI << " " << pp.name() << endl;
156 return patchIDs.shrink();
164 int main(int argc, char *argv[])
168 "Tecplot binary file format writer"
171 timeSelector::addOptions();
172 #include "addRegionOption.H"
178 "convert selected fields only. eg, '(p T U)'"
184 "restrict conversion to the specified cellSet"
190 "restrict conversion to the specified cellSet"
192 argList::addBoolOption
195 "output cell value on patches instead of patch value itself"
197 argList::addBoolOption
200 "do not generate file for mesh, only for patches"
202 argList::addBoolOption
210 "patches (wildcards) to exclude"
212 argList::addBoolOption
218 #include "setRootCase.H"
219 #include "createTime.H"
221 const bool doWriteInternal = !args.optionFound("noInternal");
222 const bool doFaceZones = !args.optionFound("noFaceZones");
223 const bool nearCellValue = args.optionFound("nearCellValue");
224 const bool noPointValues = args.optionFound("noPointValues");
228 WarningIn(args.executable())
229 << "Using neighbouring cell value instead of patch value"
235 WarningIn(args.executable())
236 << "Outputting cell values only" << nl << endl;
239 List<wordRe> excludePatches;
240 if (args.optionFound("excludePatches"))
242 args.optionLookup("excludePatches")() >> excludePatches;
244 Info<< "Not including patches " << excludePatches << nl << endl;
250 if (args.optionReadIfPresent("cellSet", cellSetName))
252 vtkName = cellSetName;
254 else if (Pstream::parRun())
256 // Strip off leading casename, leaving just processor_DDD ending.
257 vtkName = runTime.caseName();
259 string::size_type i = vtkName.rfind("processor");
261 if (i != string::npos)
263 vtkName = vtkName.substr(i);
268 vtkName = runTime.caseName();
272 instantList timeDirs = timeSelector::select0(runTime, args);
274 # include "createNamedMesh.H"
276 // TecplotData/ directory in the case
277 fileName fvPath(runTime.path()/"Tecplot360");
278 // Directory of mesh (region0 gets filtered out)
279 fileName regionPrefix = "";
281 if (regionName != polyMesh::defaultRegion)
283 fvPath = fvPath/regionName;
284 regionPrefix = regionName;
291 args.optionFound("time")
292 || args.optionFound("latestTime")
293 || cellSetName.size()
294 || regionName != polyMesh::defaultRegion
297 Info<< "Keeping old files in " << fvPath << nl << endl;
301 Info<< "Deleting old VTK files in " << fvPath << nl << endl;
310 // mesh wrapper; does subsetting and decomposition
311 vtkMesh vMesh(mesh, cellSetName);
313 forAll(timeDirs, timeI)
315 runTime.setTime(timeDirs[timeI], timeI);
317 Info<< "Time: " << runTime.timeName() << endl;
319 const word timeDesc = name(timeI); //name(runTime.timeIndex());
321 // Check for new polyMesh/ and update mesh, fvMeshSubset and cell
323 polyMesh::readUpdateState meshState = vMesh.readUpdate();
325 const fvMesh& mesh = vMesh.mesh();
327 INTEGER4 nFaceNodes = 0;
328 forAll(mesh.faces(), faceI)
330 nFaceNodes += mesh.faces()[faceI].size();
334 // Read all fields on the new mesh
335 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
337 // Search for list of objects for this time
338 IOobjectList objects(mesh, runTime.timeName());
340 HashSet<word> selectedFields;
341 if (args.optionFound("fields"))
343 args.optionLookup("fields")() >> selectedFields;
346 // Construct the vol fields (on the original mesh if subsetted)
348 PtrList<volScalarField> vsf;
349 readFields(vMesh, vMesh.baseMesh(), objects, selectedFields, vsf);
350 print(" volScalarFields :", Info, vsf);
352 PtrList<volVectorField> vvf;
353 readFields(vMesh, vMesh.baseMesh(), objects, selectedFields, vvf);
354 print(" volVectorFields :", Info, vvf);
356 PtrList<volSphericalTensorField> vSpheretf;
357 readFields(vMesh, vMesh.baseMesh(), objects, selectedFields, vSpheretf);
358 print(" volSphericalTensorFields :", Info, vSpheretf);
360 PtrList<volSymmTensorField> vSymmtf;
361 readFields(vMesh, vMesh.baseMesh(), objects, selectedFields, vSymmtf);
362 print(" volSymmTensorFields :", Info, vSymmtf);
364 PtrList<volTensorField> vtf;
365 readFields(vMesh, vMesh.baseMesh(), objects, selectedFields, vtf);
366 print(" volTensorFields :", Info, vtf);
369 // Construct pointMesh only if nessecary since constructs edge
370 // addressing (expensive on polyhedral meshes)
373 Info<< " pointScalarFields : switched off"
374 << " (\"-noPointValues\" (at your option)\n";
375 Info<< " pointVectorFields : switched off"
376 << " (\"-noPointValues\" (at your option)\n";
379 PtrList<pointScalarField> psf;
380 PtrList<pointVectorField> pvf;
381 //PtrList<pointSphericalTensorField> pSpheretf;
382 //PtrList<pointSymmTensorField> pSymmtf;
383 //PtrList<pointTensorField> ptf;
388 //// Add interpolated volFields
389 //const volPointInterpolation& pInterp = volPointInterpolation::New
394 //label nPsf = psf.size();
395 //psf.setSize(nPsf+vsf.size());
398 // Info<< "Interpolating " << vsf[i].name() << endl;
399 // tmp<pointScalarField> tvsf(pInterp.interpolate(vsf[i]));
400 // tvsf().rename(vsf[i].name() + "_point");
401 // psf.set(nPsf, tvsf);
405 //label nPvf = pvf.size();
406 //pvf.setSize(vvf.size());
409 // Info<< "Interpolating " << vvf[i].name() << endl;
410 // tmp<pointVectorField> tvvf(pInterp.interpolate(vvf[i]));
411 // tvvf().rename(vvf[i].name() + "_point");
412 // pvf.set(nPvf, tvvf);
419 pointMesh::New(vMesh.baseMesh()),
424 print(" pointScalarFields :", Info, psf);
429 pointMesh::New(vMesh.baseMesh()),
434 print(" pointVectorFields :", Info, pvf);
439 // pointMesh::New(vMesh.baseMesh()),
444 //print(" pointSphericalTensorFields :", Info, pSpheretf);
449 // pointMesh::New(vMesh.baseMesh()),
454 //print(" pointSymmTensorFields :", Info, pSymmtf);
459 // pointMesh::New(vMesh.baseMesh()),
464 //print(" pointTensorFields :", Info, ptf);
473 DynamicList<INTEGER4> varLocation;
476 DynamicList<INTEGER4> cellVarLocation;
479 tecplotWriter::getTecplotNames
482 ValueLocation_CellCentered,
486 tecplotWriter::getTecplotNames
489 ValueLocation_CellCentered,
494 tecplotWriter::getTecplotNames
497 ValueLocation_CellCentered,
501 tecplotWriter::getTecplotNames
504 ValueLocation_CellCentered,
509 tecplotWriter::getTecplotNames
512 ValueLocation_CellCentered,
516 tecplotWriter::getTecplotNames
519 ValueLocation_CellCentered,
524 tecplotWriter::getTecplotNames
527 ValueLocation_CellCentered,
531 tecplotWriter::getTecplotNames
534 ValueLocation_CellCentered,
539 tecplotWriter::getTecplotNames
542 ValueLocation_CellCentered,
546 tecplotWriter::getTecplotNames
549 ValueLocation_CellCentered,
557 tecplotWriter::getTecplotNames
565 tecplotWriter::getTecplotNames
573 // strandID (= piece id. Gets incremented for every piece of geometry
575 INTEGER4 strandID = 1;
578 if (meshState != polyMesh::UNCHANGED)
582 // Output mesh and fields
591 tecplotWriter writer(runTime);
593 string allVarNames = string("X Y Z ") + varNames;
594 DynamicList<INTEGER4> allVarLocation;
595 allVarLocation.append(ValueLocation_Nodal);
596 allVarLocation.append(ValueLocation_Nodal);
597 allVarLocation.append(ValueLocation_Nodal);
598 allVarLocation.append(varLocation);
608 writer.writePolyhedralZone
610 mesh.name(), // regionName
611 strandID++, // strandID
618 writer.writeField(mesh.points().component(0)());
619 writer.writeField(mesh.points().component(1)());
620 writer.writeField(mesh.points().component(2)());
625 writer.writeField(vsf[i]);
629 writer.writeField(vvf[i]);
633 writer.writeField(vSpheretf[i]);
637 writer.writeField(vSymmtf[i]);
641 writer.writeField(vtf[i]);
646 writer.writeField(psf[i]);
650 writer.writeField(pvf[i]);
653 writer.writeConnectivity(mesh);
663 // Output static mesh only
672 tecplotWriter writer(runTime);
682 writer.writePolyhedralZone
684 mesh.name(), // regionName
685 strandID, // strandID
687 List<INTEGER4>(3, ValueLocation_Nodal),
692 writer.writeField(mesh.points().component(0)());
693 writer.writeField(mesh.points().component(1)());
694 writer.writeField(mesh.points().component(2)());
696 writer.writeConnectivity(mesh);
700 // Output solution file
709 tecplotWriter writer(runTime);
716 DataFileType_Solution
719 writer.writePolyhedralZone
721 mesh.name(), // regionName
722 strandID++, // strandID
731 writer.writeField(vsf[i]);
735 writer.writeField(vvf[i]);
739 writer.writeField(vSpheretf[i]);
743 writer.writeField(vSymmtf[i]);
747 writer.writeField(vtf[i]);
752 writer.writeField(psf[i]);
756 writer.writeField(pvf[i]);
763 //---------------------------------------------------------------------
767 //---------------------------------------------------------------------
769 if (args.optionFound("faceSet"))
772 const word setName = args["faceSet"];
773 labelList faceLabels(faceSet(mesh, setName).toc());
775 // Filename as if patch with same name.
776 mkDir(fvPath/setName);
778 fileName patchFileName
780 fvPath/setName/setName
786 Info<< " FaceSet : " << patchFileName << endl;
788 tecplotWriter writer(runTime);
790 string allVarNames = string("X Y Z ") + cellVarNames;
791 DynamicList<INTEGER4> allVarLocation;
792 allVarLocation.append(ValueLocation_Nodal);
793 allVarLocation.append(ValueLocation_Nodal);
794 allVarLocation.append(ValueLocation_Nodal);
795 allVarLocation.append(cellVarLocation);
805 const indirectPrimitivePatch ipp
807 IndirectList<face>(mesh.faces(), faceLabels),
811 writer.writePolygonalZone
820 writer.writeField(ipp.localPoints().component(0)());
821 writer.writeField(ipp.localPoints().component(1)());
822 writer.writeField(ipp.localPoints().component(2)());
824 // Write all volfields
831 linearInterpolate(vsf[i])(),
842 linearInterpolate(vvf[i])(),
853 linearInterpolate(vSpheretf[i])(),
864 linearInterpolate(vSymmtf[i])(),
875 linearInterpolate(vtf[i])(),
880 writer.writeConnectivity(ipp);
887 //---------------------------------------------------------------------
889 // Write patches as multi-zone file
891 //---------------------------------------------------------------------
893 const polyBoundaryMesh& patches = mesh.boundaryMesh();
895 labelList patchIDs(getSelectedPatches(patches, excludePatches));
897 mkDir(fvPath/"boundaryMesh");
899 fileName patchFileName;
901 if (vMesh.useSubMesh())
904 fvPath/"boundaryMesh"/cellSetName
912 fvPath/"boundaryMesh"/"boundaryMesh"
918 Info<< " Combined patches : " << patchFileName << endl;
920 tecplotWriter writer(runTime);
922 string allVarNames = string("X Y Z ") + varNames;
923 DynamicList<INTEGER4> allVarLocation;
924 allVarLocation.append(ValueLocation_Nodal);
925 allVarLocation.append(ValueLocation_Nodal);
926 allVarLocation.append(ValueLocation_Nodal);
927 allVarLocation.append(varLocation);
939 label patchID = patchIDs[i];
940 const polyPatch& pp = patches[patchID];
941 //INTEGER4 strandID = 1 + i;
945 Info<< " Writing patch " << patchID << "\t" << pp.name()
946 << "\tstrand:" << strandID << nl << endl;
948 const indirectPrimitivePatch ipp
950 IndirectList<face>(pp, identity(pp.size())),
954 writer.writePolygonalZone
957 strandID++, //strandID,
963 writer.writeField(ipp.localPoints().component(0)());
964 writer.writeField(ipp.localPoints().component(1)());
965 writer.writeField(ipp.localPoints().component(2)());
1008 writer.getPatchField
1020 writer.getPatchField
1033 psf[i].boundaryField()[patchID].patchInternalField()()
1040 pvf[i].boundaryField()[patchID].patchInternalField()()
1044 writer.writeConnectivity(ipp);
1048 Info<< " Skipping zero sized patch " << patchID
1049 << "\t" << pp.name()
1058 //---------------------------------------------------------------------
1060 // Write faceZones as multi-zone file
1062 //---------------------------------------------------------------------
1064 const faceZoneMesh& zones = mesh.faceZones();
1066 if (doFaceZones && zones.size() > 0)
1068 mkDir(fvPath/"faceZoneMesh");
1070 fileName patchFileName;
1072 if (vMesh.useSubMesh())
1075 fvPath/"faceZoneMesh"/cellSetName
1083 fvPath/"faceZoneMesh"/"faceZoneMesh"
1089 Info<< " FaceZone : " << patchFileName << endl;
1091 tecplotWriter writer(runTime);
1093 string allVarNames = string("X Y Z ") + cellVarNames;
1094 DynamicList<INTEGER4> allVarLocation;
1095 allVarLocation.append(ValueLocation_Nodal);
1096 allVarLocation.append(ValueLocation_Nodal);
1097 allVarLocation.append(ValueLocation_Nodal);
1098 allVarLocation.append(cellVarLocation);
1108 forAll(zones, zoneI)
1110 const faceZone& pp = zones[zoneI];
1114 const indirectPrimitivePatch ipp
1116 IndirectList<face>(mesh.faces(), pp),
1120 writer.writePolygonalZone
1123 strandID++, //1+patchIDs.size()+zoneI, //strandID,
1128 // Write coordinates
1129 writer.writeField(ipp.localPoints().component(0)());
1130 writer.writeField(ipp.localPoints().component(1)());
1131 writer.writeField(ipp.localPoints().component(2)());
1133 // Write all volfields
1140 linearInterpolate(vsf[i])(),
1151 linearInterpolate(vvf[i])(),
1156 forAll(vSpheretf, i)
1162 linearInterpolate(vSpheretf[i])(),
1173 linearInterpolate(vSymmtf[i])(),
1184 linearInterpolate(vtf[i])(),
1190 writer.writeConnectivity(ipp);
1194 Info<< " Skipping zero sized faceZone " << zoneI
1195 << "\t" << pp.name()
1205 //---------------------------------------------------------------------
1207 // Write lagrangian data
1209 //---------------------------------------------------------------------
1211 fileNameList cloudDirs
1215 runTime.timePath()/regionPrefix/cloud::prefix,
1220 forAll(cloudDirs, cloudI)
1222 IOobjectList sprayObjs
1226 cloud::prefix/cloudDirs[cloudI]
1229 IOobject* positionsPtr = sprayObjs.lookup("positions");
1233 mkDir(fvPath/cloud::prefix/cloudDirs[cloudI]);
1235 fileName lagrFileName
1237 fvPath/cloud::prefix/cloudDirs[cloudI]/cloudDirs[cloudI]
1238 + "_" + timeDesc + ".plt"
1241 Info<< " Lagrangian: " << lagrFileName << endl;
1243 wordList labelNames(sprayObjs.names(labelIOField::typeName));
1245 print(Info, labelNames);
1247 wordList scalarNames(sprayObjs.names(scalarIOField::typeName));
1248 Info<< " scalars :";
1249 print(Info, scalarNames);
1251 wordList vectorNames(sprayObjs.names(vectorIOField::typeName));
1252 Info<< " vectors :";
1253 print(Info, vectorNames);
1255 //wordList sphereNames
1259 // sphericalTensorIOField::typeName
1262 //Info<< " spherical tensors :";
1263 //print(Info, sphereNames);
1265 //wordList symmNames
1269 // symmTensorIOField::typeName
1272 //Info<< " symm tensors :";
1273 //print(Info, symmNames);
1275 //wordList tensorNames
1277 // sprayObjs.names(tensorIOField::typeName)
1279 //Info<< " tensors :";
1280 //print(Info, tensorNames);
1283 // Load cloud positions
1284 passiveParticleCloud parcels(mesh, cloudDirs[cloudI]);
1286 // Get positions as pointField
1287 pointField positions(parcels.size());
1289 forAllConstIter(passiveParticleCloud, parcels, elmnt)
1291 positions[n++] = elmnt().position();
1295 string allVarNames = string("X Y Z");
1296 DynamicList<INTEGER4> allVarLocation;
1297 allVarLocation.append(ValueLocation_Nodal);
1298 allVarLocation.append(ValueLocation_Nodal);
1299 allVarLocation.append(ValueLocation_Nodal);
1301 tecplotWriter::getTecplotNames<label>
1304 ValueLocation_Nodal,
1309 tecplotWriter::getTecplotNames<scalar>
1312 ValueLocation_Nodal,
1317 tecplotWriter::getTecplotNames<vector>
1320 ValueLocation_Nodal,
1326 tecplotWriter writer(runTime);
1336 writer.writeOrderedZone
1339 strandID++, //strandID,
1344 // Write coordinates
1345 writer.writeField(positions.component(0)());
1346 writer.writeField(positions.component(1)());
1347 writer.writeField(positions.component(2)());
1350 forAll(labelNames, i)
1357 mesh.time().timeName(),
1358 cloud::prefix/cloudDirs[cloudI],
1360 IOobject::MUST_READ,
1366 scalarField sfld(fld.size());
1369 sfld[j] = scalar(fld[j]);
1371 writer.writeField(sfld);
1374 forAll(scalarNames, i)
1381 mesh.time().timeName(),
1382 cloud::prefix/cloudDirs[cloudI],
1384 IOobject::MUST_READ,
1389 writer.writeField(fld);
1392 forAll(vectorNames, i)
1399 mesh.time().timeName(),
1400 cloud::prefix/cloudDirs[cloudI],
1402 IOobject::MUST_READ,
1407 writer.writeField(fld);
1415 Info<< "End\n" << endl;
1421 // ************************************************************************* //