Merge branch 'master' of ssh://git.code.sf.net/p/foam-extend/foam-extend-3.2
[foam-extend-3.2.git] / applications / utilities / postProcessing / dataConversion / foamToEnsight / ensightMesh.C
blobd215c2da040a4eeb04561a840367ffeb294f7a07
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 \*---------------------------------------------------------------------------*/
26 #include "argList.H"
27 #include "objectRegistry.H"
28 #include "foamTime.H"
29 #include "ensightMesh.H"
30 #include "fvMesh.H"
31 #include "globalMeshData.H"
32 #include "PstreamCombineReduceOps.H"
33 #include "processorPolyPatch.H"
34 #include "cellModeller.H"
35 #include "IOmanip.H"
36 #include "itoa.H"
37 #include "ensightWriteBinary.H"
38 #include <fstream>
40 // * * * * * * * * * * * * * Private Functions * * * * * * * * * * * * * * //
42 namespace Foam
44     //- Proxy-class to hold the patch processor list combination operator
45     class concatPatchProcs
46     {
48     public:
50         void operator()
51         (
52             List<labelList>& x,
53             const List<labelList>& y
54         ) const
55         {
56             forAll(y, i)
57             {
58                 const labelList& yPatches = y[i];
60                 if (yPatches.size())
61                 {
62                     labelList& xPatches = x[i];
64                     label offset = xPatches.size();
65                     xPatches.setSize(offset + yPatches.size());
67                     forAll(yPatches, i)
68                     {
69                         xPatches[i + offset] = yPatches[i];
70                     }
71                 }
72             }
73         }
74     };
75 } // End namespace Foam
78 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
80 Foam::ensightMesh::ensightMesh
82     const fvMesh& mesh,
83     const argList& args,
84     const bool binary
87     mesh_(mesh),
88     binary_(binary),
89     patchPartOffset_(2),
90     meshCellSets_(mesh_.nCells()),
91     boundaryFaceSets_(mesh_.boundary().size()),
92     allPatchNames_(0),
93     allPatchProcs_(0),
94     patchNames_(0),
95     nPatchPrims_(0)
97     const cellShapeList& cellShapes = mesh.cellShapes();
99     const cellModel& tet = *(cellModeller::lookup("tet"));
100     const cellModel& pyr = *(cellModeller::lookup("pyr"));
101     const cellModel& prism = *(cellModeller::lookup("prism"));
102     const cellModel& wedge = *(cellModeller::lookup("wedge"));
103     const cellModel& hex = *(cellModeller::lookup("hex"));
105     if (!args.optionFound("noPatches"))
106     {
107         allPatchNames_ = wordList::subList
108         (
109             mesh_.boundaryMesh().names(), mesh_.boundary().size()
110           - mesh_.globalData().processorPatches().size()
111         );
113         allPatchProcs_.setSize(allPatchNames_.size());
115         forAll (allPatchProcs_, patchi)
116         {
117             if (mesh_.boundary()[patchi].size())
118             {
119                 allPatchProcs_[patchi].setSize(1);
120                 allPatchProcs_[patchi][0] = Pstream::myProcNo();
121             }
122         }
124         combineReduce(allPatchProcs_, concatPatchProcs());
126         if (args.optionFound("patches"))
127         {
128             wordList patchNameList(args.optionLookup("patches")());
130             if (patchNameList.empty())
131             {
132                 patchNameList = allPatchNames_;
133             }
135             forAll (patchNameList, i)
136             {
137                 patchNames_.insert(patchNameList[i]);
138             }
139         }
140     }
142     if (patchNames_.size())
143     {
144         // no internalMesh
145         patchPartOffset_ = 1;
146     }
147     else
148     {
149         // Count the shapes
150         labelList& tets = meshCellSets_.tets;
151         labelList& pyrs = meshCellSets_.pyrs;
152         labelList& prisms = meshCellSets_.prisms;
153         labelList& wedges = meshCellSets_.wedges;
154         labelList& hexes = meshCellSets_.hexes;
155         labelList& polys = meshCellSets_.polys;
157         label nTets = 0;
158         label nPyrs = 0;
159         label nPrisms = 0;
160         label nWedges = 0;
161         label nHexes = 0;
162         label nPolys = 0;
164         forAll(cellShapes, cellI)
165         {
166             const cellShape& cellShape = cellShapes[cellI];
167             const cellModel& cellModel = cellShape.model();
169             if (cellModel == tet)
170             {
171                 tets[nTets++] = cellI;
172             }
173             else if (cellModel == pyr)
174             {
175                 pyrs[nPyrs++] = cellI;
176             }
177             else if (cellModel == prism)
178             {
179                 prisms[nPrisms++] = cellI;
180             }
181             else if (cellModel == wedge)
182             {
183                 wedges[nWedges++] = cellI;
184             }
185             else if (cellModel == hex)
186             {
187                 hexes[nHexes++] = cellI;
188             }
189             else
190             {
191                 polys[nPolys++] = cellI;
192             }
193         }
195         tets.setSize(nTets);
196         pyrs.setSize(nPyrs);
197         prisms.setSize(nPrisms);
198         wedges.setSize(nWedges);
199         hexes.setSize(nHexes);
200         polys.setSize(nPolys);
202         meshCellSets_.nTets = nTets;
203         reduce(meshCellSets_.nTets, sumOp<label>());
205         meshCellSets_.nPyrs = nPyrs;
206         reduce(meshCellSets_.nPyrs, sumOp<label>());
208         meshCellSets_.nPrisms = nPrisms;
209         reduce(meshCellSets_.nPrisms, sumOp<label>());
211         meshCellSets_.nHexesWedges = nHexes + nWedges;
212         reduce(meshCellSets_.nHexesWedges, sumOp<label>());
214         meshCellSets_.nPolys = nPolys;
215         reduce(meshCellSets_.nPolys, sumOp<label>());
216     }
218     if (!args.optionFound("noPatches"))
219     {
220         forAll (mesh.boundary(), patchi)
221         {
222             if (mesh.boundary()[patchi].size())
223             {
224                 const polyPatch& p = mesh.boundaryMesh()[patchi];
226                 labelList& tris = boundaryFaceSets_[patchi].tris;
227                 labelList& quads = boundaryFaceSets_[patchi].quads;
228                 labelList& polys = boundaryFaceSets_[patchi].polys;
230                 tris.setSize(p.size());
231                 quads.setSize(p.size());
232                 polys.setSize(p.size());
234                 label nTris = 0;
235                 label nQuads = 0;
236                 label nPolys = 0;
238                 forAll(p, faceI)
239                 {
240                     const face& f = p[faceI];
242                     if (f.size() == 3)
243                     {
244                         tris[nTris++] = faceI;
245                     }
246                     else if (f.size() == 4)
247                     {
248                         quads[nQuads++] = faceI;
249                     }
250                     else
251                     {
252                         polys[nPolys++] = faceI;
253                     }
254                 }
256                 tris.setSize(nTris);
257                 quads.setSize(nQuads);
258                 polys.setSize(nPolys);
259             }
260         }
261     }
264     forAll(allPatchNames_, patchi)
265     {
266         const word& patchName = allPatchNames_[patchi];
267         nFacePrimitives nfp;
269         if (patchNames_.empty() || patchNames_.found(patchName))
270         {
271             if (mesh.boundary()[patchi].size())
272             {
273                 nfp.nPoints = mesh.boundaryMesh()[patchi].localPoints().size();
274                 nfp.nTris   = boundaryFaceSets_[patchi].tris.size();
275                 nfp.nQuads  = boundaryFaceSets_[patchi].quads.size();
276                 nfp.nPolys  = boundaryFaceSets_[patchi].polys.size();
277             }
278         }
280         reduce(nfp.nPoints, sumOp<label>());
281         reduce(nfp.nTris, sumOp<label>());
282         reduce(nfp.nQuads, sumOp<label>());
283         reduce(nfp.nPolys, sumOp<label>());
285         nPatchPrims_.insert(patchName, nfp);
286     }
290 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
292 Foam::ensightMesh::~ensightMesh()
296 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
298 void Foam::ensightMesh::writePoints
300     const scalarField& pointsComponent,
301     OFstream& ensightGeometryFile
302 ) const
304     forAll(pointsComponent, pointI)
305     {
306         ensightGeometryFile<< setw(12) << float(pointsComponent[pointI]) << nl;
307     }
311 Foam::cellShapeList Foam::ensightMesh::ensMap
313     const cellShapeList& cellShapes,
314     const labelList& prims
315 ) const
317     cellShapeList mcsl(prims.size());
319     forAll(prims, i)
320     {
321         mcsl[i] = cellShapes[prims[i]];
322     }
324     return mcsl;
328 Foam::cellShapeList Foam::ensightMesh::ensMap
330     const cellShapeList& cellShapes,
331     const labelList& hexes,
332     const labelList& wedges
333 ) const
335     cellShapeList mcsl(hexes.size() + wedges.size());
337     forAll(hexes, i)
338     {
339         mcsl[i] = cellShapes[hexes[i]];
340     }
342     label offset = hexes.size();
344     const cellModel& hex = *(cellModeller::lookup("hex"));
345     labelList hexLabels(8);
347     forAll(wedges, i)
348     {
349         const cellShape& cellPoints = cellShapes[wedges[i]];
351         hexLabels[0] = cellPoints[0];
352         hexLabels[1] = cellPoints[1];
353         hexLabels[2] = cellPoints[0];
354         hexLabels[3] = cellPoints[2];
355         hexLabels[4] = cellPoints[3];
356         hexLabels[5] = cellPoints[4];
357         hexLabels[6] = cellPoints[6];
358         hexLabels[7] = cellPoints[5];
360         mcsl[i + offset] = cellShape(hex, hexLabels);
361     }
363     return mcsl;
367 void Foam::ensightMesh::writePrims
369     const cellShapeList& cellShapes,
370     const label pointOffset,
371     OFstream& ensightGeometryFile
372 ) const
374     label po = pointOffset + 1;
376     forAll(cellShapes, i)
377     {
378         const cellShape& cellPoints = cellShapes[i];
380         forAll(cellPoints, pointI)
381         {
382             ensightGeometryFile<< setw(10) << cellPoints[pointI] + po;
383         }
384         ensightGeometryFile << nl;
385     }
389 void Foam::ensightMesh::writePrimsBinary
391     const cellShapeList& cellShapes,
392     const label pointOffset,
393     std::ofstream& ensightGeometryFile
394 ) const
396     label po = pointOffset + 1;
398     if (cellShapes.size())
399     {
400         // All the cellShapes have the same number of elements!
401         int numIntElem = cellShapes.size()*cellShapes[0].size();
402         List<int> temp(numIntElem);
404         int n = 0;
406         forAll(cellShapes, i)
407         {
408             const cellShape& cellPoints = cellShapes[i];
410             forAll(cellPoints, pointI)
411             {
412                 temp[n] = cellPoints[pointI] + po;
413                 n++;
414             }
415         }
417         ensightGeometryFile.write
418         (
419             reinterpret_cast<char*>(temp.begin()),
420             numIntElem*sizeof(int)
421         );
422     }
426 void Foam::ensightMesh::writePolysNFaces
428     const labelList& polys,
429     const cellList& cellFaces,
430     OFstream& ensightGeometryFile
431 ) const
433         forAll(polys, i)
434         {
435             ensightGeometryFile
436                 << setw(10) << cellFaces[polys[i]].size() << nl;
437         }
441 void Foam::ensightMesh::writePolysNPointsPerFace
443     const labelList& polys,
444     const cellList& cellFaces,
445     const faceList& faces,
446     OFstream& ensightGeometryFile
447 ) const
449     forAll(polys, i)
450     {
451         const labelList& cf = cellFaces[polys[i]];
453         forAll(cf, faceI)
454         {
455             ensightGeometryFile
456                 << setw(10) << faces[cf[faceI]].size() << nl;
457         }
458     }
462 void Foam::ensightMesh::writePolysPoints
464     const labelList& polys,
465     const cellList& cellFaces,
466     const faceList& faces,
467     const label pointOffset,
468     OFstream& ensightGeometryFile
469 ) const
471     label po = pointOffset + 1;
473     forAll(polys, i)
474     {
475         const labelList& cf = cellFaces[polys[i]];
477         forAll(cf, faceI)
478         {
479             const face& f = faces[cf[faceI]];
481             forAll(f, pointI)
482             {
483                 ensightGeometryFile << setw(10) << f[pointI] + po;
484             }
485             ensightGeometryFile << nl;
486         }
487     }
491 void Foam::ensightMesh::writeAllPolys
493     const labelList& pointOffsets,
494     OFstream& ensightGeometryFile
495 ) const
497     if (meshCellSets_.nPolys)
498     {
499         const cellList& cellFaces = mesh_.cells();
500         const faceList& faces = mesh_.faces();
502         if (Pstream::master())
503         {
504             ensightGeometryFile
505                 << "nfaced" << nl << setw(10) << meshCellSets_.nPolys << nl;
506         }
508         // Number of faces for each poly cell
509         if (Pstream::master())
510         {
511             // Master
512             writePolysNFaces
513             (
514                 meshCellSets_.polys,
515                 cellFaces,
516                 ensightGeometryFile
517             );
518             // Slaves
519             for (int slave=1; slave<Pstream::nProcs(); slave++)
520             {
521                 IPstream fromSlave(Pstream::scheduled, slave);
522                 labelList polys(fromSlave);
523                 cellList cellFaces(fromSlave);
525                 writePolysNFaces
526                 (
527                     polys,
528                     cellFaces,
529                     ensightGeometryFile
530                 );
531             }
532         }
533         else
534         {
535             OPstream toMaster(Pstream::scheduled, Pstream::masterNo());
536             toMaster<< meshCellSets_.polys << cellFaces;
537         }
539         // Number of points for each face of the above list
540         if (Pstream::master())
541         {
542             // Master
543             writePolysNPointsPerFace
544             (
545                 meshCellSets_.polys,
546                 cellFaces,
547                 faces,
548                 ensightGeometryFile
549             );
550             // Slaves
551             for (int slave=1; slave<Pstream::nProcs(); slave++)
552             {
553                 IPstream fromSlave(Pstream::scheduled, slave);
554                 labelList polys(fromSlave);
555                 cellList cellFaces(fromSlave);
556                 faceList faces(fromSlave);
558                 writePolysNPointsPerFace
559                 (
560                     polys,
561                     cellFaces,
562                     faces,
563                     ensightGeometryFile
564                 );
565             }
566         }
567         else
568         {
569             OPstream toMaster(Pstream::scheduled, Pstream::masterNo());
570             toMaster<< meshCellSets_.polys << cellFaces << faces;
571         }
573         // List of points id for each face of the above list
574         if (Pstream::master())
575         {
576             // Master
577             writePolysPoints
578             (
579                 meshCellSets_.polys,
580                 cellFaces,
581                 faces,
582                 0,
583                 ensightGeometryFile
584             );
585             // Slaves
586             for (int slave=1; slave<Pstream::nProcs(); slave++)
587             {
588                 IPstream fromSlave(Pstream::scheduled, slave);
589                 labelList polys(fromSlave);
590                 cellList cellFaces(fromSlave);
591                 faceList faces(fromSlave);
593                 writePolysPoints
594                 (
595                     polys,
596                     cellFaces,
597                     faces,
598                     pointOffsets[slave-1],
599                     ensightGeometryFile
600                 );
601             }
602         }
603         else
604         {
605             OPstream toMaster(Pstream::scheduled, Pstream::masterNo());
606             toMaster<< meshCellSets_.polys << cellFaces << faces;
607         }
608     }
612 void Foam::ensightMesh::writePolysNFacesBinary
614     const labelList& polys,
615     const cellList& cellFaces,
616     std::ofstream& ensightGeometryFile
617 ) const
619     forAll(polys, i)
620     {
621         writeEnsDataBinary
622         (
623             cellFaces[polys[i]].size(),
624             ensightGeometryFile
625         );
626     }
630 void Foam::ensightMesh::writePolysNPointsPerFaceBinary
632     const labelList& polys,
633     const cellList& cellFaces,
634     const faceList& faces,
635     std::ofstream& ensightGeometryFile
636 ) const
638     forAll(polys, i)
639     {
640         const labelList& cf = cellFaces[polys[i]];
642         forAll(cf, faceI)
643         {
644             writeEnsDataBinary
645             (
646                 faces[cf[faceI]].size(),
647                 ensightGeometryFile
648             );
649         }
650     }
654 void Foam::ensightMesh::writePolysPointsBinary
656     const labelList& polys,
657     const cellList& cellFaces,
658     const faceList& faces,
659     const label pointOffset,
660     std::ofstream& ensightGeometryFile
661 ) const
663     label po = pointOffset + 1;
665     forAll(polys, i)
666     {
667         const labelList& cf = cellFaces[polys[i]];
669         forAll(cf, faceI)
670         {
671             const face& f = faces[cf[faceI]];
673             forAll(f, pointI)
674             {
675                 writeEnsDataBinary(f[pointI] + po,ensightGeometryFile);
676             }
677         }
678     }
682 void Foam::ensightMesh::writeAllPolysBinary
684     const labelList& pointOffsets,
685     std::ofstream& ensightGeometryFile
686 ) const
688     if (meshCellSets_.nPolys)
689     {
690         const cellList& cellFaces = mesh_.cells();
691         const faceList& faces = mesh_.faces();
693         if (Pstream::master())
694         {
695             writeEnsDataBinary("nfaced",ensightGeometryFile);
696             writeEnsDataBinary(meshCellSets_.nPolys,ensightGeometryFile);
697         }
699         // Number of faces for each poly cell
700         if (Pstream::master())
701         {
702             // Master
703             writePolysNFacesBinary
704             (
705                 meshCellSets_.polys,
706                 cellFaces,
707                 ensightGeometryFile
708             );
709             // Slaves
710             for (int slave=1; slave<Pstream::nProcs(); slave++)
711             {
712                 IPstream fromSlave(Pstream::scheduled, slave);
713                 labelList polys(fromSlave);
714                 cellList cellFaces(fromSlave);
716                 writePolysNFacesBinary
717                 (
718                     polys,
719                     cellFaces,
720                     ensightGeometryFile
721                 );
722             }
723         }
724         else
725         {
726             OPstream toMaster(Pstream::scheduled, Pstream::masterNo());
727             toMaster<< meshCellSets_.polys << cellFaces;
728         }
730         // Number of points for each face of the above list
731         if (Pstream::master())
732         {
733             // Master
734             writePolysNPointsPerFaceBinary
735             (
736                 meshCellSets_.polys,
737                 cellFaces,
738                 faces,
739                 ensightGeometryFile
740             );
741             // Slaves
742             for (int slave=1; slave<Pstream::nProcs(); slave++)
743             {
744                 IPstream fromSlave(Pstream::scheduled, slave);
745                 labelList polys(fromSlave);
746                 cellList cellFaces(fromSlave);
747                 faceList faces(fromSlave);
749                 writePolysNPointsPerFaceBinary
750                 (
751                     polys,
752                     cellFaces,
753                     faces,
754                     ensightGeometryFile
755                 );
756             }
757         }
758         else
759         {
760             OPstream toMaster(Pstream::scheduled, Pstream::masterNo());
761             toMaster<< meshCellSets_.polys << cellFaces << faces;
762         }
764         // List of points id for each face of the above list
765         if (Pstream::master())
766         {
767             // Master
768             writePolysPointsBinary
769             (
770                 meshCellSets_.polys,
771                 cellFaces,
772                 faces,
773                 0,
774                 ensightGeometryFile
775             );
776             // Slaves
777             for (int slave=1; slave<Pstream::nProcs(); slave++)
778             {
779                 IPstream fromSlave(Pstream::scheduled, slave);
780                 labelList polys(fromSlave);
781                 cellList cellFaces(fromSlave);
782                 faceList faces(fromSlave);
784                 writePolysPointsBinary
785                 (
786                     polys,
787                     cellFaces,
788                     faces,
789                     pointOffsets[slave-1],
790                     ensightGeometryFile
791                 );
792             }
793         }
794         else
795         {
796             OPstream toMaster(Pstream::scheduled, Pstream::masterNo());
797             toMaster<< meshCellSets_.polys << cellFaces << faces;
798         }
799     }
803 void Foam::ensightMesh::writeAllPrims
805     const char* key,
806     const label nPrims,
807     const cellShapeList& cellShapes,
808     const labelList& pointOffsets,
809     OFstream& ensightGeometryFile
810 ) const
812     if (nPrims)
813     {
814         if (Pstream::master())
815         {
816             ensightGeometryFile << key << nl << setw(10) << nPrims << nl;
818             writePrims(cellShapes, 0, ensightGeometryFile);
820             for (int slave=1; slave<Pstream::nProcs(); slave++)
821             {
822                 IPstream fromSlave(Pstream::scheduled, slave);
823                 cellShapeList cellShapes(fromSlave);
825                 writePrims
826                 (
827                     cellShapes,
828                     pointOffsets[slave-1],
829                     ensightGeometryFile
830                 );
831             }
832         }
833         else
834         {
835             OPstream toMaster(Pstream::scheduled, Pstream::masterNo());
836             toMaster<< cellShapes;
837         }
838     }
842 void Foam::ensightMesh::writeAllPrimsBinary
844     const char* key,
845     const label nPrims,
846     const cellShapeList& cellShapes,
847     const labelList& pointOffsets,
848     std::ofstream& ensightGeometryFile
849 ) const
851     if (nPrims)
852     {
853         if (Pstream::master())
854         {
855             writeEnsDataBinary(key,ensightGeometryFile);
856             writeEnsDataBinary(nPrims,ensightGeometryFile);
858             writePrimsBinary(cellShapes, 0, ensightGeometryFile);
860             for (int slave=1; slave<Pstream::nProcs(); slave++)
861             {
862                 IPstream fromSlave(Pstream::scheduled, slave);
863                 cellShapeList cellShapes(fromSlave);
865                 writePrimsBinary
866                 (
867                     cellShapes,
868                     pointOffsets[slave-1],
869                     ensightGeometryFile
870                 );
871             }
872         }
873         else
874         {
875             OPstream toMaster(Pstream::scheduled, Pstream::masterNo());
876             toMaster<< cellShapes;
877         }
878     }
882 void Foam::ensightMesh::writeFacePrims
884     const faceList& patchFaces,
885     const label pointOffset,
886     OFstream& ensightGeometryFile
887 ) const
889     if (patchFaces.size())
890     {
891         label po = pointOffset + 1;
893         forAll(patchFaces, i)
894         {
895             const face& patchFace = patchFaces[i];
897             forAll(patchFace, pointI)
898             {
899                 ensightGeometryFile << setw(10) << patchFace[pointI] + po;
900             }
901             ensightGeometryFile << nl;
902         }
903     }
907 void Foam::ensightMesh::writeFacePrimsBinary
909     const faceList& patchFaces,
910     const label pointOffset,
911     std::ofstream& ensightGeometryFile
912 ) const
914     if (patchFaces.size())
915     {
916         label po = pointOffset + 1;
918         forAll(patchFaces, i)
919         {
920             const face& patchFace = patchFaces[i];
922             forAll(patchFace, pointI)
923             {
924                 writeEnsDataBinary
925                 (
926                     patchFace[pointI] + po,
927                     ensightGeometryFile
928                 );
929             }
930         }
931     }
935 Foam::faceList Foam::ensightMesh::ensMap
937     const faceList& patchFaces,
938     const labelList& prims
939 ) const
941     faceList ppf(prims.size());
943     forAll (prims, i)
944     {
945         ppf[i] = patchFaces[prims[i]];
946     }
948     return ppf;
952 void Foam::ensightMesh::writeAllFacePrims
954     const char* key,
955     const labelList& prims,
956     const label nPrims,
957     const faceList& patchFaces,
958     const labelList& pointOffsets,
959     const labelList& patchProcessors,
960     OFstream& ensightGeometryFile
961 ) const
963     if (nPrims)
964     {
965         if (Pstream::master())
966         {
967             ensightGeometryFile << key << nl << setw(10) << nPrims << nl;
969             if (&prims != NULL)
970             {
971                 writeFacePrims
972                 (
973                     ensMap(patchFaces, prims),
974                     0,
975                     ensightGeometryFile
976                 );
977             }
979             forAll (patchProcessors, i)
980             {
981                 if (patchProcessors[i] != 0)
982                 {
983                     label slave = patchProcessors[i];
984                     IPstream fromSlave(Pstream::scheduled, slave);
985                     faceList patchFaces(fromSlave);
987                     writeFacePrims
988                     (
989                         patchFaces,
990                         pointOffsets[i],
991                         ensightGeometryFile
992                     );
993                 }
994             }
995         }
996         else if (&prims != NULL)
997         {
998             OPstream toMaster(Pstream::scheduled, Pstream::masterNo());
999             toMaster<< ensMap(patchFaces, prims);
1000         }
1001     }
1005 void Foam::ensightMesh::writeNSidedNPointsPerFace
1007     const faceList& patchFaces,
1008     OFstream& ensightGeometryFile
1009 ) const
1011     forAll(patchFaces, i)
1012     {
1013         ensightGeometryFile
1014             << setw(10) << patchFaces[i].size() << nl;
1015     }
1019 void Foam::ensightMesh::writeNSidedPoints
1021     const faceList& patchFaces,
1022     const label pointOffset,
1023     OFstream& ensightGeometryFile
1024 ) const
1026     writeFacePrims
1027     (
1028         patchFaces,
1029         pointOffset,
1030         ensightGeometryFile
1031     );
1035 void Foam::ensightMesh::writeAllNSided
1037     const labelList& prims,
1038     const label nPrims,
1039     const faceList& patchFaces,
1040     const labelList& pointOffsets,
1041     const labelList& patchProcessors,
1042     OFstream& ensightGeometryFile
1043 ) const
1045     if (nPrims)
1046     {
1047         if (Pstream::master())
1048         {
1049             ensightGeometryFile
1050                 << "nsided" << nl << setw(10) << nPrims << nl;
1051         }
1053         // Number of points for each face
1054         if (Pstream::master())
1055         {
1056             if (&prims != NULL)
1057             {
1058                 writeNSidedNPointsPerFace
1059                 (
1060                     ensMap(patchFaces, prims),
1061                     ensightGeometryFile
1062                 );
1063             }
1065             forAll (patchProcessors, i)
1066             {
1067                 if (patchProcessors[i] != 0)
1068                 {
1069                     label slave = patchProcessors[i];
1070                     IPstream fromSlave(Pstream::scheduled, slave);
1071                     faceList patchFaces(fromSlave);
1073                     writeNSidedNPointsPerFace
1074                     (
1075                         patchFaces,
1076                         ensightGeometryFile
1077                     );
1078                 }
1079             }
1080         }
1081         else if (&prims != NULL)
1082         {
1083             OPstream toMaster(Pstream::scheduled, Pstream::masterNo());
1084             toMaster<< ensMap(patchFaces, prims);
1085         }
1087         // List of points id for each face
1088         if (Pstream::master())
1089         {
1090             if (&prims != NULL)
1091             {
1092                 writeNSidedPoints
1093                 (
1094                     ensMap(patchFaces, prims),
1095                     0,
1096                     ensightGeometryFile
1097                 );
1098             }
1100             forAll (patchProcessors, i)
1101             {
1102                 if (patchProcessors[i] != 0)
1103                 {
1104                     label slave = patchProcessors[i];
1105                     IPstream fromSlave(Pstream::scheduled, slave);
1106                     faceList patchFaces(fromSlave);
1108                     writeNSidedPoints
1109                     (
1110                         patchFaces,
1111                         pointOffsets[i],
1112                         ensightGeometryFile
1113                     );
1114                 }
1115             }
1116         }
1117         else if (&prims != NULL)
1118         {
1119             OPstream toMaster(Pstream::scheduled, Pstream::masterNo());
1120             toMaster<< ensMap(patchFaces, prims);
1121         }
1122     }
1126 void Foam::ensightMesh::writeNSidedPointsBinary
1128     const faceList& patchFaces,
1129     const label pointOffset,
1130     std::ofstream& ensightGeometryFile
1131 ) const
1133     writeFacePrimsBinary
1134     (
1135         patchFaces,
1136         pointOffset,
1137         ensightGeometryFile
1138     );
1142 void Foam::ensightMesh::writeNSidedNPointsPerFaceBinary
1144     const faceList& patchFaces,
1145     std::ofstream& ensightGeometryFile
1146 ) const
1148     forAll(patchFaces, i)
1149     {
1150         writeEnsDataBinary
1151         (
1152             patchFaces[i].size(),
1153             ensightGeometryFile
1154         );
1155     }
1159 void Foam::ensightMesh::writeAllNSidedBinary
1161     const labelList& prims,
1162     const label nPrims,
1163     const faceList& patchFaces,
1164     const labelList& pointOffsets,
1165     const labelList& patchProcessors,
1166     std::ofstream& ensightGeometryFile
1167 ) const
1169     if (nPrims)
1170     {
1171         if (Pstream::master())
1172         {
1173             writeEnsDataBinary("nsided",ensightGeometryFile);
1174             writeEnsDataBinary(nPrims,ensightGeometryFile);
1175         }
1177         // Number of points for each face
1178         if (Pstream::master())
1179         {
1180             if (&prims != NULL)
1181             {
1182                 writeNSidedNPointsPerFaceBinary
1183                 (
1184                     ensMap(patchFaces, prims),
1185                     ensightGeometryFile
1186                 );
1187             }
1189             forAll (patchProcessors, i)
1190             {
1191                 if (patchProcessors[i] != 0)
1192                 {
1193                     label slave = patchProcessors[i];
1194                     IPstream fromSlave(Pstream::scheduled, slave);
1195                     faceList patchFaces(fromSlave);
1197                     writeNSidedNPointsPerFaceBinary
1198                     (
1199                         patchFaces,
1200                         ensightGeometryFile
1201                     );
1202                 }
1203             }
1204         }
1205         else if (&prims != NULL)
1206         {
1207             OPstream toMaster(Pstream::scheduled, Pstream::masterNo());
1208             toMaster<< ensMap(patchFaces, prims);
1209         }
1211         // List of points id for each face
1212         if (Pstream::master())
1213         {
1214             if (&prims != NULL)
1215             {
1216                 writeNSidedPointsBinary
1217                 (
1218                     ensMap(patchFaces, prims),
1219                     0,
1220                     ensightGeometryFile
1221                 );
1222             }
1224             forAll (patchProcessors, i)
1225             {
1226                 if (patchProcessors[i] != 0)
1227                 {
1228                     label slave = patchProcessors[i];
1229                     IPstream fromSlave(Pstream::scheduled, slave);
1230                     faceList patchFaces(fromSlave);
1232                     writeNSidedPointsBinary
1233                     (
1234                         patchFaces,
1235                         pointOffsets[i],
1236                         ensightGeometryFile
1237                     );
1238                 }
1239             }
1240         }
1241         else if (&prims != NULL)
1242         {
1243             OPstream toMaster(Pstream::scheduled, Pstream::masterNo());
1244             toMaster<< ensMap(patchFaces, prims);
1245         }
1246     }
1250 void Foam::ensightMesh::writeAllFacePrimsBinary
1252     const char* key,
1253     const labelList& prims,
1254     const label nPrims,
1255     const faceList& patchFaces,
1256     const labelList& pointOffsets,
1257     const labelList& patchProcessors,
1258     std::ofstream& ensightGeometryFile
1259 ) const
1261     if (nPrims)
1262     {
1263         if (Pstream::master())
1264         {
1265             writeEnsDataBinary(key,ensightGeometryFile);
1266             writeEnsDataBinary(nPrims,ensightGeometryFile);
1268             if (&prims != NULL)
1269             {
1270                 writeFacePrimsBinary
1271                 (
1272                     ensMap(patchFaces, prims),
1273                     0,
1274                     ensightGeometryFile
1275                 );
1276             }
1278             forAll (patchProcessors, i)
1279             {
1280                 if (patchProcessors[i] != 0)
1281                 {
1282                     label slave = patchProcessors[i];
1283                     IPstream fromSlave(Pstream::scheduled, slave);
1284                     faceList patchFaces(fromSlave);
1286                     writeFacePrimsBinary
1287                     (
1288                         patchFaces,
1289                         pointOffsets[i],
1290                         ensightGeometryFile
1291                     );
1292                 }
1293             }
1294         }
1295         else if (&prims != NULL)
1296         {
1297             OPstream toMaster(Pstream::scheduled, Pstream::masterNo());
1298             toMaster<< ensMap(patchFaces, prims);
1299         }
1300     }
1304 void Foam::ensightMesh::write
1306     const fileName& postProcPath,
1307     const word& prepend,
1308     const label timeIndex,
1309     Ostream& ensightCaseFile
1310 ) const
1312     if (binary_)
1313     {
1314         writeBinary(postProcPath, prepend, timeIndex, ensightCaseFile);
1315     }
1316     else
1317     {
1318         writeAscii(postProcPath, prepend, timeIndex, ensightCaseFile);
1319     }
1323 void Foam::ensightMesh::writeAscii
1325     const fileName& postProcPath,
1326     const word& prepend,
1327     const label timeIndex,
1328     Ostream& ensightCaseFile
1329 ) const
1331     const Time& runTime = mesh_.time();
1332     const pointField& points = mesh_.points();
1333     const cellShapeList& cellShapes = mesh_.cellShapes();
1335     word timeFile = prepend;
1337     if (timeIndex == 0)
1338     {
1339         timeFile += "000.";
1340     }
1341     else if (mesh_.moving())
1342     {
1343         timeFile += itoa(timeIndex) + '.';
1344     }
1346     // set the filename of the ensight file
1347     fileName ensightGeometryFileName = timeFile + "mesh";
1349     OFstream *ensightGeometryFilePtr = NULL;
1350     if (Pstream::master())
1351     {
1352         ensightGeometryFilePtr = new OFstream
1353         (
1354             postProcPath/ensightGeometryFileName,
1355             ios_base::out|ios_base::trunc,
1356             runTime.writeFormat(),
1357             runTime.writeVersion(),
1358             IOstream::UNCOMPRESSED
1359         );
1360     }
1362     OFstream& ensightGeometryFile = *ensightGeometryFilePtr;
1364     if (Pstream::master())
1365     {
1366         // Set Format
1367         ensightGeometryFile.setf
1368         (
1369             ios_base::scientific,
1370             ios_base::floatfield
1371         );
1372         ensightGeometryFile.precision(5);
1374         ensightGeometryFile
1375             << "EnSight Geometry File" << nl
1376             << "written from FOAM-" << Foam::FOAMversion << nl
1377             << "node id assign" << nl
1378             << "element id assign" << nl;
1379     }
1381     labelList pointOffsets(Pstream::nProcs(), 0);
1383     if (patchNames_.empty())
1384     {
1385         label nPoints = points.size();
1386         Pstream::gather(nPoints, sumOp<label>());
1388         if (Pstream::master())
1389         {
1390             ensightGeometryFile
1391                 << "part" << nl
1392                 << setw(10) << 1 << nl
1393                 << "internalMesh" << nl
1394                 << "coordinates" << nl
1395                 << setw(10) << nPoints
1396                 << endl;
1398             for (direction d=0; d<vector::nComponents; d++)
1399             {
1400                 writePoints(points.component(d), ensightGeometryFile);
1401                 pointOffsets[0] = points.size();
1403                 for (int slave=1; slave<Pstream::nProcs(); slave++)
1404                 {
1405                     IPstream fromSlave(Pstream::scheduled, slave);
1406                     scalarField pointsComponent(fromSlave);
1407                     writePoints(pointsComponent, ensightGeometryFile);
1408                     pointOffsets[slave] =
1409                         pointOffsets[slave-1]
1410                       + pointsComponent.size();
1411                 }
1412             }
1413         }
1414         else
1415         {
1416             for (direction d=0; d<vector::nComponents; d++)
1417             {
1418                 OPstream toMaster(Pstream::scheduled, Pstream::masterNo());
1419                 toMaster<< points.component(d);
1420             }
1421         }
1423         writeAllPrims
1424         (
1425             "hexa8",
1426             meshCellSets_.nHexesWedges,
1427             ensMap(cellShapes, meshCellSets_.hexes, meshCellSets_.wedges),
1428             pointOffsets,
1429             ensightGeometryFile
1430         );
1432         writeAllPrims
1433         (
1434             "penta6",
1435             meshCellSets_.nPrisms,
1436             ensMap(cellShapes, meshCellSets_.prisms),
1437             pointOffsets,
1438             ensightGeometryFile
1439         );
1441         writeAllPrims
1442         (
1443             "pyramid5",
1444             meshCellSets_.nPyrs,
1445             ensMap(cellShapes, meshCellSets_.pyrs),
1446             pointOffsets,
1447             ensightGeometryFile
1448         );
1450         writeAllPrims
1451         (
1452             "tetra4",
1453             meshCellSets_.nTets,
1454             ensMap(cellShapes, meshCellSets_.tets),
1455             pointOffsets,
1456             ensightGeometryFile
1457         );
1459         writeAllPolys
1460         (
1461             pointOffsets,
1462             ensightGeometryFile
1463         );
1464     }
1467     label ensightPatchI = patchPartOffset_;
1469     forAll(allPatchNames_, patchi)
1470     {
1471         const word& patchName = allPatchNames_[patchi];
1472         const labelList& patchProcessors = allPatchProcs_[patchi];
1474         if (patchNames_.empty() || patchNames_.found(patchName))
1475         {
1476             const nFacePrimitives& nfp = nPatchPrims_.find(patchName)();
1478             const labelList *trisPtr  = NULL;
1479             const labelList *quadsPtr = NULL;
1480             const labelList *polysPtr = NULL;
1482             const pointField *patchPointsPtr = NULL;
1483             const faceList *patchFacesPtr = NULL;
1485             if (mesh_.boundary()[patchi].size())
1486             {
1487                 const polyPatch& p = mesh_.boundaryMesh()[patchi];
1489                 trisPtr  = &boundaryFaceSets_[patchi].tris;
1490                 quadsPtr = &boundaryFaceSets_[patchi].quads;
1491                 polysPtr = &boundaryFaceSets_[patchi].polys;
1493                 patchPointsPtr = &(p.localPoints());
1494                 patchFacesPtr  = &(p.localFaces());
1495             }
1497             const labelList& tris = *trisPtr;
1498             const labelList& quads = *quadsPtr;
1499             const labelList& polys = *polysPtr;
1500             const pointField& patchPoints = *patchPointsPtr;
1501             const faceList& patchFaces = *patchFacesPtr;
1503             if (nfp.nTris || nfp.nQuads || nfp.nPolys)
1504             {
1505                 labelList patchPointOffsets(Pstream::nProcs(), 0);
1507                 if (Pstream::master())
1508                 {
1509                     ensightGeometryFile
1510                         << "part" << nl
1511                         << setw(10) << ensightPatchI++ << nl
1512                         << patchName << nl
1513                         << "coordinates" << nl
1514                         << setw(10) << nfp.nPoints
1515                         << endl;
1517                     for (direction d=0; d<vector::nComponents; d++)
1518                     {
1519                         if (patchPointsPtr)
1520                         {
1521                             writePoints
1522                             (
1523                                 patchPoints.component(d),
1524                                 ensightGeometryFile
1525                             );
1526                         }
1528                         patchPointOffsets = 0;
1530                         forAll (patchProcessors, i)
1531                         {
1532                             if (patchProcessors[i] != 0)
1533                             {
1534                                 label slave = patchProcessors[i];
1535                                 IPstream fromSlave(Pstream::scheduled, slave);
1536                                 scalarField patchPointsComponent(fromSlave);
1538                                 writePoints
1539                                 (
1540                                     patchPointsComponent,
1541                                     ensightGeometryFile
1542                                 );
1544                                 if (i < Pstream::nProcs()-1)
1545                                 {
1546                                     patchPointOffsets[i+1] =
1547                                         patchPointOffsets[i]
1548                                       + patchPointsComponent.size();
1549                                 }
1550                             }
1551                             else
1552                             {
1553                                 if (i < Pstream::nProcs()-1)
1554                                 {
1555                                     patchPointOffsets[i+1] =
1556                                         patchPointOffsets[i]
1557                                       + patchPoints.size();
1558                                 }
1559                             }
1560                         }
1561                     }
1562                 }
1563                 else if (patchPointsPtr)
1564                 {
1565                     for (direction d=0; d<vector::nComponents; d++)
1566                     {
1567                         OPstream toMaster
1568                         (
1569                             Pstream::scheduled,
1570                             Pstream::masterNo()
1571                         );
1572                         toMaster<< patchPoints.component(d);
1573                     }
1574                 }
1576                 writeAllFacePrims
1577                 (
1578                     "tria3",
1579                     tris,
1580                     nfp.nTris,
1581                     patchFaces,
1582                     patchPointOffsets,
1583                     patchProcessors,
1584                     ensightGeometryFile
1585                 );
1587                 writeAllFacePrims
1588                 (
1589                     "quad4",
1590                     quads,
1591                     nfp.nQuads,
1592                     patchFaces,
1593                     patchPointOffsets,
1594                     patchProcessors,
1595                     ensightGeometryFile
1596                 );
1598                 writeAllNSided
1599                 (
1600                     polys,
1601                     nfp.nPolys,
1602                     patchFaces,
1603                     patchPointOffsets,
1604                     patchProcessors,
1605                     ensightGeometryFile
1606                 );
1607             }
1608         }
1609     }
1611     if (Pstream::master())
1612     {
1613         delete ensightGeometryFilePtr;
1614     }
1618 void Foam::ensightMesh::writeBinary
1620     const fileName& postProcPath,
1621     const word& prepend,
1622     const label timeIndex,
1623     Ostream& ensightCaseFile
1624 ) const
1626     //const Time& runTime = mesh.time();
1627     const pointField& points = mesh_.points();
1628     const cellShapeList& cellShapes = mesh_.cellShapes();
1630     word timeFile = prepend;
1632     if (timeIndex == 0)
1633     {
1634         timeFile += "000.";
1635     }
1636     else if (mesh_.moving())
1637     {
1638         timeFile += itoa(timeIndex) + '.';
1639     }
1641     // set the filename of the ensight file
1642     fileName ensightGeometryFileName = timeFile + "mesh";
1644     std::ofstream *ensightGeometryFilePtr = NULL;
1646     if (Pstream::master())
1647     {
1648         ensightGeometryFilePtr = new std::ofstream
1649         (
1650             (postProcPath/ensightGeometryFileName).c_str(),
1651             ios_base::out | ios_base::binary | ios_base::trunc
1652         );
1653         // Check on file opened?
1654     }
1656     std::ofstream& ensightGeometryFile = *ensightGeometryFilePtr;
1658     if (Pstream::master())
1659     {
1660         writeEnsDataBinary("C binary", ensightGeometryFile);
1661         writeEnsDataBinary("EnSight Geometry File", ensightGeometryFile);
1662         writeEnsDataBinary("written from FOAM", ensightGeometryFile);
1663         writeEnsDataBinary("node id assign", ensightGeometryFile);
1664         writeEnsDataBinary("element id assign", ensightGeometryFile);
1665     }
1667     labelList pointOffsets(Pstream::nProcs(), 0);
1669     if (patchNames_.empty())
1670     {
1671         label nPoints = points.size();
1672         Pstream::gather(nPoints, sumOp<label>());
1674         if (Pstream::master())
1675         {
1676             writeEnsDataBinary("part",ensightGeometryFile);
1677             writeEnsDataBinary(1,ensightGeometryFile);
1678             writeEnsDataBinary("internalMesh",ensightGeometryFile);
1679             writeEnsDataBinary("coordinates",ensightGeometryFile);
1680             writeEnsDataBinary(nPoints,ensightGeometryFile);
1682             for (direction d=0; d<vector::nComponents; d++)
1683             {
1684                 //writePointsBinary(points.component(d), ensightGeometryFile);
1685                 writeEnsDataBinary(points.component(d), ensightGeometryFile);
1686                 pointOffsets[0] = points.size();
1688                 for (int slave=1; slave<Pstream::nProcs(); slave++)
1689                 {
1690                     IPstream fromSlave(Pstream::scheduled, slave);
1691                     scalarField pointsComponent(fromSlave);
1692                     //writePointsBinary(pointsComponent, ensightGeometryFile);
1693                     writeEnsDataBinary(pointsComponent, ensightGeometryFile);
1694                     pointOffsets[slave] =
1695                         pointOffsets[slave-1]
1696                       + pointsComponent.size();
1697                 }
1698             }
1699         }
1700         else
1701         {
1702             for (direction d=0; d<vector::nComponents; d++)
1703             {
1704                 OPstream toMaster(Pstream::scheduled, Pstream::masterNo());
1705                 toMaster<< points.component(d);
1706             }
1707         }
1709         writeAllPrimsBinary
1710         (
1711             "hexa8",
1712             meshCellSets_.nHexesWedges,
1713             ensMap(cellShapes, meshCellSets_.hexes, meshCellSets_.wedges),
1714             pointOffsets,
1715             ensightGeometryFile
1716         );
1718         writeAllPrimsBinary
1719         (
1720             "penta6",
1721             meshCellSets_.nPrisms,
1722             ensMap(cellShapes, meshCellSets_.prisms),
1723             pointOffsets,
1724             ensightGeometryFile
1725         );
1727         writeAllPrimsBinary
1728         (
1729             "pyramid5",
1730             meshCellSets_.nPyrs,
1731             ensMap(cellShapes, meshCellSets_.pyrs),
1732             pointOffsets,
1733             ensightGeometryFile
1734         );
1736         writeAllPrimsBinary
1737         (
1738             "tetra4",
1739             meshCellSets_.nTets,
1740             ensMap(cellShapes, meshCellSets_.tets),
1741             pointOffsets,
1742             ensightGeometryFile
1743         );
1745         writeAllPolysBinary
1746         (
1747             pointOffsets,
1748             ensightGeometryFile
1749         );
1751     }
1753     label ensightPatchI = patchPartOffset_;
1754     label iCount = 0;
1756     forAll(allPatchNames_, patchi)
1757     {
1758         iCount ++;
1759         const word& patchName = allPatchNames_[patchi];
1760         const labelList& patchProcessors = allPatchProcs_[patchi];
1762         if (patchNames_.empty() || patchNames_.found(patchName))
1763         {
1764             const nFacePrimitives& nfp = nPatchPrims_.find(patchName)();
1766             const labelList *trisPtr = NULL;
1767             const labelList *quadsPtr = NULL;
1768             const labelList *polysPtr = NULL;
1770             const pointField *patchPointsPtr = NULL;
1771             const faceList *patchFacesPtr = NULL;
1773             if (mesh_.boundary()[patchi].size())
1774             {
1775                 const polyPatch& p = mesh_.boundaryMesh()[patchi];
1777                 trisPtr = &boundaryFaceSets_[patchi].tris;
1778                 quadsPtr = &boundaryFaceSets_[patchi].quads;
1779                 polysPtr = &boundaryFaceSets_[patchi].polys;
1781                 patchPointsPtr = &(p.localPoints());
1782                 patchFacesPtr = &(p.localFaces());
1783             }
1785             const labelList& tris = *trisPtr;
1786             const labelList& quads = *quadsPtr;
1787             const labelList& polys = *polysPtr;
1788             const pointField& patchPoints = *patchPointsPtr;
1789             const faceList& patchFaces = *patchFacesPtr;
1791             if (nfp.nTris || nfp.nQuads || nfp.nPolys)
1792             {
1793                 labelList patchPointOffsets(Pstream::nProcs(), 0);
1795                 if (Pstream::master())
1796                 {
1797                     writeEnsDataBinary("part",ensightGeometryFile);
1798                     writeEnsDataBinary(ensightPatchI++,ensightGeometryFile);
1799                     //writeEnsDataBinary(patchName.c_str(),ensightGeometryFile);
1800                     writeEnsDataBinary(patchName.c_str(),ensightGeometryFile);
1801                     writeEnsDataBinary("coordinates",ensightGeometryFile);
1802                     writeEnsDataBinary(nfp.nPoints,ensightGeometryFile);
1804                     for (direction d=0; d<vector::nComponents; d++)
1805                     {
1806                         if (patchPointsPtr)
1807                         {
1808                             //writePointsBinary
1809                             writeEnsDataBinary
1810                             (
1811                                 patchPoints.component(d),
1812                                 ensightGeometryFile
1813                             );
1814                         }
1816                         patchPointOffsets = 0;
1819                         forAll (patchProcessors, i)
1820                         {
1821                             if (patchProcessors[i] != 0)
1822                             {
1823                                 label slave = patchProcessors[i];
1824                                 IPstream fromSlave(Pstream::scheduled, slave);
1825                                 scalarField patchPointsComponent(fromSlave);
1827                                 //writePointsBinary
1828                                 writeEnsDataBinary
1829                                 (
1830                                     patchPointsComponent,
1831                                     ensightGeometryFile
1832                                 );
1834                                 if (i < Pstream::nProcs()-1)
1835                                 {
1836                                     patchPointOffsets[i+1] =
1837                                         patchPointOffsets[i]
1838                                       + patchPointsComponent.size();
1839                                 }
1840                             }
1841                             else
1842                             {
1843                                 if (i < Pstream::nProcs()-1)
1844                                 {
1845                                     patchPointOffsets[i+1] =
1846                                         patchPointOffsets[i]
1847                                       + patchPoints.size();
1848                                 }
1849                             }
1850                         }
1851                     }
1852                 }
1853                 else if (patchPointsPtr)
1854                 {
1855                     for (direction d=0; d<vector::nComponents; d++)
1856                     {
1857                         OPstream toMaster
1858                         (
1859                             Pstream::scheduled,
1860                             Pstream::masterNo()
1861                         );
1862                         toMaster<< patchPoints.component(d);
1863                     }
1864                 }
1866                 writeAllFacePrimsBinary
1867                 (
1868                     "tria3",
1869                     tris,
1870                     nfp.nTris,
1871                     patchFaces,
1872                     patchPointOffsets,
1873                     patchProcessors,
1874                     ensightGeometryFile
1875                 );
1877                 writeAllFacePrimsBinary
1878                 (
1879                     "quad4",
1880                     quads,
1881                     nfp.nQuads,
1882                     patchFaces,
1883                     patchPointOffsets,
1884                     patchProcessors,
1885                     ensightGeometryFile
1886                 );
1888                 writeAllNSidedBinary
1889                 (
1890                     polys,
1891                     nfp.nPolys,
1892                     patchFaces,
1893                     patchPointOffsets,
1894                     patchProcessors,
1895                     ensightGeometryFile
1896                 );
1897             }
1898         }
1899     }
1902     if (Pstream::master())
1903     {
1904         delete ensightGeometryFilePtr;
1905     }
1909 // ************************************************************************* //