Forward compatibility: flex
[foam-extend-3.2.git] / applications / utilities / postProcessing / dataConversion / foamToTecplot360 / foamToTecplot360.C
blobea26dfa5ec832bd2644251651224da3d45d7379c
1 /*---------------------------------------------------------------------------*\
2   =========                 |
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 -------------------------------------------------------------------------------
8 License
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/>.
24 Application
25     foamToTecplot360
27 Description
28     Tecplot binary file format writer.
30 Usage
32     - foamToTecplot360 [OPTION]
34     @param -fields \<fields\>\n
35     Convert selected fields only. For example,
36     @verbatim
37          -fields '( p T U )'
38     @endverbatim
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
49     @param -noInternal \n
50     Do not generate file for mesh, only for patches
52     @param -noPointValues \n
53     No pointFields
55     @param -noFaceZones \n
56     No faceZones
58     @param -excludePatches \<patchNames\>\n
59     Specify patches (wildcards) to exclude. For example,
60     @verbatim
61          -excludePatches '( inlet_1 inlet_2 "proc.*")'
62     @endverbatim
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
65     expression.
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"
81 #include "faceSet.H"
82 #include "stringListOps.H"
83 #include "wordRe.H"
85 #include "vtkMesh.H"
86 #include "readFields.H"
87 #include "tecplotWriter.H"
89 #include "TECIO.h"
91 // Note: needs to be after TECIO to prevent Foam::Time conflicting with
92 // Xlib Time.
93 #include "fvCFD.H"
95 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
97 template<class GeoField>
98 void print(const char* msg, Ostream& os, const PtrList<GeoField>& flds)
100     if (flds.size())
101     {
102         os << msg;
103         forAll(flds, i)
104         {
105             os<< ' ' << flds[i].name();
106         }
107         os << endl;
108     }
112 void print(Ostream& os, const wordList& flds)
114     forAll(flds, i)
115     {
116         os<< ' ' << flds[i];
117     }
118     os << endl;
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)
133     {
134         const polyPatch& pp = patches[patchI];
136         if
137         (
138             isType<emptyPolyPatch>(pp)
139             || (Pstream::parRun() && isType<processorPolyPatch>(pp))
140         )
141         {
142             Info<< "    discarding empty/processor patch " << patchI
143                 << " " << pp.name() << endl;
144         }
145         else if (findStrings(excludePatches, pp.name()))
146         {
147             Info<< "    excluding patch " << patchI
148                 << " " << pp.name() << endl;
149         }
150         else
151         {
152             patchIDs.append(patchI);
153             Info<< "    patch " << patchI << " " << pp.name() << endl;
154         }
155     }
156     return patchIDs.shrink();
162 // Main program:
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
177     (
178         "excludePatches",
179         "patches (wildcards) to exclude"
180     );
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");
191     if (nearCellValue)
192     {
193         WarningIn(args.executable())
194             << "Using neighbouring cell value instead of patch value"
195             << nl << endl;
196     }
198     bool noPointValues = args.optionFound("noPointValues");
200     if (noPointValues)
201     {
202         WarningIn(args.executable())
203             << "Outputting cell values only" << nl << endl;
204     }
206     List<wordRe> excludePatches;
207     if (args.optionFound("excludePatches"))
208     {
209         args.optionLookup("excludePatches")() >> excludePatches;
211         Info<< "Not including patches " << excludePatches << nl << endl;
212     }
214     word cellSetName;
215     string vtkName;
217     if (args.optionFound("cellSet"))
218     {
219         cellSetName = args.option("cellSet");
220         vtkName = cellSetName;
221     }
222     else if (Pstream::parRun())
223     {
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)
230         {
231             vtkName = vtkName.substr(i);
232         }
233     }
234     else
235     {
236         vtkName = runTime.caseName();
237     }
240     instantList timeDirs = timeSelector::select0(runTime, args);
242     // Get region name
243     word regionName;
245     if (args.optionReadIfPresent("region", regionName))
246     {
247         Info<< "Create mesh " << regionName << " for time = "
248             << runTime.timeName() << Foam::nl << Foam::endl;
249     }
250     else
251     {
252         regionName = fvMesh::defaultRegion;
253         Info<< "Create mesh for time = "
254             << runTime.timeName() << Foam::nl << Foam::endl;
255     }
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)
263     {
264         fvPath = fvPath/regionName;
265         regionPrefix = regionName;
266     }
268     if (isDir(fvPath))
269     {
270         if
271         (
272             args.optionFound("time")
273          || args.optionFound("latestTime")
274          || cellSetName.size()
275          || regionName != polyMesh::defaultRegion
276         )
277         {
278             Info<< "Keeping old files in " << fvPath << nl << endl;
279         }
280         else
281         {
282             Info<< "Deleting old VTK files in " << fvPath << nl << endl;
284             rmDir(fvPath);
285         }
286     }
288     mkDir(fvPath);
291     // mesh wrapper; does subsetting and decomposition
292     vtkMesh vMesh
293     (
294         Foam::IOobject
295         (
296             regionName,
297             runTime.timeName(),
298             runTime,
299             Foam::IOobject::MUST_READ
300         ),
301         cellSetName
302     );
304     forAll(timeDirs, timeI)
305     {
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
313         // decomposition.
314         polyMesh::readUpdateState meshState = vMesh.readUpdate();
316         const fvMesh& mesh = vMesh.mesh();
318         INTEGER4 nFaceNodes = 0;
319         forAll(mesh.faces(), faceI)
320         {
321             nFaceNodes += mesh.faces()[faceI].size();
322         }
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"))
333         {
334             args.optionLookup("fields")() >> selectedFields;
335         }
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)
362         if (noPointValues)
363         {
364             Info<< "    pointScalarFields : switched off"
365                 << " (\"-noPointValues\" option)\n";
366             Info<< "    pointVectorFields : switched off"
367                 << " (\"-noPointValues\" option)\n";
368         }
370         PtrList<pointScalarField> psf;
371         PtrList<pointVectorField> pvf;
372         //PtrList<pointSphericalTensorField> pSpheretf;
373         //PtrList<pointSymmTensorField> pSymmtf;
374         //PtrList<pointTensorField> ptf;
377         if (!noPointValues)
378         {
379             //// Add interpolated volFields
380             //const volPointInterpolation& pInterp = volPointInterpolation::New
381             //(
382             //    mesh
383             //);
384             //
385             //label nPsf = psf.size();
386             //psf.setSize(nPsf+vsf.size());
387             //forAll(vsf, i)
388             //{
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);
393             //    nPsf++;
394             //}
395             //
396             //label nPvf = pvf.size();
397             //pvf.setSize(vvf.size());
398             //forAll(vvf, i)
399             //{
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);
404             //    nPvf++;
405             //}
407             readFields
408             (
409                 vMesh,
410                 pointMesh::New(vMesh),
411                 objects,
412                 selectedFields,
413                 psf
414             );
415             print("    pointScalarFields          :", Info, psf);
417             readFields
418             (
419                 vMesh,
420                 pointMesh::New(vMesh),
421                 objects,
422                 selectedFields,
423                 pvf
424             );
425             print("    pointVectorFields          :", Info, pvf);
427             //readFields
428             //(
429             //    vMesh,
430             //    pointMesh::New(vMesh),
431             //    objects,
432             //    selectedFields,
433             //    pSpheretf
434             //);
435             //print("    pointSphericalTensorFields :", Info, pSpheretf);
436             //
437             //readFields
438             //(
439             //    vMesh,
440             //    pointMesh::New(vMesh),
441             //    objects,
442             //    selectedFields,
443             //    pSymmtf
444             //);
445             //print("    pointSymmTensorFields      :", Info, pSymmtf);
446             //
447             //readFields
448             //(
449             //    vMesh,
450             //    pointMesh::New(vMesh),
451             //    objects,
452             //    selectedFields,
453             //    ptf
454             //);
455             //print("    pointTensorFields          :", Info, ptf);
456         }
457         Info<< endl;
460         // Get field names
461         // ~~~~~~~~~~~~~~~
463         string varNames;
464         DynamicList<INTEGER4> varLocation;
466         string cellVarNames;
467         DynamicList<INTEGER4> cellVarLocation;
469         // volFields
470         tecplotWriter::getTecplotNames
471         (
472             vsf,
473             ValueLocation_CellCentered,
474             varNames,
475             varLocation
476         );
477         tecplotWriter::getTecplotNames
478         (
479             vsf,
480             ValueLocation_CellCentered,
481             cellVarNames,
482             cellVarLocation
483         );
485         tecplotWriter::getTecplotNames
486         (
487             vvf,
488             ValueLocation_CellCentered,
489             varNames,
490             varLocation
491         );
492         tecplotWriter::getTecplotNames
493         (
494             vvf,
495             ValueLocation_CellCentered,
496             cellVarNames,
497             cellVarLocation
498         );
500         tecplotWriter::getTecplotNames
501         (
502             vSpheretf,
503             ValueLocation_CellCentered,
504             varNames,
505             varLocation
506         );
507         tecplotWriter::getTecplotNames
508         (
509             vSpheretf,
510             ValueLocation_CellCentered,
511             cellVarNames,
512             cellVarLocation
513         );
515         tecplotWriter::getTecplotNames
516         (
517             vSymmtf,
518             ValueLocation_CellCentered,
519             varNames,
520             varLocation
521         );
522         tecplotWriter::getTecplotNames
523         (
524             vSymmtf,
525             ValueLocation_CellCentered,
526             cellVarNames,
527             cellVarLocation
528         );
530         tecplotWriter::getTecplotNames
531         (
532             vtf,
533             ValueLocation_CellCentered,
534             varNames,
535             varLocation
536         );
537         tecplotWriter::getTecplotNames
538         (
539             vtf,
540             ValueLocation_CellCentered,
541             cellVarNames,
542             cellVarLocation
543         );
547         // pointFields
548         tecplotWriter::getTecplotNames
549         (
550             psf,
551             ValueLocation_Nodal,
552             varNames,
553             varLocation
554         );
556         tecplotWriter::getTecplotNames
557         (
558             pvf,
559             ValueLocation_Nodal,
560             varNames,
561             varLocation
562         );
564         // strandID (= piece id. Gets incremented for every piece of geometry
565         // that is output)
566         INTEGER4 strandID = 1;
569         if (meshState != polyMesh::UNCHANGED)
570         {
571             if (doWriteInternal)
572             {
573                 // Output mesh and fields
574                 fileName vtkFileName
575                 (
576                     fvPath/vtkName
577                   + "_"
578                   + timeDesc
579                   + ".plt"
580                 );
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);
591                 writer.writeInit
592                 (
593                     runTime.caseName(),
594                     allVarNames,
595                     vtkFileName,
596                     DataFileType_Full
597                 );
599                 writer.writePolyhedralZone
600                 (
601                     mesh.name(),        // regionName
602                     strandID++,         // strandID
603                     mesh,
604                     allVarLocation,
605                     nFaceNodes
606                 );
608                 // Write coordinates
609                 writer.writeField(mesh.points().component(0)());
610                 writer.writeField(mesh.points().component(1)());
611                 writer.writeField(mesh.points().component(2)());
613                 // Write all fields
614                 forAll(vsf, i)
615                 {
616                     writer.writeField(vsf[i]);
617                 }
618                 forAll(vvf, i)
619                 {
620                     writer.writeField(vvf[i]);
621                 }
622                 forAll(vSpheretf, i)
623                 {
624                     writer.writeField(vSpheretf[i]);
625                 }
626                 forAll(vSymmtf, i)
627                 {
628                     writer.writeField(vSymmtf[i]);
629                 }
630                 forAll(vtf, i)
631                 {
632                     writer.writeField(vtf[i]);
633                 }
635                 forAll(psf, i)
636                 {
637                     writer.writeField(psf[i]);
638                 }
639                 forAll(pvf, i)
640                 {
641                     writer.writeField(pvf[i]);
642                 }
644                 writer.writeConnectivity(mesh);
645                 writer.writeEnd();
646             }
647         }
648         else
649         {
650             if (doWriteInternal)
651             {
652                 if (timeI == 0)
653                 {
654                     // Output static mesh only
655                     fileName vtkFileName
656                     (
657                         fvPath/vtkName
658                       + "_grid_"
659                       + timeDesc
660                       + ".plt"
661                     );
663                     tecplotWriter writer(runTime);
665                     writer.writeInit
666                     (
667                         runTime.caseName(),
668                         "X Y Z",
669                         vtkFileName,
670                         DataFileType_Grid
671                     );
673                     writer.writePolyhedralZone
674                     (
675                         mesh.name(),        // regionName
676                         strandID,           // strandID
677                         mesh,
678                         List<INTEGER4>(3, ValueLocation_Nodal),
679                         nFaceNodes
680                     );
682                     // Write coordinates
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);
688                     writer.writeEnd();
689                 }
691                 // Output solution file
692                 fileName vtkFileName
693                 (
694                     fvPath/vtkName
695                   + "_"
696                   + timeDesc
697                   + ".plt"
698                 );
700                 tecplotWriter writer(runTime);
702                 writer.writeInit
703                 (
704                     runTime.caseName(),
705                     varNames,
706                     vtkFileName,
707                     DataFileType_Solution
708                 );
710                 writer.writePolyhedralZone
711                 (
712                     mesh.name(),        // regionName
713                     strandID++,         // strandID
714                     mesh,
715                     varLocation,
716                     0
717                 );
719                 // Write all fields
720                 forAll(vsf, i)
721                 {
722                     writer.writeField(vsf[i]);
723                 }
724                 forAll(vvf, i)
725                 {
726                     writer.writeField(vvf[i]);
727                 }
728                 forAll(vSpheretf, i)
729                 {
730                     writer.writeField(vSpheretf[i]);
731                 }
732                 forAll(vSymmtf, i)
733                 {
734                     writer.writeField(vSymmtf[i]);
735                 }
736                 forAll(vtf, i)
737                 {
738                     writer.writeField(vtf[i]);
739                 }
741                 forAll(psf, i)
742                 {
743                     writer.writeField(psf[i]);
744                 }
745                 forAll(pvf, i)
746                 {
747                     writer.writeField(pvf[i]);
748                 }
749                 writer.writeEnd();
750             }
751         }
754         //---------------------------------------------------------------------
755         //
756         // Write faceSet
757         //
758         //---------------------------------------------------------------------
760         if (args.optionFound("faceSet"))
761         {
762             // Load the 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
770             (
771                 fvPath/setName/setName
772               + "_"
773               + timeDesc
774               + ".plt"
775             );
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);
788             writer.writeInit
789             (
790                 runTime.caseName(),
791                 cellVarNames,
792                 patchFileName,
793                 DataFileType_Full
794             );
796             const indirectPrimitivePatch ipp
797             (
798                 IndirectList<face>(mesh.faces(), faceLabels),
799                 mesh.points()
800             );
802             writer.writePolygonalZone
803             (
804                 setName,
805                 strandID++,
806                 ipp,
807                 allVarLocation
808             );
810             // Write coordinates
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
816             forAll(vsf, i)
817             {
818                 writer.writeField
819                 (
820                     writer.getFaceField
821                     (
822                         linearInterpolate(vsf[i])(),
823                         faceLabels
824                     )()
825                 );
826             }
827             forAll(vvf, i)
828             {
829                 writer.writeField
830                 (
831                     writer.getFaceField
832                     (
833                         linearInterpolate(vvf[i])(),
834                         faceLabels
835                     )()
836                 );
837             }
838             forAll(vSpheretf, i)
839             {
840                 writer.writeField
841                 (
842                     writer.getFaceField
843                     (
844                         linearInterpolate(vSpheretf[i])(),
845                         faceLabels
846                     )()
847                 );
848             }
849             forAll(vSymmtf, i)
850             {
851                 writer.writeField
852                 (
853                     writer.getFaceField
854                     (
855                         linearInterpolate(vSymmtf[i])(),
856                         faceLabels
857                     )()
858                 );
859             }
860             forAll(vtf, i)
861             {
862                 writer.writeField
863                 (
864                     writer.getFaceField
865                     (
866                         linearInterpolate(vtf[i])(),
867                         faceLabels
868                     )()
869                 );
870             }
871             writer.writeConnectivity(ipp);
873             continue;
874         }
878         //---------------------------------------------------------------------
879         //
880         // Write patches as multi-zone file
881         //
882         //---------------------------------------------------------------------
884         const polyBoundaryMesh& patches = mesh.boundaryMesh();
886         labelList patchIDs(getSelectedPatches(patches, excludePatches));
888         mkDir(fvPath/"boundaryMesh");
890         fileName patchFileName;
892         if (vMesh.useSubMesh())
893         {
894             patchFileName =
895                 fvPath/"boundaryMesh"/cellSetName
896               + "_"
897               + timeDesc
898               + ".plt";
899         }
900         else
901         {
902             patchFileName =
903                 fvPath/"boundaryMesh"/"boundaryMesh"
904               + "_"
905               + timeDesc
906               + ".plt";
907         }
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);
920         writer.writeInit
921         (
922             runTime.caseName(),
923             allVarNames,
924             patchFileName,
925             DataFileType_Full
926         );
928         forAll(patchIDs, i)
929         {
930             label patchID = patchIDs[i];
931             const polyPatch& pp = patches[patchID];
932             //INTEGER4 strandID = 1 + i;
934             if (pp.size() > 0)
935             {
936                 Info<< "    Writing patch " << patchID << "\t" << pp.name()
937                     << "\tstrand:" << strandID << nl << endl;
939                 const indirectPrimitivePatch ipp
940                 (
941                     IndirectList<face>(pp, identity(pp.size())),
942                     pp.points()
943                 );
945                 writer.writePolygonalZone
946                 (
947                     pp.name(),
948                     strandID++,     //strandID,
949                     ipp,
950                     allVarLocation
951                 );
953                 // Write coordinates
954                 writer.writeField(ipp.localPoints().component(0)());
955                 writer.writeField(ipp.localPoints().component(1)());
956                 writer.writeField(ipp.localPoints().component(2)());
958                 // Write all fields
959                 forAll(vsf, i)
960                 {
961                     writer.writeField
962                     (
963                         writer.getPatchField
964                         (
965                             nearCellValue,
966                             vsf[i],
967                             patchID
968                         )()
969                     );
970                 }
971                 forAll(vvf, i)
972                 {
973                     writer.writeField
974                     (
975                         writer.getPatchField
976                         (
977                             nearCellValue,
978                             vvf[i],
979                             patchID
980                         )()
981                     );
982                 }
983                 forAll(vSpheretf, i)
984                 {
985                     writer.writeField
986                     (
987                         writer.getPatchField
988                         (
989                             nearCellValue,
990                             vSpheretf[i],
991                             patchID
992                         )()
993                     );
994                 }
995                 forAll(vSymmtf, i)
996                 {
997                     writer.writeField
998                     (
999                         writer.getPatchField
1000                         (
1001                             nearCellValue,
1002                             vSymmtf[i],
1003                             patchID
1004                         )()
1005                     );
1006                 }
1007                 forAll(vtf, i)
1008                 {
1009                     writer.writeField
1010                     (
1011                         writer.getPatchField
1012                         (
1013                             nearCellValue,
1014                             vtf[i],
1015                             patchID
1016                         )()
1017                     );
1018                 }
1020                 forAll(psf, i)
1021                 {
1022                     writer.writeField
1023                     (
1024                         psf[i].boundaryField()[patchID].patchInternalField()()
1025                     );
1026                 }
1027                 forAll(pvf, i)
1028                 {
1029                     writer.writeField
1030                     (
1031                         pvf[i].boundaryField()[patchID].patchInternalField()()
1032                     );
1033                 }
1035                 writer.writeConnectivity(ipp);
1036             }
1037             else
1038             {
1039                 Info<< "    Skipping zero sized patch " << patchID
1040                     << "\t" << pp.name()
1041                     << nl << endl;
1042             }
1043         }
1044         writer.writeEnd();
1045         Info<< endl;
1049         //---------------------------------------------------------------------
1050         //
1051         // Write faceZones as multi-zone file
1052         //
1053         //---------------------------------------------------------------------
1055         const faceZoneMesh& zones = mesh.faceZones();
1057         if (doFaceZones && zones.size() > 0)
1058         {
1059             mkDir(fvPath/"faceZoneMesh");
1061             fileName patchFileName;
1063             if (vMesh.useSubMesh())
1064             {
1065                 patchFileName =
1066                     fvPath/"faceZoneMesh"/cellSetName
1067                   + "_"
1068                   + timeDesc
1069                   + ".plt";
1070             }
1071             else
1072             {
1073                 patchFileName =
1074                     fvPath/"faceZoneMesh"/"faceZoneMesh"
1075                   + "_"
1076                   + timeDesc
1077                   + ".plt";
1078             }
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);
1091             writer.writeInit
1092             (
1093                 runTime.caseName(),
1094                 allVarNames,
1095                 patchFileName,
1096                 DataFileType_Full
1097             );
1099             forAll(zones, zoneI)
1100             {
1101                 const faceZone& pp = zones[zoneI];
1103                 if (pp.size() > 0)
1104                 {
1105                     const indirectPrimitivePatch ipp
1106                     (
1107                         IndirectList<face>(mesh.faces(), pp),
1108                         mesh.points()
1109                     );
1111                     writer.writePolygonalZone
1112                     (
1113                         pp.name(),
1114                         strandID++, //1+patchIDs.size()+zoneI,    //strandID,
1115                         ipp,
1116                         allVarLocation
1117                     );
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
1125                     forAll(vsf, i)
1126                     {
1127                         writer.writeField
1128                         (
1129                             writer.getFaceField
1130                             (
1131                                 linearInterpolate(vsf[i])(),
1132                                 pp
1133                             )()
1134                         );
1135                     }
1136                     forAll(vvf, i)
1137                     {
1138                         writer.writeField
1139                         (
1140                             writer.getFaceField
1141                             (
1142                                 linearInterpolate(vvf[i])(),
1143                                 pp
1144                             )()
1145                         );
1146                     }
1147                     forAll(vSpheretf, i)
1148                     {
1149                         writer.writeField
1150                         (
1151                             writer.getFaceField
1152                             (
1153                                 linearInterpolate(vSpheretf[i])(),
1154                                 pp
1155                             )()
1156                         );
1157                     }
1158                     forAll(vSymmtf, i)
1159                     {
1160                         writer.writeField
1161                         (
1162                             writer.getFaceField
1163                             (
1164                                 linearInterpolate(vSymmtf[i])(),
1165                                 pp
1166                             )()
1167                         );
1168                     }
1169                     forAll(vtf, i)
1170                     {
1171                         writer.writeField
1172                         (
1173                             writer.getFaceField
1174                             (
1175                                 linearInterpolate(vtf[i])(),
1176                                 pp
1177                             )()
1178                         );
1179                     }
1181                     writer.writeConnectivity(ipp);
1182                 }
1183                 else
1184                 {
1185                     Info<< "    Skipping zero sized faceZone " << zoneI
1186                         << "\t" << pp.name()
1187                         << nl << endl;
1188                 }
1189             }
1190             writer.writeEnd();
1191             Info<< endl;
1192         }
1196         //---------------------------------------------------------------------
1197         //
1198         // Write lagrangian data
1199         //
1200         //---------------------------------------------------------------------
1202         fileNameList cloudDirs
1203         (
1204             readDir
1205             (
1206                 runTime.timePath()/regionPrefix/cloud::prefix,
1207                 fileName::DIRECTORY
1208             )
1209         );
1211         forAll(cloudDirs, cloudI)
1212         {
1213             IOobjectList sprayObjs
1214             (
1215                 mesh,
1216                 runTime.timeName(),
1217                 cloud::prefix/cloudDirs[cloudI]
1218             );
1220             IOobject* positionsPtr = sprayObjs.lookup("positions");
1222             if (positionsPtr)
1223             {
1224                 mkDir(fvPath/cloud::prefix/cloudDirs[cloudI]);
1226                 fileName lagrFileName
1227                 (
1228                     fvPath/cloud::prefix/cloudDirs[cloudI]/cloudDirs[cloudI]
1229                   + "_" + timeDesc + ".plt"
1230                 );
1232                 Info<< "    Lagrangian: " << lagrFileName << endl;
1234                 wordList labelNames(sprayObjs.names(labelIOField::typeName));
1235                 Info<< "        labels            :";
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
1247                 //(
1248                 //    sprayObjs.names
1249                 //    (
1250                 //        sphericalTensorIOField::typeName
1251                 //    )
1252                 //);
1253                 //Info<< "        spherical tensors :";
1254                 //print(Info, sphereNames);
1255                 //
1256                 //wordList symmNames
1257                 //(
1258                 //    sprayObjs.names
1259                 //    (
1260                 //        symmTensorIOField::typeName
1261                 //    )
1262                 //);
1263                 //Info<< "        symm tensors      :";
1264                 //print(Info, symmNames);
1265                 //
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());
1276                 label n = 0;
1277                 forAllConstIter(passiveParticleCloud, parcels, elmnt)
1278                 {
1279                     positions[n++] = elmnt().position();
1280                 }
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>
1290                 (
1291                     labelNames,
1292                     ValueLocation_Nodal,
1293                     allVarNames,
1294                     allVarLocation
1295                 );
1297                 tecplotWriter::getTecplotNames<scalar>
1298                 (
1299                     scalarNames,
1300                     ValueLocation_Nodal,
1301                     allVarNames,
1302                     allVarLocation
1303                 );
1305                 tecplotWriter::getTecplotNames<vector>
1306                 (
1307                     vectorNames,
1308                     ValueLocation_Nodal,
1309                     allVarNames,
1310                     allVarLocation
1311                 );
1314                 tecplotWriter writer(runTime);
1316                 writer.writeInit
1317                 (
1318                     runTime.caseName(),
1319                     allVarNames,
1320                     lagrFileName,
1321                     DataFileType_Full
1322                 );
1324                 writer.writeOrderedZone
1325                 (
1326                     cloudDirs[cloudI],
1327                     strandID++,     //strandID,
1328                     parcels.size(),
1329                     allVarLocation
1330                 );
1332                 // Write coordinates
1333                 writer.writeField(positions.component(0)());
1334                 writer.writeField(positions.component(1)());
1335                 writer.writeField(positions.component(2)());
1337                 // labelFields
1338                 forAll(labelNames, i)
1339                 {
1340                     IOField<label> fld
1341                     (
1342                         IOobject
1343                         (
1344                             labelNames[i],
1345                             mesh.time().timeName(),
1346                             cloud::prefix/cloudDirs[cloudI],
1347                             mesh,
1348                             IOobject::MUST_READ,
1349                             IOobject::NO_WRITE,
1350                             false
1351                         )
1352                     );
1354                     scalarField sfld(fld.size());
1355                     forAll(fld, j)
1356                     {
1357                         sfld[j] = scalar(fld[j]);
1358                     }
1359                     writer.writeField(sfld);
1360                 }
1361                 // scalarFields
1362                 forAll(scalarNames, i)
1363                 {
1364                     IOField<scalar> fld
1365                     (
1366                         IOobject
1367                         (
1368                             scalarNames[i],
1369                             mesh.time().timeName(),
1370                             cloud::prefix/cloudDirs[cloudI],
1371                             mesh,
1372                             IOobject::MUST_READ,
1373                             IOobject::NO_WRITE,
1374                             false
1375                         )
1376                     );
1377                     writer.writeField(fld);
1378                 }
1379                 // vectorFields
1380                 forAll(vectorNames, i)
1381                 {
1382                     IOField<vector> fld
1383                     (
1384                         IOobject
1385                         (
1386                             vectorNames[i],
1387                             mesh.time().timeName(),
1388                             cloud::prefix/cloudDirs[cloudI],
1389                             mesh,
1390                             IOobject::MUST_READ,
1391                             IOobject::NO_WRITE,
1392                             false
1393                         )
1394                     );
1395                     writer.writeField(fld);
1396                 }
1398                 writer.writeEnd();
1399             }
1400         }
1401     }
1403     Info<< "End\n" << endl;
1405     return 0;
1409 // ************************************************************************* //