Changed order of include files; removed extraneous whitespace
[OpenFOAM-1.6.x.git] / applications / utilities / postProcessing / dataConversion / foamToTecplot360 / foamToTecplot360.C
blobc491e5993cd1a22704248522204f75f0abe12875
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
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
19     for more details.
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
25 Application
26     foamToTecplot360
28 Description
29     Tecplot binary file format writer.
31 Usage
33     - foamToTecplot360 [OPTION]
35     @param -fields \<fields\>\n
36     Convert selected fields only. For example,
37     @verbatim
38          -fields '( p T U )'
39     @endverbatim
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
50     @param -noInternal \n
51     Do not generate file for mesh, only for patches
53     @param -noPointValues \n
54     No pointFields
56     @param -noFaceZones \n
57     No faceZones
59     @param -excludePatches \<patchNames\>\n
60     Specify patches (wildcards) to exclude. For example,
61     @verbatim
62          -excludePatches '( inlet_1 inlet_2 "proc.*")'
63     @endverbatim
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
66     expression.
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"
82 #include "faceSet.H"
83 #include "stringListOps.H"
84 #include "wordRe.H"
86 #include "vtkMesh.H"
87 #include "readFields.H"
88 #include "tecplotWriter.H"
90 #include "TECIO.h"
92 // Note: needs to be after TECIO to prevent Foam::Time conflicting with
93 // Xlib Time.
94 #include "fvCFD.H"
96 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
98 template<class GeoField>
99 void print(const char* msg, Ostream& os, const PtrList<GeoField>& flds)
101     if (flds.size())
102     {
103         os << msg;
104         forAll(flds, i)
105         {
106             os<< ' ' << flds[i].name();
107         }
108         os << endl;
109     }
113 void print(Ostream& os, const wordList& flds)
115     forAll(flds, i)
116     {
117         os<< ' ' << flds[i];
118     }
119     os << endl;
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)
134     {
135         const polyPatch& pp = patches[patchI];
137         if
138         (
139             isType<emptyPolyPatch>(pp)
140             || (Pstream::parRun() && isType<processorPolyPatch>(pp))
141         )
142         {
143             Info<< "    discarding empty/processor patch " << patchI
144                 << " " << pp.name() << endl;
145         }
146         else if (findStrings(excludePatches, pp.name()))
147         {
148             Info<< "    excluding patch " << patchI
149                 << " " << pp.name() << endl;
150         }
151         else
152         {
153             patchIDs.append(patchI);
154             Info<< "    patch " << patchI << " " << pp.name() << endl;
155         }
156     }
157     return patchIDs.shrink();
163 // Main program:
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
178     (
179         "excludePatches",
180         "patches (wildcards) to exclude"
181     );
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");
192     if (nearCellValue)
193     {
194         WarningIn(args.executable())
195             << "Using neighbouring cell value instead of patch value"
196             << nl << endl;
197     }
199     bool noPointValues = args.optionFound("noPointValues");
201     if (noPointValues)
202     {
203         WarningIn(args.executable())
204             << "Outputting cell values only" << nl << endl;
205     }
207     List<wordRe> excludePatches;
208     if (args.optionFound("excludePatches"))
209     {
210         args.optionLookup("excludePatches")() >> excludePatches;
212         Info<< "Not including patches " << excludePatches << nl << endl;
213     }
215     word cellSetName;
216     string vtkName;
218     if (args.optionFound("cellSet"))
219     {
220         cellSetName = args.option("cellSet");
221         vtkName = cellSetName;
222     }
223     else if (Pstream::parRun())
224     {
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)
231         {
232             vtkName = vtkName.substr(i);
233         }
234     }
235     else
236     {
237         vtkName = runTime.caseName();
238     }
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)
251     {
252         fvPath = fvPath/regionName;
253         regionPrefix = regionName;
254     }
256     if (isDir(fvPath))
257     {
258         if
259         (
260             args.optionFound("time")
261          || args.optionFound("latestTime")
262          || cellSetName.size()
263          || regionName != polyMesh::defaultRegion
264         )
265         {
266             Info<< "Keeping old files in " << fvPath << nl << endl;
267         }
268         else
269         {
270             Info<< "Deleting old VTK files in " << fvPath << nl << endl;
272             rmDir(fvPath);
273         }
274     }
276     mkDir(fvPath);
279     // mesh wrapper; does subsetting and decomposition
280     vtkMesh vMesh(mesh, cellSetName);
282     forAll(timeDirs, timeI)
283     {
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
291         // decomposition.
292         polyMesh::readUpdateState meshState = vMesh.readUpdate();
294         const fvMesh& mesh = vMesh.mesh();
296         INTEGER4 nFaceNodes = 0;
297         forAll(mesh.faces(), faceI)
298         {
299             nFaceNodes += mesh.faces()[faceI].size();
300         }
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"))
311         {
312             args.optionLookup("fields")() >> selectedFields;
313         }
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)
340         if (noPointValues)
341         {
342             Info<< "    pointScalarFields : switched off"
343                 << " (\"-noPointValues\" option)\n";
344             Info<< "    pointVectorFields : switched off"
345                 << " (\"-noPointValues\" option)\n";
346         }
348         PtrList<pointScalarField> psf;
349         PtrList<pointVectorField> pvf;
350         //PtrList<pointSphericalTensorField> pSpheretf;
351         //PtrList<pointSymmTensorField> pSymmtf;
352         //PtrList<pointTensorField> ptf;
355         if (!noPointValues)
356         {
357             //// Add interpolated volFields
358             //const volPointInterpolation& pInterp = volPointInterpolation::New
359             //(
360             //    mesh
361             //);
362             //
363             //label nPsf = psf.size();
364             //psf.setSize(nPsf+vsf.size());
365             //forAll(vsf, i)
366             //{
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);
371             //    nPsf++;
372             //}
373             //
374             //label nPvf = pvf.size();
375             //pvf.setSize(vvf.size());
376             //forAll(vvf, i)
377             //{
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);
382             //    nPvf++;
383             //}
385             readFields
386             (
387                 vMesh,
388                 pointMesh::New(vMesh.baseMesh()),
389                 objects,
390                 selectedFields,
391                 psf
392             );
393             print("    pointScalarFields          :", Info, psf);
395             readFields
396             (
397                 vMesh,
398                 pointMesh::New(vMesh.baseMesh()),
399                 objects,
400                 selectedFields,
401                 pvf
402             );
403             print("    pointVectorFields          :", Info, pvf);
405             //readFields
406             //(
407             //    vMesh,
408             //    pointMesh::New(vMesh.baseMesh()),
409             //    objects,
410             //    selectedFields,
411             //    pSpheretf
412             //);
413             //print("    pointSphericalTensorFields :", Info, pSpheretf);
414             //
415             //readFields
416             //(
417             //    vMesh,
418             //    pointMesh::New(vMesh.baseMesh()),
419             //    objects,
420             //    selectedFields,
421             //    pSymmtf
422             //);
423             //print("    pointSymmTensorFields      :", Info, pSymmtf);
424             //
425             //readFields
426             //(
427             //    vMesh,
428             //    pointMesh::New(vMesh.baseMesh()),
429             //    objects,
430             //    selectedFields,
431             //    ptf
432             //);
433             //print("    pointTensorFields          :", Info, ptf);
434         }
435         Info<< endl;
438         // Get field names
439         // ~~~~~~~~~~~~~~~
441         string varNames;
442         DynamicList<INTEGER4> varLocation;
444         string cellVarNames;
445         DynamicList<INTEGER4> cellVarLocation;
447         // volFields
448         tecplotWriter::getTecplotNames
449         (
450             vsf,
451             ValueLocation_CellCentered,
452             varNames,
453             varLocation
454         );
455         tecplotWriter::getTecplotNames
456         (
457             vsf,
458             ValueLocation_CellCentered,
459             cellVarNames,
460             cellVarLocation
461         );
463         tecplotWriter::getTecplotNames
464         (
465             vvf,
466             ValueLocation_CellCentered,
467             varNames,
468             varLocation
469         );
470         tecplotWriter::getTecplotNames
471         (
472             vvf,
473             ValueLocation_CellCentered,
474             cellVarNames,
475             cellVarLocation
476         );
478         tecplotWriter::getTecplotNames
479         (
480             vSpheretf,
481             ValueLocation_CellCentered,
482             varNames,
483             varLocation
484         );
485         tecplotWriter::getTecplotNames
486         (
487             vSpheretf,
488             ValueLocation_CellCentered,
489             cellVarNames,
490             cellVarLocation
491         );
493         tecplotWriter::getTecplotNames
494         (
495             vSymmtf,
496             ValueLocation_CellCentered,
497             varNames,
498             varLocation
499         );
500         tecplotWriter::getTecplotNames
501         (
502             vSymmtf,
503             ValueLocation_CellCentered,
504             cellVarNames,
505             cellVarLocation
506         );
508         tecplotWriter::getTecplotNames
509         (
510             vtf,
511             ValueLocation_CellCentered,
512             varNames,
513             varLocation
514         );
515         tecplotWriter::getTecplotNames
516         (
517             vtf,
518             ValueLocation_CellCentered,
519             cellVarNames,
520             cellVarLocation
521         );
525         // pointFields
526         tecplotWriter::getTecplotNames
527         (
528             psf,
529             ValueLocation_Nodal,
530             varNames,
531             varLocation
532         );
534         tecplotWriter::getTecplotNames
535         (
536             pvf,
537             ValueLocation_Nodal,
538             varNames,
539             varLocation
540         );
542         // strandID (= piece id. Gets incremented for every piece of geometry
543         // that is output)
544         INTEGER4 strandID = 1;
547         if (meshState != polyMesh::UNCHANGED)
548         {
549             if (doWriteInternal)
550             {
551                 // Output mesh and fields
552                 fileName vtkFileName
553                 (
554                     fvPath/vtkName
555                   + "_"
556                   + timeDesc
557                   + ".plt"
558                 );
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);
569                 writer.writeInit
570                 (
571                     runTime.caseName(),
572                     allVarNames,
573                     vtkFileName,
574                     DataFileType_Full
575                 );
577                 writer.writePolyhedralZone
578                 (
579                     mesh.name(),        // regionName
580                     strandID++,         // strandID
581                     mesh,
582                     allVarLocation,
583                     nFaceNodes
584                 );
586                 // Write coordinates
587                 writer.writeField(mesh.points().component(0)());
588                 writer.writeField(mesh.points().component(1)());
589                 writer.writeField(mesh.points().component(2)());
591                 // Write all fields
592                 forAll(vsf, i)
593                 {
594                     writer.writeField(vsf[i]);
595                 }
596                 forAll(vvf, i)
597                 {
598                     writer.writeField(vvf[i]);
599                 }
600                 forAll(vSpheretf, i)
601                 {
602                     writer.writeField(vSpheretf[i]);
603                 }
604                 forAll(vSymmtf, i)
605                 {
606                     writer.writeField(vSymmtf[i]);
607                 }
608                 forAll(vtf, i)
609                 {
610                     writer.writeField(vtf[i]);
611                 }
613                 forAll(psf, i)
614                 {
615                     writer.writeField(psf[i]);
616                 }
617                 forAll(pvf, i)
618                 {
619                     writer.writeField(pvf[i]);
620                 }
622                 writer.writeConnectivity(mesh);
623                 writer.writeEnd();
624             }
625         }
626         else
627         {
628             if (doWriteInternal)
629             {
630                 if (timeI == 0)
631                 {
632                     // Output static mesh only
633                     fileName vtkFileName
634                     (
635                         fvPath/vtkName
636                       + "_grid_"
637                       + timeDesc
638                       + ".plt"
639                     );
641                     tecplotWriter writer(runTime);
643                     writer.writeInit
644                     (
645                         runTime.caseName(),
646                         "X Y Z",
647                         vtkFileName,
648                         DataFileType_Grid
649                     );
651                     writer.writePolyhedralZone
652                     (
653                         mesh.name(),        // regionName
654                         strandID,           // strandID
655                         mesh,
656                         List<INTEGER4>(3, ValueLocation_Nodal),
657                         nFaceNodes
658                     );
660                     // Write coordinates
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);
666                     writer.writeEnd();
667                 }
669                 // Output solution file
670                 fileName vtkFileName
671                 (
672                     fvPath/vtkName
673                   + "_"
674                   + timeDesc
675                   + ".plt"
676                 );
678                 tecplotWriter writer(runTime);
680                 writer.writeInit
681                 (
682                     runTime.caseName(),
683                     varNames,
684                     vtkFileName,
685                     DataFileType_Solution
686                 );
688                 writer.writePolyhedralZone
689                 (
690                     mesh.name(),        // regionName
691                     strandID++,         // strandID
692                     mesh,
693                     varLocation,
694                     0
695                 );
697                 // Write all fields
698                 forAll(vsf, i)
699                 {
700                     writer.writeField(vsf[i]);
701                 }
702                 forAll(vvf, i)
703                 {
704                     writer.writeField(vvf[i]);
705                 }
706                 forAll(vSpheretf, i)
707                 {
708                     writer.writeField(vSpheretf[i]);
709                 }
710                 forAll(vSymmtf, i)
711                 {
712                     writer.writeField(vSymmtf[i]);
713                 }
714                 forAll(vtf, i)
715                 {
716                     writer.writeField(vtf[i]);
717                 }
719                 forAll(psf, i)
720                 {
721                     writer.writeField(psf[i]);
722                 }
723                 forAll(pvf, i)
724                 {
725                     writer.writeField(pvf[i]);
726                 }
727                 writer.writeEnd();
728             }
729         }
732         //---------------------------------------------------------------------
733         //
734         // Write faceSet
735         //
736         //---------------------------------------------------------------------
738         if (args.optionFound("faceSet"))
739         {
740             // Load the 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
748             (
749                 fvPath/setName/setName
750               + "_"
751               + timeDesc
752               + ".plt"
753             );
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);
766             writer.writeInit
767             (
768                 runTime.caseName(),
769                 cellVarNames,
770                 patchFileName,
771                 DataFileType_Full
772             );
774             const indirectPrimitivePatch ipp
775             (
776                 IndirectList<face>(mesh.faces(), faceLabels),
777                 mesh.points()
778             );
780             writer.writePolygonalZone
781             (
782                 setName,
783                 strandID++,
784                 ipp,
785                 allVarLocation
786             );
788             // Write coordinates
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
794             forAll(vsf, i)
795             {
796                 writer.writeField
797                 (
798                     writer.getFaceField
799                     (
800                         linearInterpolate(vsf[i])(),
801                         faceLabels
802                     )()
803                 );
804             }
805             forAll(vvf, i)
806             {
807                 writer.writeField
808                 (
809                     writer.getFaceField
810                     (
811                         linearInterpolate(vvf[i])(),
812                         faceLabels
813                     )()
814                 );
815             }
816             forAll(vSpheretf, i)
817             {
818                 writer.writeField
819                 (
820                     writer.getFaceField
821                     (
822                         linearInterpolate(vSpheretf[i])(),
823                         faceLabels
824                     )()
825                 );
826             }
827             forAll(vSymmtf, i)
828             {
829                 writer.writeField
830                 (
831                     writer.getFaceField
832                     (
833                         linearInterpolate(vSymmtf[i])(),
834                         faceLabels
835                     )()
836                 );
837             }
838             forAll(vtf, i)
839             {
840                 writer.writeField
841                 (
842                     writer.getFaceField
843                     (
844                         linearInterpolate(vtf[i])(),
845                         faceLabels
846                     )()
847                 );
848             }
849             writer.writeConnectivity(ipp);
851             continue;
852         }
856         //---------------------------------------------------------------------
857         //
858         // Write patches as multi-zone file
859         //
860         //---------------------------------------------------------------------
862         const polyBoundaryMesh& patches = mesh.boundaryMesh();
864         labelList patchIDs(getSelectedPatches(patches, excludePatches));
866         mkDir(fvPath/"boundaryMesh");
868         fileName patchFileName;
870         if (vMesh.useSubMesh())
871         {
872             patchFileName =
873                 fvPath/"boundaryMesh"/cellSetName
874               + "_"
875               + timeDesc
876               + ".plt";
877         }
878         else
879         {
880             patchFileName =
881                 fvPath/"boundaryMesh"/"boundaryMesh"
882               + "_"
883               + timeDesc
884               + ".plt";
885         }
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);
898         writer.writeInit
899         (
900             runTime.caseName(),
901             allVarNames,
902             patchFileName,
903             DataFileType_Full
904         );
906         forAll(patchIDs, i)
907         {
908             label patchID = patchIDs[i];
909             const polyPatch& pp = patches[patchID];
910             //INTEGER4 strandID = 1 + i;
912             if (pp.size() > 0)
913             {
914                 Info<< "    Writing patch " << patchID << "\t" << pp.name()
915                     << "\tstrand:" << strandID << nl << endl;
917                 const indirectPrimitivePatch ipp
918                 (
919                     IndirectList<face>(pp, identity(pp.size())),
920                     pp.points()
921                 );
923                 writer.writePolygonalZone
924                 (
925                     pp.name(),
926                     strandID++,     //strandID,
927                     ipp,
928                     allVarLocation
929                 );
931                 // Write coordinates
932                 writer.writeField(ipp.localPoints().component(0)());
933                 writer.writeField(ipp.localPoints().component(1)());
934                 writer.writeField(ipp.localPoints().component(2)());
936                 // Write all fields
937                 forAll(vsf, i)
938                 {
939                     writer.writeField
940                     (
941                         writer.getPatchField
942                         (
943                             nearCellValue,
944                             vsf[i],
945                             patchID
946                         )()
947                     );
948                 }
949                 forAll(vvf, i)
950                 {
951                     writer.writeField
952                     (
953                         writer.getPatchField
954                         (
955                             nearCellValue,
956                             vvf[i],
957                             patchID
958                         )()
959                     );
960                 }
961                 forAll(vSpheretf, i)
962                 {
963                     writer.writeField
964                     (
965                         writer.getPatchField
966                         (
967                             nearCellValue,
968                             vSpheretf[i],
969                             patchID
970                         )()
971                     );
972                 }
973                 forAll(vSymmtf, i)
974                 {
975                     writer.writeField
976                     (
977                         writer.getPatchField
978                         (
979                             nearCellValue,
980                             vSymmtf[i],
981                             patchID
982                         )()
983                     );
984                 }
985                 forAll(vtf, i)
986                 {
987                     writer.writeField
988                     (
989                         writer.getPatchField
990                         (
991                             nearCellValue,
992                             vtf[i],
993                             patchID
994                         )()
995                     );
996                 }
998                 forAll(psf, i)
999                 {
1000                     writer.writeField
1001                     (
1002                         psf[i].boundaryField()[patchID].patchInternalField()()
1003                     );
1004                 }
1005                 forAll(pvf, i)
1006                 {
1007                     writer.writeField
1008                     (
1009                         pvf[i].boundaryField()[patchID].patchInternalField()()
1010                     );
1011                 }
1013                 writer.writeConnectivity(ipp);
1014             }
1015             else
1016             {
1017                 Info<< "    Skipping zero sized patch " << patchID
1018                     << "\t" << pp.name()
1019                     << nl << endl;
1020             }
1021         }
1022         writer.writeEnd();
1026         //---------------------------------------------------------------------
1027         //
1028         // Write faceZones as multi-zone file
1029         //
1030         //---------------------------------------------------------------------
1032         const faceZoneMesh& zones = mesh.faceZones();
1034         if (doFaceZones && zones.size() > 0)
1035         {
1036             mkDir(fvPath/"faceZoneMesh");
1038             fileName patchFileName;
1040             if (vMesh.useSubMesh())
1041             {
1042                 patchFileName =
1043                     fvPath/"faceZoneMesh"/cellSetName
1044                   + "_"
1045                   + timeDesc
1046                   + ".plt";
1047             }
1048             else
1049             {
1050                 patchFileName =
1051                     fvPath/"faceZoneMesh"/"faceZoneMesh"
1052                   + "_"
1053                   + timeDesc
1054                   + ".plt";
1055             }
1057             Info<< "    FaceZone  : " << patchFileName << endl;
1059             tecplotWriter writer(runTime);
1061             string allVarNames = string("X Y Z ") + cellVarNames;
1062             DynamicList<INTEGER4> allVarLocation;
1063             allVarLocation.append(ValueLocation_Nodal);
1064             allVarLocation.append(ValueLocation_Nodal);
1065             allVarLocation.append(ValueLocation_Nodal);
1066             allVarLocation.append(cellVarLocation);
1068             writer.writeInit
1069             (
1070                 runTime.caseName(),
1071                 cellVarNames,
1072                 patchFileName,
1073                 DataFileType_Full
1074             );
1076             forAll(zones, zoneI)
1077             {
1078                 const faceZone& pp = zones[zoneI];
1080                 if (pp.size() > 0)
1081                 {
1082                     const indirectPrimitivePatch ipp
1083                     (
1084                         IndirectList<face>(mesh.faces(), pp),
1085                         mesh.points()
1086                     );
1088                     writer.writePolygonalZone
1089                     (
1090                         pp.name(),
1091                         strandID++, //1+patchIDs.size()+zoneI,    //strandID,
1092                         ipp,
1093                         allVarLocation
1094                     );
1096                     // Write coordinates
1097                     writer.writeField(ipp.localPoints().component(0)());
1098                     writer.writeField(ipp.localPoints().component(1)());
1099                     writer.writeField(ipp.localPoints().component(2)());
1101                     // Write all volfields
1102                     forAll(vsf, i)
1103                     {
1104                         writer.writeField
1105                         (
1106                             writer.getFaceField
1107                             (
1108                                 linearInterpolate(vsf[i])(),
1109                                 pp
1110                             )()
1111                         );
1112                     }
1113                     forAll(vvf, i)
1114                     {
1115                         writer.writeField
1116                         (
1117                             writer.getFaceField
1118                             (
1119                                 linearInterpolate(vvf[i])(),
1120                                 pp
1121                             )()
1122                         );
1123                     }
1124                     forAll(vSpheretf, i)
1125                     {
1126                         writer.writeField
1127                         (
1128                             writer.getFaceField
1129                             (
1130                                 linearInterpolate(vSpheretf[i])(),
1131                                 pp
1132                             )()
1133                         );
1134                     }
1135                     forAll(vSymmtf, i)
1136                     {
1137                         writer.writeField
1138                         (
1139                             writer.getFaceField
1140                             (
1141                                 linearInterpolate(vSymmtf[i])(),
1142                                 pp
1143                             )()
1144                         );
1145                     }
1146                     forAll(vtf, i)
1147                     {
1148                         writer.writeField
1149                         (
1150                             writer.getFaceField
1151                             (
1152                                 linearInterpolate(vtf[i])(),
1153                                 pp
1154                             )()
1155                         );
1156                     }
1158                     writer.writeConnectivity(ipp);
1159                 }
1160                 else
1161                 {
1162                     Info<< "    Skipping zero sized faceZone " << zoneI
1163                         << "\t" << pp.name()
1164                         << nl << endl;
1165                 }
1166             }
1167         }
1171         //---------------------------------------------------------------------
1172         //
1173         // Write lagrangian data
1174         //
1175         //---------------------------------------------------------------------
1177         fileNameList cloudDirs
1178         (
1179             readDir
1180             (
1181                 runTime.timePath()/regionPrefix/cloud::prefix,
1182                 fileName::DIRECTORY
1183             )
1184         );
1186         forAll(cloudDirs, cloudI)
1187         {
1188             IOobjectList sprayObjs
1189             (
1190                 mesh,
1191                 runTime.timeName(),
1192                 cloud::prefix/cloudDirs[cloudI]
1193             );
1195             IOobject* positionsPtr = sprayObjs.lookup("positions");
1197             if (positionsPtr)
1198             {
1199                 mkDir(fvPath/cloud::prefix/cloudDirs[cloudI]);
1201                 fileName lagrFileName
1202                 (
1203                     fvPath/cloud::prefix/cloudDirs[cloudI]/cloudDirs[cloudI]
1204                   + "_" + timeDesc + ".plt"
1205                 );
1207                 Info<< "    Lagrangian: " << lagrFileName << endl;
1209                 wordList labelNames(sprayObjs.names(labelIOField::typeName));
1210                 Info<< "        labels            :";
1211                 print(Info, labelNames);
1213                 wordList scalarNames(sprayObjs.names(scalarIOField::typeName));
1214                 Info<< "        scalars           :";
1215                 print(Info, scalarNames);
1217                 wordList vectorNames(sprayObjs.names(vectorIOField::typeName));
1218                 Info<< "        vectors           :";
1219                 print(Info, vectorNames);
1221                 //wordList sphereNames
1222                 //(
1223                 //    sprayObjs.names
1224                 //    (
1225                 //        sphericalTensorIOField::typeName
1226                 //    )
1227                 //);
1228                 //Info<< "        spherical tensors :";
1229                 //print(Info, sphereNames);
1230                 //
1231                 //wordList symmNames
1232                 //(
1233                 //    sprayObjs.names
1234                 //    (
1235                 //        symmTensorIOField::typeName
1236                 //    )
1237                 //);
1238                 //Info<< "        symm tensors      :";
1239                 //print(Info, symmNames);
1240                 //
1241                 //wordList tensorNames(sprayObjs.names(tensorIOField::typeName));
1242                 //Info<< "        tensors           :";
1243                 //print(Info, tensorNames);
1246                 // Load cloud positions
1247                 passiveParticleCloud parcels(mesh, cloudDirs[cloudI]);
1249                 // Get positions as pointField
1250                 pointField positions(parcels.size());
1251                 label n = 0;
1252                 forAllConstIter(passiveParticleCloud, parcels, elmnt)
1253                 {
1254                     positions[n++] = elmnt().position();
1255                 }
1258                 string allVarNames = string("X Y Z");
1259                 DynamicList<INTEGER4> allVarLocation;
1260                 allVarLocation.append(ValueLocation_Nodal);
1261                 allVarLocation.append(ValueLocation_Nodal);
1262                 allVarLocation.append(ValueLocation_Nodal);
1264                 tecplotWriter::getTecplotNames<label>
1265                 (
1266                     labelNames,
1267                     ValueLocation_Nodal,
1268                     allVarNames,
1269                     allVarLocation
1270                 );
1272                 tecplotWriter::getTecplotNames<scalar>
1273                 (
1274                     scalarNames,
1275                     ValueLocation_Nodal,
1276                     allVarNames,
1277                     allVarLocation
1278                 );
1280                 tecplotWriter::getTecplotNames<vector>
1281                 (
1282                     vectorNames,
1283                     ValueLocation_Nodal,
1284                     allVarNames,
1285                     allVarLocation
1286                 );
1289                 tecplotWriter writer(runTime);
1291                 writer.writeInit
1292                 (
1293                     runTime.caseName(),
1294                     allVarNames,
1295                     lagrFileName,
1296                     DataFileType_Full
1297                 );
1299                 writer.writeOrderedZone
1300                 (
1301                     cloudDirs[cloudI],
1302                     strandID++,     //strandID,
1303                     parcels.size(),
1304                     allVarLocation
1305                 );
1307                 // Write coordinates
1308                 writer.writeField(positions.component(0)());
1309                 writer.writeField(positions.component(1)());
1310                 writer.writeField(positions.component(2)());
1312                 // labelFields
1313                 forAll(labelNames, i)
1314                 {
1315                     IOField<label> fld
1316                     (
1317                         IOobject
1318                         (
1319                             labelNames[i],
1320                             mesh.time().timeName(),
1321                             cloud::prefix/cloudDirs[cloudI],
1322                             mesh,
1323                             IOobject::MUST_READ,
1324                             IOobject::NO_WRITE,
1325                             false
1326                         )
1327                     );
1329                     scalarField sfld(fld.size());
1330                     forAll(fld, j)
1331                     {
1332                         sfld[j] = scalar(fld[j]);
1333                     }
1334                     writer.writeField(sfld);
1335                 }
1336                 // scalarFields
1337                 forAll(scalarNames, i)
1338                 {
1339                     IOField<scalar> fld
1340                     (
1341                         IOobject
1342                         (
1343                             scalarNames[i],
1344                             mesh.time().timeName(),
1345                             cloud::prefix/cloudDirs[cloudI],
1346                             mesh,
1347                             IOobject::MUST_READ,
1348                             IOobject::NO_WRITE,
1349                             false
1350                         )
1351                     );
1352                     writer.writeField(fld);
1353                 }
1354                 // vectorFields
1355                 forAll(vectorNames, i)
1356                 {
1357                     IOField<vector> fld
1358                     (
1359                         IOobject
1360                         (
1361                             vectorNames[i],
1362                             mesh.time().timeName(),
1363                             cloud::prefix/cloudDirs[cloudI],
1364                             mesh,
1365                             IOobject::MUST_READ,
1366                             IOobject::NO_WRITE,
1367                             false
1368                         )
1369                     );
1370                     writer.writeField(fld);
1371                 }
1373                 writer.writeEnd();
1374             }
1375         }
1376     }
1378     Info<< "End\n" << endl;
1380     return 0;
1384 // ************************************************************************* //