BUGFIX: Uninitialised member variables
[foam-extend-3.2.git] / applications / utilities / postProcessing / dataConversion / foamToEnsight / ensightMesh.C
bloba79b3bfc530fe9ee208e8ae4222537b5c722a9a9
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright held by original author
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 \*---------------------------------------------------------------------------*/
27 #include "argList.H"
28 #include "objectRegistry.H"
29 #include "Time.H"
30 #include "ensightMesh.H"
31 #include "fvMesh.H"
32 #include "globalMeshData.H"
33 #include "PstreamCombineReduceOps.H"
34 #include "processorPolyPatch.H"
35 #include "cellModeller.H"
36 #include "IOmanip.H"
37 #include "itoa.H"
38 #include "ensightWriteBinary.H"
39 #include <fstream>
41 // * * * * * * * * * * * * * Private Functions * * * * * * * * * * * * * * //
43 namespace Foam
45     //- Proxy-class to hold the patch processor list combination operator
46     class concatPatchProcs
47     {
49     public:
51         void operator()
52         (
53             List<labelList>& x,
54             const List<labelList>& y
55         ) const
56         {
57             forAll(y, i)
58             {
59                 const labelList& yPatches = y[i];
61                 if (yPatches.size())
62                 {
63                     labelList& xPatches = x[i];
65                     label offset = xPatches.size();
66                     xPatches.setSize(offset + yPatches.size());
68                     forAll(yPatches, i)
69                     {
70                         xPatches[i + offset] = yPatches[i];
71                     }
72                 }
73             }
74         }
75     };
76 } // End namespace Foam
79 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
81 Foam::ensightMesh::ensightMesh
83     const fvMesh& mesh,
84     const argList& args,
85     const bool binary
88     mesh_(mesh),
89     binary_(binary),
90     patchPartOffset_(2),
91     meshCellSets_(mesh_.nCells()),
92     boundaryFaceSets_(mesh_.boundary().size()),
93     allPatchNames_(0),
94     allPatchProcs_(0),
95     patchNames_(0),
96     nPatchPrims_(0)
98     const cellShapeList& cellShapes = mesh.cellShapes();
100     const cellModel& tet = *(cellModeller::lookup("tet"));
101     const cellModel& pyr = *(cellModeller::lookup("pyr"));
102     const cellModel& prism = *(cellModeller::lookup("prism"));
103     const cellModel& wedge = *(cellModeller::lookup("wedge"));
104     const cellModel& hex = *(cellModeller::lookup("hex"));
106     if (!args.optionFound("noPatches"))
107     {
108         allPatchNames_ = wordList::subList
109         (
110             mesh_.boundaryMesh().names(), mesh_.boundary().size()
111           - mesh_.globalData().processorPatches().size()
112         );
114         allPatchProcs_.setSize(allPatchNames_.size());
116         forAll (allPatchProcs_, patchi)
117         {
118             if (mesh_.boundary()[patchi].size())
119             {
120                 allPatchProcs_[patchi].setSize(1);
121                 allPatchProcs_[patchi][0] = Pstream::myProcNo();
122             }
123         }
125         combineReduce(allPatchProcs_, concatPatchProcs());
127         if (args.optionFound("patches"))
128         {
129             wordList patchNameList(args.optionLookup("patches")());
131             if (patchNameList.empty())
132             {
133                 patchNameList = allPatchNames_;
134             }
136             forAll (patchNameList, i)
137             {
138                 patchNames_.insert(patchNameList[i]);
139             }
140         }
141     }
143     if (patchNames_.size())
144     {
145         // no internalMesh
146         patchPartOffset_ = 1;
147     }
148     else
149     {
150         // Count the shapes
151         labelList& tets = meshCellSets_.tets;
152         labelList& pyrs = meshCellSets_.pyrs;
153         labelList& prisms = meshCellSets_.prisms;
154         labelList& wedges = meshCellSets_.wedges;
155         labelList& hexes = meshCellSets_.hexes;
156         labelList& polys = meshCellSets_.polys;
158         label nTets = 0;
159         label nPyrs = 0;
160         label nPrisms = 0;
161         label nWedges = 0;
162         label nHexes = 0;
163         label nPolys = 0;
165         forAll(cellShapes, cellI)
166         {
167             const cellShape& cellShape = cellShapes[cellI];
168             const cellModel& cellModel = cellShape.model();
170             if (cellModel == tet)
171             {
172                 tets[nTets++] = cellI;
173             }
174             else if (cellModel == pyr)
175             {
176                 pyrs[nPyrs++] = cellI;
177             }
178             else if (cellModel == prism)
179             {
180                 prisms[nPrisms++] = cellI;
181             }
182             else if (cellModel == wedge)
183             {
184                 wedges[nWedges++] = cellI;
185             }
186             else if (cellModel == hex)
187             {
188                 hexes[nHexes++] = cellI;
189             }
190             else
191             {
192                 polys[nPolys++] = cellI;
193             }
194         }
196         tets.setSize(nTets);
197         pyrs.setSize(nPyrs);
198         prisms.setSize(nPrisms);
199         wedges.setSize(nWedges);
200         hexes.setSize(nHexes);
201         polys.setSize(nPolys);
203         meshCellSets_.nTets = nTets;
204         reduce(meshCellSets_.nTets, sumOp<label>());
206         meshCellSets_.nPyrs = nPyrs;
207         reduce(meshCellSets_.nPyrs, sumOp<label>());
209         meshCellSets_.nPrisms = nPrisms;
210         reduce(meshCellSets_.nPrisms, sumOp<label>());
212         meshCellSets_.nHexesWedges = nHexes + nWedges;
213         reduce(meshCellSets_.nHexesWedges, sumOp<label>());
215         meshCellSets_.nPolys = nPolys;
216         reduce(meshCellSets_.nPolys, sumOp<label>());
217     }
219     if (!args.optionFound("noPatches"))
220     {
221         forAll (mesh.boundary(), patchi)
222         {
223             if (mesh.boundary()[patchi].size())
224             {
225                 const polyPatch& p = mesh.boundaryMesh()[patchi];
227                 labelList& tris = boundaryFaceSets_[patchi].tris;
228                 labelList& quads = boundaryFaceSets_[patchi].quads;
229                 labelList& polys = boundaryFaceSets_[patchi].polys;
231                 tris.setSize(p.size());
232                 quads.setSize(p.size());
233                 polys.setSize(p.size());
235                 label nTris = 0;
236                 label nQuads = 0;
237                 label nPolys = 0;
239                 forAll(p, faceI)
240                 {
241                     const face& f = p[faceI];
243                     if (f.size() == 3)
244                     {
245                         tris[nTris++] = faceI;
246                     }
247                     else if (f.size() == 4)
248                     {
249                         quads[nQuads++] = faceI;
250                     }
251                     else
252                     {
253                         polys[nPolys++] = faceI;
254                     }
255                 }
257                 tris.setSize(nTris);
258                 quads.setSize(nQuads);
259                 polys.setSize(nPolys);
260             }
261         }
262     }
265     forAll(allPatchNames_, patchi)
266     {
267         const word& patchName = allPatchNames_[patchi];
268         nFacePrimitives nfp;
270         if (patchNames_.empty() || patchNames_.found(patchName))
271         {
272             if (mesh.boundary()[patchi].size())
273             {
274                 nfp.nPoints = mesh.boundaryMesh()[patchi].localPoints().size();
275                 nfp.nTris   = boundaryFaceSets_[patchi].tris.size();
276                 nfp.nQuads  = boundaryFaceSets_[patchi].quads.size();
277                 nfp.nPolys  = boundaryFaceSets_[patchi].polys.size();
278             }
279         }
281         reduce(nfp.nPoints, sumOp<label>());
282         reduce(nfp.nTris, sumOp<label>());
283         reduce(nfp.nQuads, sumOp<label>());
284         reduce(nfp.nPolys, sumOp<label>());
286         nPatchPrims_.insert(patchName, nfp);
287     }
291 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
293 Foam::ensightMesh::~ensightMesh()
297 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
299 void Foam::ensightMesh::writePoints
301     const scalarField& pointsComponent,
302     OFstream& ensightGeometryFile
303 ) const
305     forAll(pointsComponent, pointI)
306     {
307         ensightGeometryFile<< setw(12) << float(pointsComponent[pointI]) << nl;
308     }
312 Foam::cellShapeList Foam::ensightMesh::map
314     const cellShapeList& cellShapes,
315     const labelList& prims
316 ) const
318     cellShapeList mcsl(prims.size());
320     forAll(prims, i)
321     {
322         mcsl[i] = cellShapes[prims[i]];
323     }
325     return mcsl;
329 Foam::cellShapeList Foam::ensightMesh::map
331     const cellShapeList& cellShapes,
332     const labelList& hexes,
333     const labelList& wedges
334 ) const
336     cellShapeList mcsl(hexes.size() + wedges.size());
338     forAll(hexes, i)
339     {
340         mcsl[i] = cellShapes[hexes[i]];
341     }
343     label offset = hexes.size();
345     const cellModel& hex = *(cellModeller::lookup("hex"));
346     labelList hexLabels(8);
348     forAll(wedges, i)
349     {
350         const cellShape& cellPoints = cellShapes[wedges[i]];
352         hexLabels[0] = cellPoints[0];
353         hexLabels[1] = cellPoints[1];
354         hexLabels[2] = cellPoints[0];
355         hexLabels[3] = cellPoints[2];
356         hexLabels[4] = cellPoints[3];
357         hexLabels[5] = cellPoints[4];
358         hexLabels[6] = cellPoints[6];
359         hexLabels[7] = cellPoints[5];
361         mcsl[i + offset] = cellShape(hex, hexLabels);
362     }
364     return mcsl;
368 void Foam::ensightMesh::writePrims
370     const cellShapeList& cellShapes,
371     const label pointOffset,
372     OFstream& ensightGeometryFile
373 ) const
375     label po = pointOffset + 1;
377     forAll(cellShapes, i)
378     {
379         const cellShape& cellPoints = cellShapes[i];
381         forAll(cellPoints, pointI)
382         {
383             ensightGeometryFile<< setw(10) << cellPoints[pointI] + po;
384         }
385         ensightGeometryFile << nl;
386     }
390 void Foam::ensightMesh::writePrimsBinary
392     const cellShapeList& cellShapes,
393     const label pointOffset,
394     std::ofstream& ensightGeometryFile
395 ) const
397     label po = pointOffset + 1;
399     if (cellShapes.size())
400     {
401         // All the cellShapes have the same number of elements!
402         int numIntElem = cellShapes.size()*cellShapes[0].size();
403         List<int> temp(numIntElem);
405         int n = 0;
407         forAll(cellShapes, i)
408         {
409             const cellShape& cellPoints = cellShapes[i];
411             forAll(cellPoints, pointI)
412             {
413                 temp[n] = cellPoints[pointI] + po;
414                 n++;
415             }
416         }
418         ensightGeometryFile.write
419         (
420             reinterpret_cast<char*>(temp.begin()),
421             numIntElem*sizeof(int)
422         );
423     }
427 void Foam::ensightMesh::writePolysNFaces
429     const labelList& polys,
430     const cellList& cellFaces,
431     OFstream& ensightGeometryFile
432 ) const
434         forAll(polys, i)
435         {
436             ensightGeometryFile
437                 << setw(10) << cellFaces[polys[i]].size() << nl;
438         }
442 void Foam::ensightMesh::writePolysNPointsPerFace
444     const labelList& polys,
445     const cellList& cellFaces,
446     const faceList& faces,
447     OFstream& ensightGeometryFile
448 ) const
450     forAll(polys, i)
451     {
452         const labelList& cf = cellFaces[polys[i]];
454         forAll(cf, faceI)
455         {
456             ensightGeometryFile
457                 << setw(10) << faces[cf[faceI]].size() << nl;
458         }
459     }
463 void Foam::ensightMesh::writePolysPoints
465     const labelList& polys,
466     const cellList& cellFaces,
467     const faceList& faces,
468     const label pointOffset,
469     OFstream& ensightGeometryFile
470 ) const
472     label po = pointOffset + 1;
474     forAll(polys, i)
475     {
476         const labelList& cf = cellFaces[polys[i]];
478         forAll(cf, faceI)
479         {
480             const face& f = faces[cf[faceI]];
482             forAll(f, pointI)
483             {
484                 ensightGeometryFile << setw(10) << f[pointI] + po;
485             }
486             ensightGeometryFile << nl;
487         }
488     }
492 void Foam::ensightMesh::writeAllPolys
494     const labelList& pointOffsets,
495     OFstream& ensightGeometryFile
496 ) const
498     if (meshCellSets_.nPolys)
499     {
500         const cellList& cellFaces = mesh_.cells();
501         const faceList& faces = mesh_.faces();
503         if (Pstream::master())
504         {
505             ensightGeometryFile
506                 << "nfaced" << nl << setw(10) << meshCellSets_.nPolys << nl;
507         }
509         // Number of faces for each poly cell
510         if (Pstream::master())
511         {
512             // Master
513             writePolysNFaces
514             (
515                 meshCellSets_.polys,
516                 cellFaces,
517                 ensightGeometryFile
518             );
519             // Slaves
520             for (int slave=1; slave<Pstream::nProcs(); slave++)
521             {
522                 IPstream fromSlave(Pstream::scheduled, slave);
523                 labelList polys(fromSlave);
524                 cellList cellFaces(fromSlave);
526                 writePolysNFaces
527                 (
528                     polys,
529                     cellFaces,
530                     ensightGeometryFile
531                 );
532             }
533         }
534         else
535         {
536             OPstream toMaster(Pstream::scheduled, Pstream::masterNo());
537             toMaster<< meshCellSets_.polys << cellFaces;
538         }
540         // Number of points for each face of the above list
541         if (Pstream::master())
542         {
543             // Master
544             writePolysNPointsPerFace
545             (
546                 meshCellSets_.polys,
547                 cellFaces,
548                 faces,
549                 ensightGeometryFile
550             );
551             // Slaves
552             for (int slave=1; slave<Pstream::nProcs(); slave++)
553             {
554                 IPstream fromSlave(Pstream::scheduled, slave);
555                 labelList polys(fromSlave);
556                 cellList cellFaces(fromSlave);
557                 faceList faces(fromSlave);
559                 writePolysNPointsPerFace
560                 (
561                     polys,
562                     cellFaces,
563                     faces,
564                     ensightGeometryFile
565                 );
566             }
567         }
568         else
569         {
570             OPstream toMaster(Pstream::scheduled, Pstream::masterNo());
571             toMaster<< meshCellSets_.polys << cellFaces << faces;
572         }
574         // List of points id for each face of the above list
575         if (Pstream::master())
576         {
577             // Master
578             writePolysPoints
579             (
580                 meshCellSets_.polys,
581                 cellFaces,
582                 faces,
583                 0,
584                 ensightGeometryFile
585             );
586             // Slaves
587             for (int slave=1; slave<Pstream::nProcs(); slave++)
588             {
589                 IPstream fromSlave(Pstream::scheduled, slave);
590                 labelList polys(fromSlave);
591                 cellList cellFaces(fromSlave);
592                 faceList faces(fromSlave);
594                 writePolysPoints
595                 (
596                     polys,
597                     cellFaces,
598                     faces,
599                     pointOffsets[slave-1],
600                     ensightGeometryFile
601                 );
602             }
603         }
604         else
605         {
606             OPstream toMaster(Pstream::scheduled, Pstream::masterNo());
607             toMaster<< meshCellSets_.polys << cellFaces << faces;
608         }
609     }
613 void Foam::ensightMesh::writePolysNFacesBinary
615     const labelList& polys,
616     const cellList& cellFaces,
617     std::ofstream& ensightGeometryFile
618 ) const
620     forAll(polys, i)
621     {
622         writeEnsDataBinary
623         (
624             cellFaces[polys[i]].size(),
625             ensightGeometryFile
626         );
627     }
631 void Foam::ensightMesh::writePolysNPointsPerFaceBinary
633     const labelList& polys,
634     const cellList& cellFaces,
635     const faceList& faces,
636     std::ofstream& ensightGeometryFile
637 ) const
639     forAll(polys, i)
640     {
641         const labelList& cf = cellFaces[polys[i]];
643         forAll(cf, faceI)
644         {
645             writeEnsDataBinary
646             (
647                 faces[cf[faceI]].size(),
648                 ensightGeometryFile
649             );
650         }
651     }
655 void Foam::ensightMesh::writePolysPointsBinary
657     const labelList& polys,
658     const cellList& cellFaces,
659     const faceList& faces,
660     const label pointOffset,
661     std::ofstream& ensightGeometryFile
662 ) const
664     label po = pointOffset + 1;
666     forAll(polys, i)
667     {
668         const labelList& cf = cellFaces[polys[i]];
670         forAll(cf, faceI)
671         {
672             const face& f = faces[cf[faceI]];
674             forAll(f, pointI)
675             {
676                 writeEnsDataBinary(f[pointI] + po,ensightGeometryFile);
677             }
678         }
679     }
683 void Foam::ensightMesh::writeAllPolysBinary
685     const labelList& pointOffsets,
686     std::ofstream& ensightGeometryFile
687 ) const
689     if (meshCellSets_.nPolys)
690     {
691         const cellList& cellFaces = mesh_.cells();
692         const faceList& faces = mesh_.faces();
694         if (Pstream::master())
695         {
696             writeEnsDataBinary("nfaced",ensightGeometryFile);
697             writeEnsDataBinary(meshCellSets_.nPolys,ensightGeometryFile);
698         }
700         // Number of faces for each poly cell
701         if (Pstream::master())
702         {
703             // Master
704             writePolysNFacesBinary
705             (
706                 meshCellSets_.polys,
707                 cellFaces,
708                 ensightGeometryFile
709             );
710             // Slaves
711             for (int slave=1; slave<Pstream::nProcs(); slave++)
712             {
713                 IPstream fromSlave(Pstream::scheduled, slave);
714                 labelList polys(fromSlave);
715                 cellList cellFaces(fromSlave);
717                 writePolysNFacesBinary
718                 (
719                     polys,
720                     cellFaces,
721                     ensightGeometryFile
722                 );
723             }
724         }
725         else
726         {
727             OPstream toMaster(Pstream::scheduled, Pstream::masterNo());
728             toMaster<< meshCellSets_.polys << cellFaces;
729         }
731         // Number of points for each face of the above list
732         if (Pstream::master())
733         {
734             // Master
735             writePolysNPointsPerFaceBinary
736             (
737                 meshCellSets_.polys,
738                 cellFaces,
739                 faces,
740                 ensightGeometryFile
741             );
742             // Slaves
743             for (int slave=1; slave<Pstream::nProcs(); slave++)
744             {
745                 IPstream fromSlave(Pstream::scheduled, slave);
746                 labelList polys(fromSlave);
747                 cellList cellFaces(fromSlave);
748                 faceList faces(fromSlave);
750                 writePolysNPointsPerFaceBinary
751                 (
752                     polys,
753                     cellFaces,
754                     faces,
755                     ensightGeometryFile
756                 );
757             }
758         }
759         else
760         {
761             OPstream toMaster(Pstream::scheduled, Pstream::masterNo());
762             toMaster<< meshCellSets_.polys << cellFaces << faces;
763         }
765         // List of points id for each face of the above list
766         if (Pstream::master())
767         {
768             // Master
769             writePolysPointsBinary
770             (
771                 meshCellSets_.polys,
772                 cellFaces,
773                 faces,
774                 0,
775                 ensightGeometryFile
776             );
777             // Slaves
778             for (int slave=1; slave<Pstream::nProcs(); slave++)
779             {
780                 IPstream fromSlave(Pstream::scheduled, slave);
781                 labelList polys(fromSlave);
782                 cellList cellFaces(fromSlave);
783                 faceList faces(fromSlave);
785                 writePolysPointsBinary
786                 (
787                     polys,
788                     cellFaces,
789                     faces,
790                     pointOffsets[slave-1],
791                     ensightGeometryFile
792                 );
793             }
794         }
795         else
796         {
797             OPstream toMaster(Pstream::scheduled, Pstream::masterNo());
798             toMaster<< meshCellSets_.polys << cellFaces << faces;
799         }
800     }
804 void Foam::ensightMesh::writeAllPrims
806     const char* key,
807     const label nPrims,
808     const cellShapeList& cellShapes,
809     const labelList& pointOffsets,
810     OFstream& ensightGeometryFile
811 ) const
813     if (nPrims)
814     {
815         if (Pstream::master())
816         {
817             ensightGeometryFile << key << nl << setw(10) << nPrims << nl;
819             writePrims(cellShapes, 0, ensightGeometryFile);
821             for (int slave=1; slave<Pstream::nProcs(); slave++)
822             {
823                 IPstream fromSlave(Pstream::scheduled, slave);
824                 cellShapeList cellShapes(fromSlave);
826                 writePrims
827                 (
828                     cellShapes,
829                     pointOffsets[slave-1],
830                     ensightGeometryFile
831                 );
832             }
833         }
834         else
835         {
836             OPstream toMaster(Pstream::scheduled, Pstream::masterNo());
837             toMaster<< cellShapes;
838         }
839     }
843 void Foam::ensightMesh::writeAllPrimsBinary
845     const char* key,
846     const label nPrims,
847     const cellShapeList& cellShapes,
848     const labelList& pointOffsets,
849     std::ofstream& ensightGeometryFile
850 ) const
852     if (nPrims)
853     {
854         if (Pstream::master())
855         {
856             writeEnsDataBinary(key,ensightGeometryFile);
857             writeEnsDataBinary(nPrims,ensightGeometryFile);
859             writePrimsBinary(cellShapes, 0, ensightGeometryFile);
861             for (int slave=1; slave<Pstream::nProcs(); slave++)
862             {
863                 IPstream fromSlave(Pstream::scheduled, slave);
864                 cellShapeList cellShapes(fromSlave);
866                 writePrimsBinary
867                 (
868                     cellShapes,
869                     pointOffsets[slave-1],
870                     ensightGeometryFile
871                 );
872             }
873         }
874         else
875         {
876             OPstream toMaster(Pstream::scheduled, Pstream::masterNo());
877             toMaster<< cellShapes;
878         }
879     }
883 void Foam::ensightMesh::writeFacePrims
885     const faceList& patchFaces,
886     const label pointOffset,
887     OFstream& ensightGeometryFile
888 ) const
890     if (patchFaces.size())
891     {
892         label po = pointOffset + 1;
894         forAll(patchFaces, i)
895         {
896             const face& patchFace = patchFaces[i];
898             forAll(patchFace, pointI)
899             {
900                 ensightGeometryFile << setw(10) << patchFace[pointI] + po;
901             }
902             ensightGeometryFile << nl;
903         }
904     }
908 void Foam::ensightMesh::writeFacePrimsBinary
910     const faceList& patchFaces,
911     const label pointOffset,
912     std::ofstream& ensightGeometryFile
913 ) const
915     if (patchFaces.size())
916     {
917         label po = pointOffset + 1;
919         forAll(patchFaces, i)
920         {
921             const face& patchFace = patchFaces[i];
923             forAll(patchFace, pointI)
924             {
925                 writeEnsDataBinary
926                 (
927                     patchFace[pointI] + po,
928                     ensightGeometryFile
929                 );
930             }
931         }
932     }
936 Foam::faceList Foam::ensightMesh::map
938     const faceList& patchFaces,
939     const labelList& prims
940 ) const
942     faceList ppf(prims.size());
944     forAll (prims, i)
945     {
946         ppf[i] = patchFaces[prims[i]];
947     }
949     return ppf;
953 void Foam::ensightMesh::writeAllFacePrims
955     const char* key,
956     const labelList& prims,
957     const label nPrims,
958     const faceList& patchFaces,
959     const labelList& pointOffsets,
960     const labelList& patchProcessors,
961     OFstream& ensightGeometryFile
962 ) const
964     if (nPrims)
965     {
966         if (Pstream::master())
967         {
968             ensightGeometryFile << key << nl << setw(10) << nPrims << nl;
970             if (&prims != NULL)
971             {
972                 writeFacePrims
973                 (
974                     map(patchFaces, prims),
975                     0,
976                     ensightGeometryFile
977                 );
978             }
980             forAll (patchProcessors, i)
981             {
982                 if (patchProcessors[i] != 0)
983                 {
984                     label slave = patchProcessors[i];
985                     IPstream fromSlave(Pstream::scheduled, slave);
986                     faceList patchFaces(fromSlave);
988                     writeFacePrims
989                     (
990                         patchFaces,
991                         pointOffsets[i],
992                         ensightGeometryFile
993                     );
994                 }
995             }
996         }
997         else if (&prims != NULL)
998         {
999             OPstream toMaster(Pstream::scheduled, Pstream::masterNo());
1000             toMaster<< map(patchFaces, prims);
1001         }
1002     }
1006 void Foam::ensightMesh::writeNSidedNPointsPerFace
1008     const faceList& patchFaces,
1009     OFstream& ensightGeometryFile
1010 ) const
1012     forAll(patchFaces, i)
1013     {
1014         ensightGeometryFile
1015             << setw(10) << patchFaces[i].size() << nl;
1016     }
1020 void Foam::ensightMesh::writeNSidedPoints
1022     const faceList& patchFaces,
1023     const label pointOffset,
1024     OFstream& ensightGeometryFile
1025 ) const
1027     writeFacePrims
1028     (
1029         patchFaces,
1030         pointOffset,
1031         ensightGeometryFile
1032     );
1036 void Foam::ensightMesh::writeAllNSided
1038     const labelList& prims,
1039     const label nPrims,
1040     const faceList& patchFaces,
1041     const labelList& pointOffsets,
1042     const labelList& patchProcessors,
1043     OFstream& ensightGeometryFile
1044 ) const
1046     if (nPrims)
1047     {
1048         if (Pstream::master())
1049         {
1050             ensightGeometryFile
1051                 << "nsided" << nl << setw(10) << nPrims << nl;
1052         }
1054         // Number of points for each face
1055         if (Pstream::master())
1056         {
1057             if (&prims != NULL)
1058             {
1059                 writeNSidedNPointsPerFace
1060                 (
1061                     map(patchFaces, prims),
1062                     ensightGeometryFile
1063                 );
1064             }
1066             forAll (patchProcessors, i)
1067             {
1068                 if (patchProcessors[i] != 0)
1069                 {
1070                     label slave = patchProcessors[i];
1071                     IPstream fromSlave(Pstream::scheduled, slave);
1072                     faceList patchFaces(fromSlave);
1074                     writeNSidedNPointsPerFace
1075                     (
1076                         patchFaces,
1077                         ensightGeometryFile
1078                     );
1079                 }
1080             }
1081         }
1082         else if (&prims != NULL)
1083         {
1084             OPstream toMaster(Pstream::scheduled, Pstream::masterNo());
1085             toMaster<< map(patchFaces, prims);
1086         }
1088         // List of points id for each face
1089         if (Pstream::master())
1090         {
1091             if (&prims != NULL)
1092             {
1093                 writeNSidedPoints
1094                 (
1095                     map(patchFaces, prims),
1096                     0,
1097                     ensightGeometryFile
1098                 );
1099             }
1101             forAll (patchProcessors, i)
1102             {
1103                 if (patchProcessors[i] != 0)
1104                 {
1105                     label slave = patchProcessors[i];
1106                     IPstream fromSlave(Pstream::scheduled, slave);
1107                     faceList patchFaces(fromSlave);
1109                     writeNSidedPoints
1110                     (
1111                         patchFaces,
1112                         pointOffsets[i],
1113                         ensightGeometryFile
1114                     );
1115                 }
1116             }
1117         }
1118         else if (&prims != NULL)
1119         {
1120             OPstream toMaster(Pstream::scheduled, Pstream::masterNo());
1121             toMaster<< map(patchFaces, prims);
1122         }
1123     }
1127 void Foam::ensightMesh::writeNSidedPointsBinary
1129     const faceList& patchFaces,
1130     const label pointOffset,
1131     std::ofstream& ensightGeometryFile
1132 ) const
1134     writeFacePrimsBinary
1135     (
1136         patchFaces,
1137         pointOffset,
1138         ensightGeometryFile
1139     );
1143 void Foam::ensightMesh::writeNSidedNPointsPerFaceBinary
1145     const faceList& patchFaces,
1146     std::ofstream& ensightGeometryFile
1147 ) const
1149     forAll(patchFaces, i)
1150     {
1151         writeEnsDataBinary
1152         (
1153             patchFaces[i].size(),
1154             ensightGeometryFile
1155         );
1156     }
1160 void Foam::ensightMesh::writeAllNSidedBinary
1162     const labelList& prims,
1163     const label nPrims,
1164     const faceList& patchFaces,
1165     const labelList& pointOffsets,
1166     const labelList& patchProcessors,
1167     std::ofstream& ensightGeometryFile
1168 ) const
1170     if (nPrims)
1171     {
1172         if (Pstream::master())
1173         {
1174             writeEnsDataBinary("nsided",ensightGeometryFile);
1175             writeEnsDataBinary(nPrims,ensightGeometryFile);
1176         }
1178         // Number of points for each face
1179         if (Pstream::master())
1180         {
1181             if (&prims != NULL)
1182             {
1183                 writeNSidedNPointsPerFaceBinary
1184                 (
1185                     map(patchFaces, prims),
1186                     ensightGeometryFile
1187                 );
1188             }
1190             forAll (patchProcessors, i)
1191             {
1192                 if (patchProcessors[i] != 0)
1193                 {
1194                     label slave = patchProcessors[i];
1195                     IPstream fromSlave(Pstream::scheduled, slave);
1196                     faceList patchFaces(fromSlave);
1198                     writeNSidedNPointsPerFaceBinary
1199                     (
1200                         patchFaces,
1201                         ensightGeometryFile
1202                     );
1203                 }
1204             }
1205         }
1206         else if (&prims != NULL)
1207         {
1208             OPstream toMaster(Pstream::scheduled, Pstream::masterNo());
1209             toMaster<< map(patchFaces, prims);
1210         }
1212         // List of points id for each face
1213         if (Pstream::master())
1214         {
1215             if (&prims != NULL)
1216             {
1217                 writeNSidedPointsBinary
1218                 (
1219                     map(patchFaces, prims),
1220                     0,
1221                     ensightGeometryFile
1222                 );
1223             }
1225             forAll (patchProcessors, i)
1226             {
1227                 if (patchProcessors[i] != 0)
1228                 {
1229                     label slave = patchProcessors[i];
1230                     IPstream fromSlave(Pstream::scheduled, slave);
1231                     faceList patchFaces(fromSlave);
1233                     writeNSidedPointsBinary
1234                     (
1235                         patchFaces,
1236                         pointOffsets[i],
1237                         ensightGeometryFile
1238                     );
1239                 }
1240             }
1241         }
1242         else if (&prims != NULL)
1243         {
1244             OPstream toMaster(Pstream::scheduled, Pstream::masterNo());
1245             toMaster<< map(patchFaces, prims);
1246         }
1247     }
1251 void Foam::ensightMesh::writeAllFacePrimsBinary
1253     const char* key,
1254     const labelList& prims,
1255     const label nPrims,
1256     const faceList& patchFaces,
1257     const labelList& pointOffsets,
1258     const labelList& patchProcessors,
1259     std::ofstream& ensightGeometryFile
1260 ) const
1262     if (nPrims)
1263     {
1264         if (Pstream::master())
1265         {
1266             writeEnsDataBinary(key,ensightGeometryFile);
1267             writeEnsDataBinary(nPrims,ensightGeometryFile);
1269             if (&prims != NULL)
1270             {
1271                 writeFacePrimsBinary
1272                 (
1273                     map(patchFaces, prims),
1274                     0,
1275                     ensightGeometryFile
1276                 );
1277             }
1279             forAll (patchProcessors, i)
1280             {
1281                 if (patchProcessors[i] != 0)
1282                 {
1283                     label slave = patchProcessors[i];
1284                     IPstream fromSlave(Pstream::scheduled, slave);
1285                     faceList patchFaces(fromSlave);
1287                     writeFacePrimsBinary
1288                     (
1289                         patchFaces,
1290                         pointOffsets[i],
1291                         ensightGeometryFile
1292                     );
1293                 }
1294             }
1295         }
1296         else if (&prims != NULL)
1297         {
1298             OPstream toMaster(Pstream::scheduled, Pstream::masterNo());
1299             toMaster<< map(patchFaces, prims);
1300         }
1301     }
1305 void Foam::ensightMesh::write
1307     const fileName& postProcPath,
1308     const word& prepend,
1309     const label timeIndex,
1310     Ostream& ensightCaseFile
1311 ) const
1313     if (binary_)
1314     {
1315         writeBinary(postProcPath, prepend, timeIndex, ensightCaseFile);
1316     }
1317     else
1318     {
1319         writeAscii(postProcPath, prepend, timeIndex, ensightCaseFile);
1320     }
1324 void Foam::ensightMesh::writeAscii
1326     const fileName& postProcPath,
1327     const word& prepend,
1328     const label timeIndex,
1329     Ostream& ensightCaseFile
1330 ) const
1332     const Time& runTime = mesh_.time();
1333     const pointField& points = mesh_.points();
1334     const cellShapeList& cellShapes = mesh_.cellShapes();
1336     word timeFile = prepend;
1338     if (timeIndex == 0)
1339     {
1340         timeFile += "000.";
1341     }
1342     else if (mesh_.moving())
1343     {
1344         timeFile += itoa(timeIndex) + '.';
1345     }
1347     // set the filename of the ensight file
1348     fileName ensightGeometryFileName = timeFile + "mesh";
1350     OFstream *ensightGeometryFilePtr = NULL;
1351     if (Pstream::master())
1352     {
1353         ensightGeometryFilePtr = new OFstream
1354         (
1355             postProcPath/ensightGeometryFileName,
1356             ios_base::out|ios_base::trunc,
1357             runTime.writeFormat(),
1358             runTime.writeVersion(),
1359             IOstream::UNCOMPRESSED
1360         );
1361     }
1363     OFstream& ensightGeometryFile = *ensightGeometryFilePtr;
1365     if (Pstream::master())
1366     {
1367         // Set Format
1368         ensightGeometryFile.setf
1369         (
1370             ios_base::scientific,
1371             ios_base::floatfield
1372         );
1373         ensightGeometryFile.precision(5);
1375         ensightGeometryFile
1376             << "EnSight Geometry File" << nl
1377             << "written from OpenFOAM-" << Foam::FOAMversion << nl
1378             << "node id assign" << nl
1379             << "element id assign" << nl;
1380     }
1382     labelList pointOffsets(Pstream::nProcs(), 0);
1384     if (patchNames_.empty())
1385     {
1386         label nPoints = points.size();
1387         Pstream::gather(nPoints, sumOp<label>());
1389         if (Pstream::master())
1390         {
1391             ensightGeometryFile
1392                 << "part" << nl
1393                 << setw(10) << 1 << nl
1394                 << "internalMesh" << nl
1395                 << "coordinates" << nl
1396                 << setw(10) << nPoints
1397                 << endl;
1399             for (direction d=0; d<vector::nComponents; d++)
1400             {
1401                 writePoints(points.component(d), ensightGeometryFile);
1402                 pointOffsets[0] = points.size();
1404                 for (int slave=1; slave<Pstream::nProcs(); slave++)
1405                 {
1406                     IPstream fromSlave(Pstream::scheduled, slave);
1407                     scalarField pointsComponent(fromSlave);
1408                     writePoints(pointsComponent, ensightGeometryFile);
1409                     pointOffsets[slave] =
1410                         pointOffsets[slave-1]
1411                       + pointsComponent.size();
1412                 }
1413             }
1414         }
1415         else
1416         {
1417             for (direction d=0; d<vector::nComponents; d++)
1418             {
1419                 OPstream toMaster(Pstream::scheduled, Pstream::masterNo());
1420                 toMaster<< points.component(d);
1421             }
1422         }
1424         writeAllPrims
1425         (
1426             "hexa8",
1427             meshCellSets_.nHexesWedges,
1428             map(cellShapes, meshCellSets_.hexes, meshCellSets_.wedges),
1429             pointOffsets,
1430             ensightGeometryFile
1431         );
1433         writeAllPrims
1434         (
1435             "penta6",
1436             meshCellSets_.nPrisms,
1437             map(cellShapes, meshCellSets_.prisms),
1438             pointOffsets,
1439             ensightGeometryFile
1440         );
1442         writeAllPrims
1443         (
1444             "pyramid5",
1445             meshCellSets_.nPyrs,
1446             map(cellShapes, meshCellSets_.pyrs),
1447             pointOffsets,
1448             ensightGeometryFile
1449         );
1451         writeAllPrims
1452         (
1453             "tetra4",
1454             meshCellSets_.nTets,
1455             map(cellShapes, meshCellSets_.tets),
1456             pointOffsets,
1457             ensightGeometryFile
1458         );
1460         writeAllPolys
1461         (
1462             pointOffsets,
1463             ensightGeometryFile
1464         );
1465     }
1468     label ensightPatchI = patchPartOffset_;
1470     forAll(allPatchNames_, patchi)
1471     {
1472         const word& patchName = allPatchNames_[patchi];
1473         const labelList& patchProcessors = allPatchProcs_[patchi];
1475         if (patchNames_.empty() || patchNames_.found(patchName))
1476         {
1477             const nFacePrimitives& nfp = nPatchPrims_.find(patchName)();
1479             const labelList *trisPtr  = NULL;
1480             const labelList *quadsPtr = NULL;
1481             const labelList *polysPtr = NULL;
1483             const pointField *patchPointsPtr = NULL;
1484             const faceList *patchFacesPtr = NULL;
1486             if (mesh_.boundary()[patchi].size())
1487             {
1488                 const polyPatch& p = mesh_.boundaryMesh()[patchi];
1490                 trisPtr  = &boundaryFaceSets_[patchi].tris;
1491                 quadsPtr = &boundaryFaceSets_[patchi].quads;
1492                 polysPtr = &boundaryFaceSets_[patchi].polys;
1494                 patchPointsPtr = &(p.localPoints());
1495                 patchFacesPtr  = &(p.localFaces());
1496             }
1498             const labelList& tris = *trisPtr;
1499             const labelList& quads = *quadsPtr;
1500             const labelList& polys = *polysPtr;
1501             const pointField& patchPoints = *patchPointsPtr;
1502             const faceList& patchFaces = *patchFacesPtr;
1504             if (nfp.nTris || nfp.nQuads || nfp.nPolys)
1505             {
1506                 labelList patchPointOffsets(Pstream::nProcs(), 0);
1508                 if (Pstream::master())
1509                 {
1510                     ensightGeometryFile
1511                         << "part" << nl
1512                         << setw(10) << ensightPatchI++ << nl
1513                         << patchName << nl
1514                         << "coordinates" << nl
1515                         << setw(10) << nfp.nPoints
1516                         << endl;
1518                     for (direction d=0; d<vector::nComponents; d++)
1519                     {
1520                         if (patchPointsPtr)
1521                         {
1522                             writePoints
1523                             (
1524                                 patchPoints.component(d),
1525                                 ensightGeometryFile
1526                             );
1527                         }
1529                         patchPointOffsets = 0;
1531                         forAll (patchProcessors, i)
1532                         {
1533                             if (patchProcessors[i] != 0)
1534                             {
1535                                 label slave = patchProcessors[i];
1536                                 IPstream fromSlave(Pstream::scheduled, slave);
1537                                 scalarField patchPointsComponent(fromSlave);
1539                                 writePoints
1540                                 (
1541                                     patchPointsComponent,
1542                                     ensightGeometryFile
1543                                 );
1545                                 if (i < Pstream::nProcs()-1)
1546                                 {
1547                                     patchPointOffsets[i+1] =
1548                                         patchPointOffsets[i]
1549                                       + patchPointsComponent.size();
1550                                 }
1551                             }
1552                             else
1553                             {
1554                                 if (i < Pstream::nProcs()-1)
1555                                 {
1556                                     patchPointOffsets[i+1] =
1557                                         patchPointOffsets[i]
1558                                       + patchPoints.size();
1559                                 }
1560                             }
1561                         }
1562                     }
1563                 }
1564                 else if (patchPointsPtr)
1565                 {
1566                     for (direction d=0; d<vector::nComponents; d++)
1567                     {
1568                         OPstream toMaster
1569                         (
1570                             Pstream::scheduled,
1571                             Pstream::masterNo()
1572                         );
1573                         toMaster<< patchPoints.component(d);
1574                     }
1575                 }
1577                 writeAllFacePrims
1578                 (
1579                     "tria3",
1580                     tris,
1581                     nfp.nTris,
1582                     patchFaces,
1583                     patchPointOffsets,
1584                     patchProcessors,
1585                     ensightGeometryFile
1586                 );
1588                 writeAllFacePrims
1589                 (
1590                     "quad4",
1591                     quads,
1592                     nfp.nQuads,
1593                     patchFaces,
1594                     patchPointOffsets,
1595                     patchProcessors,
1596                     ensightGeometryFile
1597                 );
1599                 writeAllNSided
1600                 (
1601                     polys,
1602                     nfp.nPolys,
1603                     patchFaces,
1604                     patchPointOffsets,
1605                     patchProcessors,
1606                     ensightGeometryFile
1607                 );
1608             }
1609         }
1610     }
1612     if (Pstream::master())
1613     {
1614         delete ensightGeometryFilePtr;
1615     }
1619 void Foam::ensightMesh::writeBinary
1621     const fileName& postProcPath,
1622     const word& prepend,
1623     const label timeIndex,
1624     Ostream& ensightCaseFile
1625 ) const
1627     //const Time& runTime = mesh.time();
1628     const pointField& points = mesh_.points();
1629     const cellShapeList& cellShapes = mesh_.cellShapes();
1631     word timeFile = prepend;
1633     if (timeIndex == 0)
1634     {
1635         timeFile += "000.";
1636     }
1637     else if (mesh_.moving())
1638     {
1639         timeFile += itoa(timeIndex) + '.';
1640     }
1642     // set the filename of the ensight file
1643     fileName ensightGeometryFileName = timeFile + "mesh";
1645     std::ofstream *ensightGeometryFilePtr = NULL;
1647     if (Pstream::master())
1648     {
1649         ensightGeometryFilePtr = new std::ofstream
1650         (
1651             (postProcPath/ensightGeometryFileName).c_str(),
1652             ios_base::out | ios_base::binary | ios_base::trunc
1653         );
1654         // Check on file opened?
1655     }
1657     std::ofstream& ensightGeometryFile = *ensightGeometryFilePtr;
1659     if (Pstream::master())
1660     {
1661         writeEnsDataBinary("C binary", ensightGeometryFile);
1662         writeEnsDataBinary("EnSight Geometry File", ensightGeometryFile);
1663         writeEnsDataBinary("written from OpenFOAM", ensightGeometryFile);
1664         writeEnsDataBinary("node id assign", ensightGeometryFile);
1665         writeEnsDataBinary("element id assign", ensightGeometryFile);
1666     }
1668     labelList pointOffsets(Pstream::nProcs(), 0);
1670     if (patchNames_.empty())
1671     {
1672         label nPoints = points.size();
1673         Pstream::gather(nPoints, sumOp<label>());
1675         if (Pstream::master())
1676         {
1677             writeEnsDataBinary("part",ensightGeometryFile);
1678             writeEnsDataBinary(1,ensightGeometryFile);
1679             writeEnsDataBinary("internalMesh",ensightGeometryFile);
1680             writeEnsDataBinary("coordinates",ensightGeometryFile);
1681             writeEnsDataBinary(nPoints,ensightGeometryFile);
1683             for (direction d=0; d<vector::nComponents; d++)
1684             {
1685                 //writePointsBinary(points.component(d), ensightGeometryFile);
1686                 writeEnsDataBinary(points.component(d), ensightGeometryFile);
1687                 pointOffsets[0] = points.size();
1689                 for (int slave=1; slave<Pstream::nProcs(); slave++)
1690                 {
1691                     IPstream fromSlave(Pstream::scheduled, slave);
1692                     scalarField pointsComponent(fromSlave);
1693                     //writePointsBinary(pointsComponent, ensightGeometryFile);
1694                     writeEnsDataBinary(pointsComponent, ensightGeometryFile);
1695                     pointOffsets[slave] =
1696                         pointOffsets[slave-1]
1697                       + pointsComponent.size();
1698                 }
1699             }
1700         }
1701         else
1702         {
1703             for (direction d=0; d<vector::nComponents; d++)
1704             {
1705                 OPstream toMaster(Pstream::scheduled, Pstream::masterNo());
1706                 toMaster<< points.component(d);
1707             }
1708         }
1710         writeAllPrimsBinary
1711         (
1712             "hexa8",
1713             meshCellSets_.nHexesWedges,
1714             map(cellShapes, meshCellSets_.hexes, meshCellSets_.wedges),
1715             pointOffsets,
1716             ensightGeometryFile
1717         );
1719         writeAllPrimsBinary
1720         (
1721             "penta6",
1722             meshCellSets_.nPrisms,
1723             map(cellShapes, meshCellSets_.prisms),
1724             pointOffsets,
1725             ensightGeometryFile
1726         );
1728         writeAllPrimsBinary
1729         (
1730             "pyramid5",
1731             meshCellSets_.nPyrs,
1732             map(cellShapes, meshCellSets_.pyrs),
1733             pointOffsets,
1734             ensightGeometryFile
1735         );
1737         writeAllPrimsBinary
1738         (
1739             "tetra4",
1740             meshCellSets_.nTets,
1741             map(cellShapes, meshCellSets_.tets),
1742             pointOffsets,
1743             ensightGeometryFile
1744         );
1746         writeAllPolysBinary
1747         (
1748             pointOffsets,
1749             ensightGeometryFile
1750         );
1752     }
1754     label ensightPatchI = patchPartOffset_;
1755     label iCount = 0;
1757     forAll(allPatchNames_, patchi)
1758     {
1759         iCount ++;
1760         const word& patchName = allPatchNames_[patchi];
1761         const labelList& patchProcessors = allPatchProcs_[patchi];
1763         if (patchNames_.empty() || patchNames_.found(patchName))
1764         {
1765             const nFacePrimitives& nfp = nPatchPrims_.find(patchName)();
1767             const labelList *trisPtr = NULL;
1768             const labelList *quadsPtr = NULL;
1769             const labelList *polysPtr = NULL;
1771             const pointField *patchPointsPtr = NULL;
1772             const faceList *patchFacesPtr = NULL;
1774             if (mesh_.boundary()[patchi].size())
1775             {
1776                 const polyPatch& p = mesh_.boundaryMesh()[patchi];
1778                 trisPtr = &boundaryFaceSets_[patchi].tris;
1779                 quadsPtr = &boundaryFaceSets_[patchi].quads;
1780                 polysPtr = &boundaryFaceSets_[patchi].polys;
1782                 patchPointsPtr = &(p.localPoints());
1783                 patchFacesPtr = &(p.localFaces());
1784             }
1786             const labelList& tris = *trisPtr;
1787             const labelList& quads = *quadsPtr;
1788             const labelList& polys = *polysPtr;
1789             const pointField& patchPoints = *patchPointsPtr;
1790             const faceList& patchFaces = *patchFacesPtr;
1792             if (nfp.nTris || nfp.nQuads || nfp.nPolys)
1793             {
1794                 labelList patchPointOffsets(Pstream::nProcs(), 0);
1796                 if (Pstream::master())
1797                 {
1798                     writeEnsDataBinary("part",ensightGeometryFile);
1799                     writeEnsDataBinary(ensightPatchI++,ensightGeometryFile);
1800                     //writeEnsDataBinary(patchName.c_str(),ensightGeometryFile);
1801                     writeEnsDataBinary(patchName.c_str(),ensightGeometryFile);
1802                     writeEnsDataBinary("coordinates",ensightGeometryFile);
1803                     writeEnsDataBinary(nfp.nPoints,ensightGeometryFile);
1805                     for (direction d=0; d<vector::nComponents; d++)
1806                     {
1807                         if (patchPointsPtr)
1808                         {
1809                             //writePointsBinary
1810                             writeEnsDataBinary
1811                             (
1812                                 patchPoints.component(d),
1813                                 ensightGeometryFile
1814                             );
1815                         }
1817                         patchPointOffsets = 0;
1820                         forAll (patchProcessors, i)
1821                         {
1822                             if (patchProcessors[i] != 0)
1823                             {
1824                                 label slave = patchProcessors[i];
1825                                 IPstream fromSlave(Pstream::scheduled, slave);
1826                                 scalarField patchPointsComponent(fromSlave);
1828                                 //writePointsBinary
1829                                 writeEnsDataBinary
1830                                 (
1831                                     patchPointsComponent,
1832                                     ensightGeometryFile
1833                                 );
1835                                 if (i < Pstream::nProcs()-1)
1836                                 {
1837                                     patchPointOffsets[i+1] =
1838                                         patchPointOffsets[i]
1839                                       + patchPointsComponent.size();
1840                                 }
1841                             }
1842                             else
1843                             {
1844                                 if (i < Pstream::nProcs()-1)
1845                                 {
1846                                     patchPointOffsets[i+1] =
1847                                         patchPointOffsets[i]
1848                                       + patchPoints.size();
1849                                 }
1850                             }
1851                         }
1852                     }
1853                 }
1854                 else if (patchPointsPtr)
1855                 {
1856                     for (direction d=0; d<vector::nComponents; d++)
1857                     {
1858                         OPstream toMaster
1859                         (
1860                             Pstream::scheduled,
1861                             Pstream::masterNo()
1862                         );
1863                         toMaster<< patchPoints.component(d);
1864                     }
1865                 }
1867                 writeAllFacePrimsBinary
1868                 (
1869                     "tria3",
1870                     tris,
1871                     nfp.nTris,
1872                     patchFaces,
1873                     patchPointOffsets,
1874                     patchProcessors,
1875                     ensightGeometryFile
1876                 );
1878                 writeAllFacePrimsBinary
1879                 (
1880                     "quad4",
1881                     quads,
1882                     nfp.nQuads,
1883                     patchFaces,
1884                     patchPointOffsets,
1885                     patchProcessors,
1886                     ensightGeometryFile
1887                 );
1889                 writeAllNSidedBinary
1890                 (
1891                     polys,
1892                     nfp.nPolys,
1893                     patchFaces,
1894                     patchPointOffsets,
1895                     patchProcessors,
1896                     ensightGeometryFile
1897                 );
1898             }
1899         }
1900     }
1903     if (Pstream::master())
1904     {
1905         delete ensightGeometryFilePtr;
1906     }
1910 // ************************************************************************* //