1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | foam-extend: Open Source CFD
4 \\ / O peration | Version: 3.2
5 \\ / A nd | Web: http://www.foam-extend.org
6 \\/ M anipulation | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
9 This file is part of foam-extend.
11 foam-extend 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 3 of the License, or (at your
14 option) any later version.
16 foam-extend is distributed in the hope that it will be useful, but
17 WITHOUT ANY WARRANTY; without even the implied warranty of
18 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 General Public License for more details.
21 You should have received a copy of the GNU General Public License
22 along with foam-extend. If not, see <http://www.gnu.org/licenses/>.
28 Tecplot binary file format writer.
32 - foamToTecplot360 [OPTION]
34 @param -fields \<fields\>\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)
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[])
166 timeSelector::addOptions();
168 # include "addRegionOption.H"
170 argList::validOptions.insert("fields", "fields");
171 argList::validOptions.insert("cellSet", "cellSet name");
172 argList::validOptions.insert("faceSet", "faceSet name");
173 argList::validOptions.insert("nearCellValue","");
174 argList::validOptions.insert("noInternal","");
175 argList::validOptions.insert("noPointValues","");
176 argList::validOptions.insert
179 "patches (wildcards) to exclude"
181 argList::validOptions.insert("noFaceZones","");
183 # include "setRootCase.H"
184 # include "createTime.H"
186 bool doWriteInternal = !args.optionFound("noInternal");
187 bool doFaceZones = !args.optionFound("noFaceZones");
189 bool nearCellValue = args.optionFound("nearCellValue");
193 WarningIn(args.executable())
194 << "Using neighbouring cell value instead of patch value"
198 bool noPointValues = args.optionFound("noPointValues");
202 WarningIn(args.executable())
203 << "Outputting cell values only" << nl << endl;
206 List<wordRe> excludePatches;
207 if (args.optionFound("excludePatches"))
209 args.optionLookup("excludePatches")() >> excludePatches;
211 Info<< "Not including patches " << excludePatches << nl << endl;
217 if (args.optionFound("cellSet"))
219 cellSetName = args.option("cellSet");
220 vtkName = cellSetName;
222 else if (Pstream::parRun())
224 // Strip off leading casename, leaving just processor_DDD ending.
225 vtkName = runTime.caseName();
227 string::size_type i = vtkName.rfind("processor");
229 if (i != string::npos)
231 vtkName = vtkName.substr(i);
236 vtkName = runTime.caseName();
240 instantList timeDirs = timeSelector::select0(runTime, args);
245 if (args.optionReadIfPresent("region", regionName))
247 Info<< "Create mesh " << regionName << " for time = "
248 << runTime.timeName() << Foam::nl << Foam::endl;
252 regionName = fvMesh::defaultRegion;
253 Info<< "Create mesh for time = "
254 << runTime.timeName() << Foam::nl << Foam::endl;
257 // TecplotData/ directory in the case
258 fileName fvPath(runTime.path()/"Tecplot360");
259 // Directory of mesh (region0 gets filtered out)
260 fileName regionPrefix = "";
262 if (regionName != polyMesh::defaultRegion)
264 fvPath = fvPath/regionName;
265 regionPrefix = regionName;
272 args.optionFound("time")
273 || args.optionFound("latestTime")
274 || cellSetName.size()
275 || regionName != polyMesh::defaultRegion
278 Info<< "Keeping old files in " << fvPath << nl << endl;
282 Info<< "Deleting old VTK files in " << fvPath << nl << endl;
291 // mesh wrapper; does subsetting and decomposition
299 Foam::IOobject::MUST_READ
304 forAll(timeDirs, timeI)
306 runTime.setTime(timeDirs[timeI], timeI);
308 Info<< "Time: " << runTime.timeName() << endl;
310 const word timeDesc = name(timeI); //name(runTime.timeIndex());
312 // Check for new polyMesh/ and update mesh, fvMeshSubset and cell
314 polyMesh::readUpdateState meshState = vMesh.readUpdate();
316 const fvMesh& mesh = vMesh.mesh();
318 INTEGER4 nFaceNodes = 0;
319 forAll(mesh.faces(), faceI)
321 nFaceNodes += mesh.faces()[faceI].size();
325 // Read all fields on the new mesh
326 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
328 // Search for list of objects for this time
329 IOobjectList objects(mesh, runTime.timeName());
331 HashSet<word> selectedFields;
332 if (args.optionFound("fields"))
334 args.optionLookup("fields")() >> selectedFields;
337 // Construct the vol fields (on the original mesh if subsetted)
339 PtrList<volScalarField> vsf;
340 readFields(vMesh, vMesh, objects, selectedFields, vsf);
341 print(" volScalarFields :", Info, vsf);
343 PtrList<volVectorField> vvf;
344 readFields(vMesh, vMesh, objects, selectedFields, vvf);
345 print(" volVectorFields :", Info, vvf);
347 PtrList<volSphericalTensorField> vSpheretf;
348 readFields(vMesh, vMesh, objects, selectedFields, vSpheretf);
349 print(" volSphericalTensorFields :", Info, vSpheretf);
351 PtrList<volSymmTensorField> vSymmtf;
352 readFields(vMesh, vMesh, objects, selectedFields, vSymmtf);
353 print(" volSymmTensorFields :", Info, vSymmtf);
355 PtrList<volTensorField> vtf;
356 readFields(vMesh, vMesh, objects, selectedFields, vtf);
357 print(" volTensorFields :", Info, vtf);
360 // Construct pointMesh only if nessecary since constructs edge
361 // addressing (expensive on polyhedral meshes)
364 Info<< " pointScalarFields : switched off"
365 << " (\"-noPointValues\" option)\n";
366 Info<< " pointVectorFields : switched off"
367 << " (\"-noPointValues\" option)\n";
370 PtrList<pointScalarField> psf;
371 PtrList<pointVectorField> pvf;
372 //PtrList<pointSphericalTensorField> pSpheretf;
373 //PtrList<pointSymmTensorField> pSymmtf;
374 //PtrList<pointTensorField> ptf;
379 //// Add interpolated volFields
380 //const volPointInterpolation& pInterp = volPointInterpolation::New
385 //label nPsf = psf.size();
386 //psf.setSize(nPsf+vsf.size());
389 // Info<< "Interpolating " << vsf[i].name() << endl;
390 // tmp<pointScalarField> tvsf(pInterp.interpolate(vsf[i]));
391 // tvsf().rename(vsf[i].name() + "_point");
392 // psf.set(nPsf, tvsf);
396 //label nPvf = pvf.size();
397 //pvf.setSize(vvf.size());
400 // Info<< "Interpolating " << vvf[i].name() << endl;
401 // tmp<pointVectorField> tvvf(pInterp.interpolate(vvf[i]));
402 // tvvf().rename(vvf[i].name() + "_point");
403 // pvf.set(nPvf, tvvf);
410 pointMesh::New(vMesh),
415 print(" pointScalarFields :", Info, psf);
420 pointMesh::New(vMesh),
425 print(" pointVectorFields :", Info, pvf);
430 // pointMesh::New(vMesh),
435 //print(" pointSphericalTensorFields :", Info, pSpheretf);
440 // pointMesh::New(vMesh),
445 //print(" pointSymmTensorFields :", Info, pSymmtf);
450 // pointMesh::New(vMesh),
455 //print(" pointTensorFields :", Info, ptf);
464 DynamicList<INTEGER4> varLocation;
467 DynamicList<INTEGER4> cellVarLocation;
470 tecplotWriter::getTecplotNames
473 ValueLocation_CellCentered,
477 tecplotWriter::getTecplotNames
480 ValueLocation_CellCentered,
485 tecplotWriter::getTecplotNames
488 ValueLocation_CellCentered,
492 tecplotWriter::getTecplotNames
495 ValueLocation_CellCentered,
500 tecplotWriter::getTecplotNames
503 ValueLocation_CellCentered,
507 tecplotWriter::getTecplotNames
510 ValueLocation_CellCentered,
515 tecplotWriter::getTecplotNames
518 ValueLocation_CellCentered,
522 tecplotWriter::getTecplotNames
525 ValueLocation_CellCentered,
530 tecplotWriter::getTecplotNames
533 ValueLocation_CellCentered,
537 tecplotWriter::getTecplotNames
540 ValueLocation_CellCentered,
548 tecplotWriter::getTecplotNames
556 tecplotWriter::getTecplotNames
564 // strandID (= piece id. Gets incremented for every piece of geometry
566 INTEGER4 strandID = 1;
569 if (meshState != polyMesh::UNCHANGED)
573 // Output mesh and fields
582 tecplotWriter writer(runTime);
584 string allVarNames = string("X Y Z ") + varNames;
585 DynamicList<INTEGER4> allVarLocation;
586 allVarLocation.append(ValueLocation_Nodal);
587 allVarLocation.append(ValueLocation_Nodal);
588 allVarLocation.append(ValueLocation_Nodal);
589 allVarLocation.append(varLocation);
599 writer.writePolyhedralZone
601 mesh.name(), // regionName
602 strandID++, // strandID
609 writer.writeField(mesh.points().component(0)());
610 writer.writeField(mesh.points().component(1)());
611 writer.writeField(mesh.points().component(2)());
616 writer.writeField(vsf[i]);
620 writer.writeField(vvf[i]);
624 writer.writeField(vSpheretf[i]);
628 writer.writeField(vSymmtf[i]);
632 writer.writeField(vtf[i]);
637 writer.writeField(psf[i]);
641 writer.writeField(pvf[i]);
644 writer.writeConnectivity(mesh);
654 // Output static mesh only
663 tecplotWriter writer(runTime);
673 writer.writePolyhedralZone
675 mesh.name(), // regionName
676 strandID, // strandID
678 List<INTEGER4>(3, ValueLocation_Nodal),
683 writer.writeField(mesh.points().component(0)());
684 writer.writeField(mesh.points().component(1)());
685 writer.writeField(mesh.points().component(2)());
687 writer.writeConnectivity(mesh);
691 // Output solution file
700 tecplotWriter writer(runTime);
707 DataFileType_Solution
710 writer.writePolyhedralZone
712 mesh.name(), // regionName
713 strandID++, // strandID
722 writer.writeField(vsf[i]);
726 writer.writeField(vvf[i]);
730 writer.writeField(vSpheretf[i]);
734 writer.writeField(vSymmtf[i]);
738 writer.writeField(vtf[i]);
743 writer.writeField(psf[i]);
747 writer.writeField(pvf[i]);
754 //---------------------------------------------------------------------
758 //---------------------------------------------------------------------
760 if (args.optionFound("faceSet"))
763 word setName(args.option("faceSet"));
764 labelList faceLabels(faceSet(mesh, setName).toc());
766 // Filename as if patch with same name.
767 mkDir(fvPath/setName);
769 fileName patchFileName
771 fvPath/setName/setName
777 Info<< " FaceSet : " << patchFileName << endl;
779 tecplotWriter writer(runTime);
781 string allVarNames = string("X Y Z ") + cellVarNames;
782 DynamicList<INTEGER4> allVarLocation;
783 allVarLocation.append(ValueLocation_Nodal);
784 allVarLocation.append(ValueLocation_Nodal);
785 allVarLocation.append(ValueLocation_Nodal);
786 allVarLocation.append(cellVarLocation);
796 const indirectPrimitivePatch ipp
798 IndirectList<face>(mesh.faces(), faceLabels),
802 writer.writePolygonalZone
811 writer.writeField(ipp.localPoints().component(0)());
812 writer.writeField(ipp.localPoints().component(1)());
813 writer.writeField(ipp.localPoints().component(2)());
815 // Write all volfields
822 linearInterpolate(vsf[i])(),
833 linearInterpolate(vvf[i])(),
844 linearInterpolate(vSpheretf[i])(),
855 linearInterpolate(vSymmtf[i])(),
866 linearInterpolate(vtf[i])(),
871 writer.writeConnectivity(ipp);
878 //---------------------------------------------------------------------
880 // Write patches as multi-zone file
882 //---------------------------------------------------------------------
884 const polyBoundaryMesh& patches = mesh.boundaryMesh();
886 labelList patchIDs(getSelectedPatches(patches, excludePatches));
888 mkDir(fvPath/"boundaryMesh");
890 fileName patchFileName;
892 if (vMesh.useSubMesh())
895 fvPath/"boundaryMesh"/cellSetName
903 fvPath/"boundaryMesh"/"boundaryMesh"
909 Info<< " Combined patches : " << patchFileName << endl;
911 tecplotWriter writer(runTime);
913 string allVarNames = string("X Y Z ") + varNames;
914 DynamicList<INTEGER4> allVarLocation;
915 allVarLocation.append(ValueLocation_Nodal);
916 allVarLocation.append(ValueLocation_Nodal);
917 allVarLocation.append(ValueLocation_Nodal);
918 allVarLocation.append(varLocation);
930 label patchID = patchIDs[i];
931 const polyPatch& pp = patches[patchID];
932 //INTEGER4 strandID = 1 + i;
936 Info<< " Writing patch " << patchID << "\t" << pp.name()
937 << "\tstrand:" << strandID << nl << endl;
939 const indirectPrimitivePatch ipp
941 IndirectList<face>(pp, identity(pp.size())),
945 writer.writePolygonalZone
948 strandID++, //strandID,
954 writer.writeField(ipp.localPoints().component(0)());
955 writer.writeField(ipp.localPoints().component(1)());
956 writer.writeField(ipp.localPoints().component(2)());
1011 writer.getPatchField
1024 psf[i].boundaryField()[patchID].patchInternalField()()
1031 pvf[i].boundaryField()[patchID].patchInternalField()()
1035 writer.writeConnectivity(ipp);
1039 Info<< " Skipping zero sized patch " << patchID
1040 << "\t" << pp.name()
1049 //---------------------------------------------------------------------
1051 // Write faceZones as multi-zone file
1053 //---------------------------------------------------------------------
1055 const faceZoneMesh& zones = mesh.faceZones();
1057 if (doFaceZones && zones.size() > 0)
1059 mkDir(fvPath/"faceZoneMesh");
1061 fileName patchFileName;
1063 if (vMesh.useSubMesh())
1066 fvPath/"faceZoneMesh"/cellSetName
1074 fvPath/"faceZoneMesh"/"faceZoneMesh"
1080 Info<< " FaceZone : " << patchFileName << endl;
1082 tecplotWriter writer(runTime);
1084 string allVarNames = string("X Y Z ") + cellVarNames;
1085 DynamicList<INTEGER4> allVarLocation;
1086 allVarLocation.append(ValueLocation_Nodal);
1087 allVarLocation.append(ValueLocation_Nodal);
1088 allVarLocation.append(ValueLocation_Nodal);
1089 allVarLocation.append(cellVarLocation);
1099 forAll(zones, zoneI)
1101 const faceZone& pp = zones[zoneI];
1105 const indirectPrimitivePatch ipp
1107 IndirectList<face>(mesh.faces(), pp),
1111 writer.writePolygonalZone
1114 strandID++, //1+patchIDs.size()+zoneI, //strandID,
1119 // Write coordinates
1120 writer.writeField(ipp.localPoints().component(0)());
1121 writer.writeField(ipp.localPoints().component(1)());
1122 writer.writeField(ipp.localPoints().component(2)());
1124 // Write all volfields
1131 linearInterpolate(vsf[i])(),
1142 linearInterpolate(vvf[i])(),
1147 forAll(vSpheretf, i)
1153 linearInterpolate(vSpheretf[i])(),
1164 linearInterpolate(vSymmtf[i])(),
1175 linearInterpolate(vtf[i])(),
1181 writer.writeConnectivity(ipp);
1185 Info<< " Skipping zero sized faceZone " << zoneI
1186 << "\t" << pp.name()
1196 //---------------------------------------------------------------------
1198 // Write lagrangian data
1200 //---------------------------------------------------------------------
1202 fileNameList cloudDirs
1206 runTime.timePath()/regionPrefix/cloud::prefix,
1211 forAll(cloudDirs, cloudI)
1213 IOobjectList sprayObjs
1217 cloud::prefix/cloudDirs[cloudI]
1220 IOobject* positionsPtr = sprayObjs.lookup("positions");
1224 mkDir(fvPath/cloud::prefix/cloudDirs[cloudI]);
1226 fileName lagrFileName
1228 fvPath/cloud::prefix/cloudDirs[cloudI]/cloudDirs[cloudI]
1229 + "_" + timeDesc + ".plt"
1232 Info<< " Lagrangian: " << lagrFileName << endl;
1234 wordList labelNames(sprayObjs.names(labelIOField::typeName));
1236 print(Info, labelNames);
1238 wordList scalarNames(sprayObjs.names(scalarIOField::typeName));
1239 Info<< " scalars :";
1240 print(Info, scalarNames);
1242 wordList vectorNames(sprayObjs.names(vectorIOField::typeName));
1243 Info<< " vectors :";
1244 print(Info, vectorNames);
1246 //wordList sphereNames
1250 // sphericalTensorIOField::typeName
1253 //Info<< " spherical tensors :";
1254 //print(Info, sphereNames);
1256 //wordList symmNames
1260 // symmTensorIOField::typeName
1263 //Info<< " symm tensors :";
1264 //print(Info, symmNames);
1266 //wordList tensorNames(sprayObjs.names(tensorIOField::typeName));
1267 //Info<< " tensors :";
1268 //print(Info, tensorNames);
1271 // Load cloud positions
1272 passiveParticleCloud parcels(mesh, cloudDirs[cloudI]);
1274 // Get positions as pointField
1275 pointField positions(parcels.size());
1277 forAllConstIter(passiveParticleCloud, parcels, elmnt)
1279 positions[n++] = elmnt().position();
1283 string allVarNames = string("X Y Z");
1284 DynamicList<INTEGER4> allVarLocation;
1285 allVarLocation.append(ValueLocation_Nodal);
1286 allVarLocation.append(ValueLocation_Nodal);
1287 allVarLocation.append(ValueLocation_Nodal);
1289 tecplotWriter::getTecplotNames<label>
1292 ValueLocation_Nodal,
1297 tecplotWriter::getTecplotNames<scalar>
1300 ValueLocation_Nodal,
1305 tecplotWriter::getTecplotNames<vector>
1308 ValueLocation_Nodal,
1314 tecplotWriter writer(runTime);
1324 writer.writeOrderedZone
1327 strandID++, //strandID,
1332 // Write coordinates
1333 writer.writeField(positions.component(0)());
1334 writer.writeField(positions.component(1)());
1335 writer.writeField(positions.component(2)());
1338 forAll(labelNames, i)
1345 mesh.time().timeName(),
1346 cloud::prefix/cloudDirs[cloudI],
1348 IOobject::MUST_READ,
1354 scalarField sfld(fld.size());
1357 sfld[j] = scalar(fld[j]);
1359 writer.writeField(sfld);
1362 forAll(scalarNames, i)
1369 mesh.time().timeName(),
1370 cloud::prefix/cloudDirs[cloudI],
1372 IOobject::MUST_READ,
1377 writer.writeField(fld);
1380 forAll(vectorNames, i)
1387 mesh.time().timeName(),
1388 cloud::prefix/cloudDirs[cloudI],
1390 IOobject::MUST_READ,
1395 writer.writeField(fld);
1403 Info<< "End\n" << endl;
1409 // ************************************************************************* //