initial commit for version 1.5.x patch release
[OpenFOAM-1.5.x.git] / applications / utilities / postProcessing / dataConversion / foamToEnsight / ensightMesh.C
blob14cd4c53ef145d4c03aada917c0e77b7d95b6b1a
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2008 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
11     OpenFOAM is free software; you can redistribute it and/or modify it
12     under the terms of the GNU General Public License as published by the
13     Free Software Foundation; either version 2 of the License, or (at your
14     option) any later version.
16     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
19     for more details.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 #include "argList.H"
28 #include "Time.H"
29 #include "ensightMesh.H"
30 #include "fvMesh.H"
31 #include "PstreamCombineReduceOps.H"
32 #include "processorPolyPatch.H"
33 #include "cellModeller.H"
34 #include "IOmanip.H"
35 #include "itoa.H"
36 #include "ensightWriteBinary.H"
37 #include <fstream>
39 // * * * * * * * * * * * * * Private Functions * * * * * * * * * * * * * * //
41 namespace Foam
44 class concatPatchNames
47 public:
49     void operator()
50     (
51         HashTable<labelList>& x,
52         const HashTable<labelList>& y
53     ) const
54     {
55         forAllConstIter(HashTable<labelList>, y, iter)
56         {
57             HashTable<labelList>::iterator xiter = x.find(iter.key());
59             if (xiter == x.end())
60             {
61                 x.insert(iter.key(), iter());
62             }
63             else
64             {
65                 labelList& xPatches = xiter();
66                 const labelList& yPatches = iter();
68                 label offset = xPatches.size();
69                 xPatches.setSize(offset + yPatches.size());
71                 forAll(yPatches, i)
72                 {
73                     xPatches[i + offset] = yPatches[i];
74                 }
75             }
76         }
77     }
80 } // End namespace Foam
83 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
85 Foam::ensightMesh::ensightMesh
87     const fvMesh& mesh,
88     const argList& args,
89     const bool binary
92     binary_(binary),
93     mesh_(mesh),
94     meshCellSets_(mesh_.nCells()),
95     boundaryFaceSets_(mesh_.boundary().size()),
96     allPatchNames_(0),
97     patchIndices_(0),
98     patchNames_(0),
99     nPatchPrims_(0)
101     forAll (mesh_.boundaryMesh(), patchi)
102     {
103         if
104         (
105             typeid(mesh_.boundaryMesh()[patchi])
106          != typeid(processorPolyPatch)
107         )
108         {
109             if (!allPatchNames_.found(mesh_.boundaryMesh()[patchi].name()))
110             {
111                 allPatchNames_.insert
112                 (
113                     mesh_.boundaryMesh()[patchi].name(),
114                     labelList(1, Pstream::myProcNo())
115                 );
117                 patchIndices_.insert
118                 (
119                     mesh_.boundaryMesh()[patchi].name(),
120                     patchi
121                 );
122             }
123         }
124     }
126     combineReduce(allPatchNames_, concatPatchNames());
128     if (args.options().found("patches"))
129     {
130         wordList patchNameList(IStringStream(args.options()["patches"])());
132         if (!patchNameList.size())
133         {
134             patchNameList = allPatchNames_.toc();
135         }
137         forAll (patchNameList, i)
138         {
139             patchNames_.insert(patchNameList[i]);
140         }
141     }
143     const cellShapeList& cellShapes = mesh.cellShapes();
145     const cellModel& tet = *(cellModeller::lookup("tet"));
146     const cellModel& pyr = *(cellModeller::lookup("pyr"));
147     const cellModel& prism = *(cellModeller::lookup("prism"));
148     const cellModel& wedge = *(cellModeller::lookup("wedge"));
149     const cellModel& hex = *(cellModeller::lookup("hex"));
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     // Count the shapes
159     label nTets = 0;
160     label nPyrs = 0;
161     label nPrisms = 0;
162     label nWedges = 0;
163     label nHexes = 0;
164     label nPolys = 0;
166     if (!patchNames_.size())
167     {
168         forAll(cellShapes, celli)
169         {
170             const cellShape& cellShape = cellShapes[celli];
171             const cellModel& cellModel = cellShape.model();
173             if (cellModel == tet)
174             {
175                 tets[nTets++] = celli;
176             }
177             else if (cellModel == pyr)
178             {
179                 pyrs[nPyrs++] = celli;
180             }
181             else if (cellModel == prism)
182             {
183                 prisms[nPrisms++] = celli;
184             }
185             else if (cellModel == wedge)
186             {
187                 wedges[nWedges++] = celli;
188             }
189             else if (cellModel == hex)
190             {
191                 hexes[nHexes++] = celli;
192             }
193             else
194             {
195                 polys[nPolys++] = celli;
196             }
197         }
199         tets.setSize(nTets);
200         pyrs.setSize(nPyrs);
201         prisms.setSize(nPrisms);
202         wedges.setSize(nWedges);
203         hexes.setSize(nHexes);
204         polys.setSize(nPolys);
206         meshCellSets_.nTets = nTets;
207         reduce(meshCellSets_.nTets, sumOp<label>());
209         meshCellSets_.nPyrs = nPyrs;
210         reduce(meshCellSets_.nPyrs, sumOp<label>());
212         meshCellSets_.nPrisms = nPrisms;
213         reduce(meshCellSets_.nPrisms, sumOp<label>());
215         meshCellSets_.nHexesWedges = nHexes + nWedges;
216         reduce(meshCellSets_.nHexesWedges, sumOp<label>());
218         meshCellSets_.nPolys = nPolys;
219         reduce(meshCellSets_.nPolys, sumOp<label>());
220     }
223     forAll (mesh.boundary(), patchi)
224     {
225         if (mesh.boundary()[patchi].size())
226         {
227             const polyPatch& p = mesh.boundaryMesh()[patchi];
229             labelList& tris = boundaryFaceSets_[patchi].tris;
230             labelList& quads = boundaryFaceSets_[patchi].quads;
231             labelList& polys = boundaryFaceSets_[patchi].polys;
233             tris.setSize(p.size());
234             quads.setSize(p.size());
235             polys.setSize(p.size());
237             label nTris = 0;
238             label nQuads = 0;
239             label nPolys = 0;
241             forAll(p, facei)
242             {
243                 const face& f = p[facei];
245                 if (f.size() == 3)
246                 {
247                     tris[nTris++] = facei;
248                 }
249                 else if (f.size() == 4)
250                 {
251                     quads[nQuads++] = facei;
252                 }
253                 else
254                 {
255                     polys[nPolys++] = facei;
256                 }
257             }
259             tris.setSize(nTris);
260             quads.setSize(nQuads);
261             polys.setSize(nPolys);
262         }
263     }
265     forAllConstIter(HashTable<labelList>, allPatchNames_, iter)
266     {
267         const word& patchName = iter.key();
268         nFacePrimitives nfp;
270         if (!patchNames_.size() || patchNames_.found(patchName))
271         {
272             if (patchIndices_.found(patchName))
273             {
274                 label patchi = patchIndices_.find(patchName)();
276                 nfp.nPoints = mesh.boundaryMesh()[patchi].localPoints().size();
277                 nfp.nTris = boundaryFaceSets_[patchi].tris.size();
278                 nfp.nQuads = boundaryFaceSets_[patchi].quads.size();
279                 nfp.nPolys = boundaryFaceSets_[patchi].polys.size();
280             }
281         }
283         reduce(nfp.nPoints, sumOp<label>());
284         reduce(nfp.nTris, sumOp<label>());
285         reduce(nfp.nQuads, sumOp<label>());
286         reduce(nfp.nPolys, sumOp<label>());
288         nPatchPrims_.insert(patchName, nfp);
289     }
293 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
295 Foam::ensightMesh::~ensightMesh()
299 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
301 void Foam::ensightMesh::writePoints
303     const scalarField& pointsComponent,
304     OFstream& ensightGeometryFile
305 ) const
307     forAll(pointsComponent, pointi)
308     {
309         ensightGeometryFile<< setw(12) << float(pointsComponent[pointi]) << nl;
310     }
314 Foam::cellShapeList Foam::ensightMesh::map
316     const cellShapeList& cellShapes,
317     const labelList& prims
318 ) const
320     cellShapeList mcsl(prims.size());
322     forAll(prims, i)
323     {
324         mcsl[i] = cellShapes[prims[i]];
325     }
327     return mcsl;
331 Foam::cellShapeList Foam::ensightMesh::map
333     const cellShapeList& cellShapes,
334     const labelList& hexes,
335     const labelList& wedges
336 ) const
338     cellShapeList mcsl(hexes.size() + wedges.size());
340     forAll(hexes, i)
341     {
342         mcsl[i] = cellShapes[hexes[i]];
343     }
345     label offset = hexes.size();
347     const cellModel& hex = *(cellModeller::lookup("hex"));
348     labelList hexLabels(8);
350     forAll(wedges, i)
351     {
352         const cellShape& cellPoints = cellShapes[wedges[i]];
354         hexLabels[0] = cellPoints[0];
355         hexLabels[1] = cellPoints[1];
356         hexLabels[2] = cellPoints[0];
357         hexLabels[3] = cellPoints[2];
358         hexLabels[4] = cellPoints[3];
359         hexLabels[5] = cellPoints[4];
360         hexLabels[6] = cellPoints[6];
361         hexLabels[7] = cellPoints[5];
363         mcsl[i + offset] = cellShape(hex, hexLabels);
364     }
366     return mcsl;
370 void Foam::ensightMesh::writePrims
372     const cellShapeList& cellShapes,
373     const label pointOffset,
374     OFstream& ensightGeometryFile
375 ) const
377     label po = pointOffset + 1;
379     forAll(cellShapes, i)
380     {
381         const cellShape& cellPoints = cellShapes[i];
383         forAll(cellPoints, pointi)
384         {
385             ensightGeometryFile<< setw(10) << cellPoints[pointi] + po;
386         }
387         ensightGeometryFile << nl;
388     }
392 void Foam::ensightMesh::writePrimsBinary
394     const cellShapeList& cellShapes,
395     const label pointOffset,
396     std::ofstream& ensightGeometryFile
397 ) const
399     label po = pointOffset + 1;
401     // Create a temp int array
402     int numElem;
404     numElem = cellShapes.size();
406     if (cellShapes.size() > 0)
407     {
408         // All the cellShapes have the same number of elements!
409         int numIntElem = cellShapes.size()*cellShapes[0].size();
410         List<int> temp(numIntElem);
412         int n = 0;
414         forAll(cellShapes, i)
415         {
416             const cellShape& cellPoints = cellShapes[i];
418             forAll(cellPoints, pointi)
419             {
420                 temp[n] = cellPoints[pointi] + po;
421                 n++;
422             }
423         }
425         ensightGeometryFile.write
426         (
427             reinterpret_cast<char*>(temp.begin()),
428             numIntElem*sizeof(int)
429         );
430     }
434 void Foam::ensightMesh::writePolys
436     const labelList& polys,
437     const cellList& cellFaces,
438     const faceList& faces,
439     const label pointOffset,
440     OFstream& ensightGeometryFile
441 ) const
443     if (polys.size())
444     {
445         ensightGeometryFile
446             << "nfaced" << nl << setw(10) << polys.size() << nl;
448         label po = pointOffset + 1;
450         forAll(polys, i)
451         {
452             ensightGeometryFile
453                 << setw(10) << cellFaces[polys[i]].size() << nl;
454         }
456         forAll(polys, i)
457         {
458             const labelList& cf = cellFaces[polys[i]];
460             forAll(cf, facei)
461             {
462                 ensightGeometryFile
463                     << setw(10) << faces[cf[facei]].size() << nl;
464             }
465         }
467         forAll(polys, i)
468         {
469             const labelList& cf = cellFaces[polys[i]];
471             forAll(cf, facei)
472             {
473                 const face& f = faces[cf[facei]];
475                 forAll(f, pointi)
476                 {
477                     ensightGeometryFile << setw(10) << f[pointi] + po;
478                 }
479                 ensightGeometryFile << nl;
480             }
481         }
482     }
486 void Foam::ensightMesh::writePolysBinary
488     const labelList& polys,
489     const cellList& cellFaces,
490     const faceList& faces,
491     const label pointOffset,
492     std::ofstream& ensightGeometryFile
493 ) const
495     if (polys.size())
496     {
497         writeEnsDataBinary("nfaced",ensightGeometryFile);
498         writeEnsDataBinary(polys.size(),ensightGeometryFile);
500         label po = pointOffset + 1;
502         //TODO No buffer at the moment. To be done for speed purposes!
503         forAll(polys, i)
504         {
505             writeEnsDataBinary
506             (
507                 cellFaces[polys[i]].size(),
508                 ensightGeometryFile
509             );
510         }
512         forAll(polys, i)
513         {
514             const labelList& cf = cellFaces[polys[i]];
516             forAll(cf, facei)
517             {
518                 writeEnsDataBinary
519                 (
520                     faces[cf[facei]].size(),
521                     ensightGeometryFile
522                 );
523             }
524         }
526         forAll(polys, i)
527         {
528             const labelList& cf = cellFaces[polys[i]];
530             forAll(cf, facei)
531             {
532                 const face& f = faces[cf[facei]];
534                 forAll(f, pointi)
535                 {
536                     writeEnsDataBinary(f[pointi] + po,ensightGeometryFile);
537                 }
538             }
539         }
540     }
544 void Foam::ensightMesh::writeAllPrims
546     const char* key,
547     const label nPrims,
548     const cellShapeList& cellShapes,
549     const labelList& pointOffsets,
550     OFstream& ensightGeometryFile
551 ) const
553     if (nPrims)
554     {
555         if (Pstream::master())
556         {
557             ensightGeometryFile << key << nl << setw(10) << nPrims << nl;
559             writePrims(cellShapes, 0, ensightGeometryFile);
561             for (int slave=1; slave<Pstream::nProcs(); slave++)
562             {
563                 IPstream fromSlave(Pstream::scheduled, slave);
564                 cellShapeList cellShapes(fromSlave);
566                 writePrims
567                 (
568                     cellShapes,
569                     pointOffsets[slave-1],
570                     ensightGeometryFile
571                 );
572             }
573         }
574         else
575         {
576             OPstream toMaster(Pstream::scheduled, Pstream::masterNo());
577             toMaster<< cellShapes;
578         }
579     }
583 void Foam::ensightMesh::writeAllPrimsBinary
585     const char* key,
586     const label nPrims,
587     const cellShapeList& cellShapes,
588     const labelList& pointOffsets,
589     std::ofstream& ensightGeometryFile
590 ) const
592     if (nPrims)
593     {
594         if (Pstream::master())
595         {
596             writeEnsDataBinary(key,ensightGeometryFile);
597             writeEnsDataBinary(nPrims,ensightGeometryFile);
599             writePrimsBinary(cellShapes, 0, ensightGeometryFile);
601             for (int slave=1; slave<Pstream::nProcs(); slave++)
602             {
603                 IPstream fromSlave(Pstream::scheduled, slave);
604                 cellShapeList cellShapes(fromSlave);
606                 writePrimsBinary
607                 (
608                     cellShapes,
609                     pointOffsets[slave-1],
610                     ensightGeometryFile
611                 );
612             }
613         }
614         else
615         {
616             OPstream toMaster(Pstream::scheduled, Pstream::masterNo());
617             toMaster<< cellShapes;
618         }
619     }
623 void Foam::ensightMesh::writeFacePrims
625     const char* key,
626     const faceList& patchFaces,
627     const label pointOffset,
628     OFstream& ensightGeometryFile
629 ) const
631     if (patchFaces.size())
632     {
633         if (word(key) == "nsided")
634         {
635             ensightGeometryFile
636                 << key << nl << setw(10) << patchFaces.size() << nl;
638             forAll(patchFaces, i)
639             {
640                 ensightGeometryFile
641                     << setw(10) << patchFaces[i].size() << nl;
642             }
643         }
645         label po = pointOffset + 1;
647         forAll(patchFaces, i)
648         {
649             const face& patchFace = patchFaces[i];
651             forAll(patchFace, pointi)
652             {
653                 ensightGeometryFile << setw(10) << patchFace[pointi] + po;
654             }
655             ensightGeometryFile << nl;
656         }
657     }
661 void Foam::ensightMesh::writeFacePrimsBinary
663     const char* key,
664     const faceList& patchFaces,
665     const label pointOffset,
666     std::ofstream& ensightGeometryFile
667 ) const
669     if (patchFaces.size())
670     {
671         //TODO No buffer at the moment. To be done for speed purposes!
672         if (word(key) == "nsided")
673         {
674             writeEnsDataBinary(key,ensightGeometryFile);
675             writeEnsDataBinary(patchFaces.size(),ensightGeometryFile);
677             forAll(patchFaces, i)
678             {
679                 writeEnsDataBinary
680                 (
681                     patchFaces[i].size(),
682                     ensightGeometryFile
683                 );
684             }
685         }
687         label po = pointOffset + 1;
689         forAll(patchFaces, i)
690         {
691             const face& patchFace = patchFaces[i];
693             forAll(patchFace, pointi)
694             {
695                 writeEnsDataBinary
696                 (
697                     patchFace[pointi] + po,
698                     ensightGeometryFile
699                 );
700             }
701         }
702     }
706 Foam::faceList Foam::ensightMesh::map
708     const faceList& patchFaces,
709     const labelList& prims
710 ) const
712     faceList ppf(prims.size());
714     forAll (prims, i)
715     {
716         ppf[i] = patchFaces[prims[i]];
717     }
719     return ppf;
723 void Foam::ensightMesh::writeAllFacePrims
725     const char* key,
726     const labelList& prims,
727     const label nPrims,
728     const faceList& patchFaces,
729     const labelList& pointOffsets,
730     const labelList& patchProcessors,
731     OFstream& ensightGeometryFile
732 ) const
734     if (nPrims)
735     {
736         if (Pstream::master())
737         {
738             if (word(key) != "nsided")
739             {
740                 ensightGeometryFile << key << nl << setw(10) << nPrims << nl;
741             }
743             if (&prims != NULL)
744             {
745                 writeFacePrims
746                 (
747                     key,
748                     map(patchFaces, prims),
749                     0,
750                     ensightGeometryFile
751                 );
752             }
754             forAll (patchProcessors, i)
755             {
756                 if (patchProcessors[i] != 0)
757                 {
758                     label slave = patchProcessors[i];
759                     IPstream fromSlave(Pstream::scheduled, slave);
760                     faceList patchFaces(fromSlave);
762                     writeFacePrims
763                     (
764                         key,
765                         patchFaces,
766                         pointOffsets[i],
767                         ensightGeometryFile
768                     );
769                 }
770             }
771         }
772         else if (&prims != NULL)
773         {
774             OPstream toMaster(Pstream::scheduled, Pstream::masterNo());
775             toMaster<< map(patchFaces, prims);
776         }
777     }
781 void Foam::ensightMesh::writeAllFacePrimsBinary
783     const char* key,
784     const labelList& prims,
785     const label nPrims,
786     const faceList& patchFaces,
787     const labelList& pointOffsets,
788     const labelList& patchProcessors,
789     std::ofstream& ensightGeometryFile
790 ) const
792     if (nPrims)
793     {
794         if (Pstream::master())
795         {
796             if (word(key) != "nsided")
797             {
798                 writeEnsDataBinary(key,ensightGeometryFile);
799                 writeEnsDataBinary(nPrims,ensightGeometryFile);
800             }
802             if (&prims != NULL)
803             {
804                 writeFacePrimsBinary
805                 (
806                     key,
807                     map(patchFaces, prims),
808                     0,
809                     ensightGeometryFile
810                 );
811             }
813             forAll (patchProcessors, i)
814             {
815                 if (patchProcessors[i] != 0)
816                 {
817                     label slave = patchProcessors[i];
818                     IPstream fromSlave(Pstream::scheduled, slave);
819                     faceList patchFaces(fromSlave);
821                     writeFacePrimsBinary
822                     (
823                         key,
824                         patchFaces,
825                         pointOffsets[i],
826                         ensightGeometryFile
827                     );
828                 }
829             }
830         }
831         else if (&prims != NULL)
832         {
833             OPstream toMaster(Pstream::scheduled, Pstream::masterNo());
834             toMaster<< map(patchFaces, prims);
835         }
836     }
840 void Foam::ensightMesh::write
842     const fileName& postProcPath,
843     const word& prepend,
844     const label timeIndex,
845     Ostream& ensightCaseFile
846 ) const
848     if (binary_)
849     {
850         writeBinary(postProcPath, prepend, timeIndex, ensightCaseFile);
851     }
852     else
853     {
854         writeAscii(postProcPath, prepend, timeIndex, ensightCaseFile);
855     }
859 void Foam::ensightMesh::writeAscii
861     const fileName& postProcPath,
862     const word& prepend,
863     const label timeIndex,
864     Ostream& ensightCaseFile
865 ) const
867     const Time& runTime = mesh_.time();
868     const pointField& points = mesh_.points();
869     const cellList& cellFaces = mesh_.cells();
870     const faceList& faces = mesh_.faces();
871     const cellShapeList& cellShapes = mesh_.cellShapes();
873     word timeFile = prepend;
875     if (timeIndex == 0)
876     {
877         timeFile += "000.";
878     }
879     else if (mesh_.moving())
880     {
881         timeFile += itoa(timeIndex) + '.';
882     }
884     // set the filename of the ensight file
885     fileName ensightGeometryFileName = timeFile + "mesh";
887     OFstream *ensightGeometryFilePtr = NULL;
888     if (Pstream::master())
889     {
890         ensightGeometryFilePtr = new OFstream
891         (
892             postProcPath/ensightGeometryFileName,
893             runTime.writeFormat(),
894             runTime.writeVersion(),
895             runTime.writeCompression()
896         );
897     }
899     OFstream& ensightGeometryFile = *ensightGeometryFilePtr;
901     if (Pstream::master())
902     {
903         // Set Format
904         ensightGeometryFile.setf
905         (
906             ios_base::scientific,
907             ios_base::floatfield
908         );
909         ensightGeometryFile.precision(5);
911         ensightGeometryFile
912             << "OpenFOAM Geometry File " << nl
913             << "OpenFOAM version " << Foam::FOAMversion << nl
914             << "node id assign" << nl
915             << "element id assign" << nl;
916     }
918     labelList pointOffsets(Pstream::nProcs(), 0);
920     if (!patchNames_.size())
921     {
922         label nPoints = points.size();
923         Pstream::gather(nPoints, sumOp<label>());
925         if (Pstream::master())
926         {
927             ensightGeometryFile
928                 << "part" << nl
929                 << setw(10) << 1 << nl
930                 << "FOAM cells" << nl
931                 << "coordinates" << nl
932                 << setw(10) << nPoints
933                 << endl;
935             for (direction d=0; d<vector::nComponents; d++)
936             {
937                 writePoints(points.component(d), ensightGeometryFile);
938                 pointOffsets[0] = points.size();
940                 for (int slave=1; slave<Pstream::nProcs(); slave++)
941                 {
942                     IPstream fromSlave(Pstream::scheduled, slave);
943                     scalarField pointsComponent(fromSlave);
944                     writePoints(pointsComponent, ensightGeometryFile);
945                     pointOffsets[slave] =
946                         pointOffsets[slave-1]
947                       + pointsComponent.size();
948                 }
949             }
950         }
951         else
952         {
953             for (direction d=0; d<vector::nComponents; d++)
954             {
955                 OPstream toMaster(Pstream::scheduled, Pstream::masterNo());
956                 toMaster<< points.component(d);
957             }
958         }
960         writeAllPrims
961         (
962             "hexa8",
963             meshCellSets_.nHexesWedges,
964             map(cellShapes, meshCellSets_.hexes, meshCellSets_.wedges),
965             pointOffsets,
966             ensightGeometryFile
967         );
969         writeAllPrims
970         (
971             "penta6",
972             meshCellSets_.nPrisms,
973             map(cellShapes, meshCellSets_.prisms),
974             pointOffsets,
975             ensightGeometryFile
976         );
978         writeAllPrims
979         (
980             "pyramid5",
981             meshCellSets_.nPyrs,
982             map(cellShapes, meshCellSets_.pyrs),
983             pointOffsets,
984             ensightGeometryFile
985         );
987         writeAllPrims
988         (
989             "tetra4",
990             meshCellSets_.nTets,
991             map(cellShapes, meshCellSets_.tets),
992             pointOffsets,
993             ensightGeometryFile
994         );
997         if (meshCellSets_.nPolys)
998         {
999             if (Pstream::master())
1000             {
1001                 /*
1002                 ensightGeometryFile
1003                     << "nfaced" << nl
1004                     << setw(10) << meshCellSets_.nPolys << nl;
1005                 */
1006                 writePolys
1007                 (
1008                     meshCellSets_.polys,
1009                     cellFaces,
1010                     faces,
1011                     0,
1012                     ensightGeometryFile
1013                 );
1015                 for (int slave=1; slave<Pstream::nProcs(); slave++)
1016                 {
1017                     IPstream fromSlave(Pstream::scheduled, slave);
1018                     labelList polys(fromSlave);
1019                     cellList cellFaces(fromSlave);
1020                     faceList faces(fromSlave);
1022                     writePolys
1023                     (
1024                         polys,
1025                         cellFaces,
1026                         faces,
1027                         pointOffsets[slave-1],
1028                         ensightGeometryFile
1029                     );
1030                 }
1031             }
1032             else
1033             {
1034                 OPstream toMaster(Pstream::scheduled, Pstream::masterNo());
1035                 toMaster<< meshCellSets_.polys << cellFaces << faces;
1036             }
1037         }
1038     }
1041     label ensightPatchi = 2;
1043     forAllConstIter(HashTable<labelList>, allPatchNames_, iter)
1044     {
1045         const labelList& patchProcessors = iter();
1047         if (!patchNames_.size() || patchNames_.found(iter.key()))
1048         {
1049             const word& patchName = iter.key();
1050             const nFacePrimitives& nfp = nPatchPrims_.find(patchName)();
1052             const labelList *trisPtr = NULL;
1053             const labelList *quadsPtr = NULL;
1054             const labelList *polysPtr = NULL;
1056             const pointField *patchPointsPtr = NULL;
1057             const faceList *patchFacesPtr = NULL;
1059             if (patchIndices_.found(iter.key()))
1060             {
1061                 label patchi = patchIndices_.find(iter.key())();
1062                 const polyPatch& p = mesh_.boundaryMesh()[patchi];
1064                 trisPtr = &boundaryFaceSets_[patchi].tris;
1065                 quadsPtr = &boundaryFaceSets_[patchi].quads;
1066                 polysPtr = &boundaryFaceSets_[patchi].polys;
1068                 patchPointsPtr = &(p.localPoints());
1069                 patchFacesPtr = &(p.localFaces());
1070             }
1072             const labelList& tris = *trisPtr;
1073             const labelList& quads = *quadsPtr;
1074             const labelList& polys = *polysPtr;
1075             const pointField& patchPoints = *patchPointsPtr;
1076             const faceList& patchFaces = *patchFacesPtr;
1078             if (nfp.nTris || nfp.nQuads || nfp.nPolys)
1079             {
1080                 labelList patchPointOffsets(Pstream::nProcs(), 0);
1082                 if (Pstream::master())
1083                 {
1084                     ensightGeometryFile
1085                         << "part" << nl
1086                         << setw(10) << ensightPatchi++ << nl
1087                         << patchName << nl
1088                         << "coordinates" << nl
1089                         << setw(10) << nfp.nPoints
1090                         << endl;
1092                     for (direction d=0; d<vector::nComponents; d++)
1093                     {
1094                         if (patchPointsPtr)
1095                         {
1096                             writePoints
1097                             (
1098                                 patchPoints.component(d),
1099                                 ensightGeometryFile
1100                             );
1101                         }
1103                         patchPointOffsets = 0;
1105                         forAll (patchProcessors, i)
1106                         {
1107                             if (patchProcessors[i] != 0)
1108                             {
1109                                 label slave = patchProcessors[i];
1110                                 IPstream fromSlave(Pstream::scheduled, slave);
1111                                 scalarField patchPointsComponent(fromSlave);
1113                                 writePoints
1114                                 (
1115                                     patchPointsComponent,
1116                                     ensightGeometryFile
1117                                 );
1119                                 if (i < Pstream::nProcs()-1)
1120                                 {
1121                                     patchPointOffsets[i+1] =
1122                                         patchPointOffsets[i]
1123                                       + patchPointsComponent.size();
1124                                 }
1125                             }
1126                             else
1127                             {
1128                                 if (i < Pstream::nProcs()-1)
1129                                 {
1130                                     patchPointOffsets[i+1] =
1131                                         patchPointOffsets[i]
1132                                       + patchPoints.size();
1133                                 }
1134                             }
1135                         }
1136                     }
1137                 }
1138                 else if (patchPointsPtr)
1139                 {
1140                     for (direction d=0; d<vector::nComponents; d++)
1141                     {
1142                         OPstream toMaster
1143                         (
1144                             Pstream::scheduled,
1145                             Pstream::masterNo()
1146                         );
1147                         toMaster<< patchPoints.component(d);
1148                     }
1149                 }
1151                 writeAllFacePrims
1152                 (
1153                     "tria3",
1154                     tris,
1155                     nfp.nTris,
1156                     patchFaces,
1157                     patchPointOffsets,
1158                     patchProcessors,
1159                     ensightGeometryFile
1160                 );
1162                 writeAllFacePrims
1163                 (
1164                     "quad4",
1165                     quads,
1166                     nfp.nQuads,
1167                     patchFaces,
1168                     patchPointOffsets,
1169                     patchProcessors,
1170                     ensightGeometryFile
1171                 );
1173                 writeAllFacePrims
1174                 (
1175                     "nsided",
1176                     polys,
1177                     nfp.nPolys,
1178                     patchFaces,
1179                     patchPointOffsets,
1180                     patchProcessors,
1181                     ensightGeometryFile
1182                 );
1183             }
1184         }
1185     }
1187     if (Pstream::master())
1188     {
1189         delete ensightGeometryFilePtr;
1190     }
1194 void Foam::ensightMesh::writeBinary
1196     const fileName& postProcPath,
1197     const word& prepend,
1198     const label timeIndex,
1199     Ostream& ensightCaseFile
1200 ) const
1202     //const Time& runTime = mesh.time();
1203     const pointField& points = mesh_.points();
1204     const cellList& cellFaces = mesh_.cells();
1205     const faceList& faces = mesh_.faces();
1206     const cellShapeList& cellShapes = mesh_.cellShapes();
1208     word timeFile = prepend;
1210     if (timeIndex == 0)
1211     {
1212         timeFile += "000.";
1213     }
1214     else if (mesh_.moving())
1215     {
1216         timeFile += itoa(timeIndex) + '.';
1217     }
1219     // set the filename of the ensight file
1220     fileName ensightGeometryFileName = timeFile + "mesh";
1222     std::ofstream *ensightGeometryFilePtr = NULL;
1224     if (Pstream::master())
1225     {
1226         ensightGeometryFilePtr = new std::ofstream
1227         (
1228             (postProcPath/ensightGeometryFileName).c_str(),
1229             ios_base::out | ios_base::binary | ios_base::trunc
1230         );
1231         // Check on file opened?
1232     }
1234     std::ofstream& ensightGeometryFile = *ensightGeometryFilePtr;
1236     if (Pstream::master())
1237     {
1238         writeEnsDataBinary("C binary",ensightGeometryFile);
1239         writeEnsDataBinary("OpenFOAM Geometry File",ensightGeometryFile);
1240         writeEnsDataBinary("Binary format",ensightGeometryFile);
1241         writeEnsDataBinary("node id assign",ensightGeometryFile);
1242         writeEnsDataBinary("element id assign",ensightGeometryFile);
1243     }
1245     labelList pointOffsets(Pstream::nProcs(), 0);
1247     if (!patchNames_.size())
1248     {
1249         label nPoints = points.size();
1250         Pstream::gather(nPoints, sumOp<label>());
1252         if (Pstream::master())
1253         {
1254             writeEnsDataBinary("part",ensightGeometryFile);
1255             writeEnsDataBinary(1,ensightGeometryFile);
1256             writeEnsDataBinary("FOAM cells",ensightGeometryFile);
1257             writeEnsDataBinary("coordinates",ensightGeometryFile);
1258             writeEnsDataBinary(nPoints,ensightGeometryFile);
1260             for (direction d=0; d<vector::nComponents; d++)
1261             {
1262                 //writePointsBinary(points.component(d), ensightGeometryFile);
1263                 writeEnsDataBinary(points.component(d), ensightGeometryFile);
1264                 pointOffsets[0] = points.size();
1266                 for (int slave=1; slave<Pstream::nProcs(); slave++)
1267                 {
1268                     IPstream fromSlave(Pstream::scheduled, slave);
1269                     scalarField pointsComponent(fromSlave);
1270                     //writePointsBinary(pointsComponent, ensightGeometryFile);
1271                     writeEnsDataBinary(pointsComponent, ensightGeometryFile);
1272                     pointOffsets[slave] =
1273                         pointOffsets[slave-1]
1274                       + pointsComponent.size();
1275                 }
1276             }
1277         }
1278         else
1279         {
1280             for (direction d=0; d<vector::nComponents; d++)
1281             {
1282                 OPstream toMaster(Pstream::scheduled, Pstream::masterNo());
1283                 toMaster<< points.component(d);
1284             }
1285         }
1287         writeAllPrimsBinary
1288         (
1289             "hexa8",
1290             meshCellSets_.nHexesWedges,
1291             map(cellShapes, meshCellSets_.hexes, meshCellSets_.wedges),
1292             pointOffsets,
1293             ensightGeometryFile
1294         );
1296         writeAllPrimsBinary
1297         (
1298             "penta6",
1299             meshCellSets_.nPrisms,
1300             map(cellShapes, meshCellSets_.prisms),
1301             pointOffsets,
1302             ensightGeometryFile
1303         );
1305         writeAllPrimsBinary
1306         (
1307             "pyramid5",
1308             meshCellSets_.nPyrs,
1309             map(cellShapes, meshCellSets_.pyrs),
1310             pointOffsets,
1311             ensightGeometryFile
1312         );
1314         writeAllPrimsBinary
1315         (
1316             "tetra4",
1317             meshCellSets_.nTets,
1318             map(cellShapes, meshCellSets_.tets),
1319             pointOffsets,
1320             ensightGeometryFile
1321         );
1323         if (meshCellSets_.nPolys)
1324         {
1325             if (Pstream::master())
1326             {
1327                 /*
1328                 ensightGeometryFile
1329                     << "nfaced" << nl
1330                     << setw(10) << meshCellSets_.nPolys << nl;
1331                 */
1332                 writePolysBinary
1333                 (
1334                     meshCellSets_.polys,
1335                     cellFaces,
1336                     faces,
1337                     0,
1338                     ensightGeometryFile
1339                 );
1341                 for (int slave=1; slave<Pstream::nProcs(); slave++)
1342                 {
1343                     IPstream fromSlave(Pstream::scheduled, slave);
1344                     labelList polys(fromSlave);
1345                     cellList cellFaces(fromSlave);
1346                     faceList faces(fromSlave);
1348                     writePolysBinary
1349                     (
1350                         polys,
1351                         cellFaces,
1352                         faces,
1353                         pointOffsets[slave-1],
1354                         ensightGeometryFile
1355                     );
1356                 }
1357             }
1358             else
1359             {
1360                 OPstream toMaster(Pstream::scheduled, Pstream::masterNo());
1361                 toMaster<< meshCellSets_.polys << cellFaces << faces;
1362             }
1363         }
1365     }
1367     label ensightPatchi = 2;
1369     label iCount = 0;
1371     forAllConstIter(HashTable<labelList>, allPatchNames_, iter)
1372     {
1373         iCount ++;
1374         const labelList& patchProcessors = iter();
1376         if (!patchNames_.size() || patchNames_.found(iter.key()))
1377         {
1378             const word& patchName = iter.key();
1379             const nFacePrimitives& nfp = nPatchPrims_.find(patchName)();
1381             const labelList *trisPtr = NULL;
1382             const labelList *quadsPtr = NULL;
1383             const labelList *polysPtr = NULL;
1385             const pointField *patchPointsPtr = NULL;
1386             const faceList *patchFacesPtr = NULL;
1388             if (patchIndices_.found(iter.key()))
1389             {
1390                 label patchi = patchIndices_.find(iter.key())();
1391                 const polyPatch& p = mesh_.boundaryMesh()[patchi];
1393                 trisPtr = &boundaryFaceSets_[patchi].tris;
1394                 quadsPtr = &boundaryFaceSets_[patchi].quads;
1395                 polysPtr = &boundaryFaceSets_[patchi].polys;
1397                 patchPointsPtr = &(p.localPoints());
1398                 patchFacesPtr = &(p.localFaces());
1399             }
1401             const labelList& tris = *trisPtr;
1402             const labelList& quads = *quadsPtr;
1403             const labelList& polys = *polysPtr;
1404             const pointField& patchPoints = *patchPointsPtr;
1405             const faceList& patchFaces = *patchFacesPtr;
1407             if (nfp.nTris || nfp.nQuads || nfp.nPolys)
1408             {
1409                 labelList patchPointOffsets(Pstream::nProcs(), 0);
1411                 if (Pstream::master())
1412                 {
1413                     writeEnsDataBinary("part",ensightGeometryFile);
1414                     writeEnsDataBinary(ensightPatchi++,ensightGeometryFile);
1415                     //writeEnsDataBinary(patchName.c_str(),ensightGeometryFile);
1416                     writeEnsDataBinary(iter.key().c_str(),ensightGeometryFile);
1417                     writeEnsDataBinary("coordinates",ensightGeometryFile);
1418                     writeEnsDataBinary(nfp.nPoints,ensightGeometryFile);
1420                     for (direction d=0; d<vector::nComponents; d++)
1421                     {
1422                         if (patchPointsPtr)
1423                         {
1424                             //writePointsBinary
1425                             writeEnsDataBinary
1426                             (
1427                                 patchPoints.component(d),
1428                                 ensightGeometryFile
1429                             );
1430                         }
1432                         patchPointOffsets = 0;
1435                         forAll (patchProcessors, i)
1436                         {
1437                             if (patchProcessors[i] != 0)
1438                             {
1439                                 label slave = patchProcessors[i];
1440                                 IPstream fromSlave(Pstream::scheduled, slave);
1441                                 scalarField patchPointsComponent(fromSlave);
1443                                 //writePointsBinary
1444                                 writeEnsDataBinary
1445                                 (
1446                                     patchPointsComponent,
1447                                     ensightGeometryFile
1448                                 );
1450                                 if (i < Pstream::nProcs()-1)
1451                                 {
1452                                     patchPointOffsets[i+1] =
1453                                         patchPointOffsets[i]
1454                                       + patchPointsComponent.size();
1455                                 }
1456                             }
1457                             else
1458                             {
1459                                 if (i < Pstream::nProcs()-1)
1460                                 {
1461                                     patchPointOffsets[i+1] =
1462                                         patchPointOffsets[i]
1463                                       + patchPoints.size();
1464                                 }
1465                             }
1466                         }
1467                     }
1468                 }
1469                 else if (patchPointsPtr)
1470                 {
1471                     for (direction d=0; d<vector::nComponents; d++)
1472                     {
1473                         OPstream toMaster
1474                         (
1475                             Pstream::scheduled,
1476                             Pstream::masterNo()
1477                         );
1478                         toMaster<< patchPoints.component(d);
1479                     }
1480                 }
1482                 writeAllFacePrimsBinary
1483                 (
1484                     "tria3",
1485                     tris,
1486                     nfp.nTris,
1487                     patchFaces,
1488                     patchPointOffsets,
1489                     patchProcessors,
1490                     ensightGeometryFile
1491                 );
1493                 writeAllFacePrimsBinary
1494                 (
1495                     "quad4",
1496                     quads,
1497                     nfp.nQuads,
1498                     patchFaces,
1499                     patchPointOffsets,
1500                     patchProcessors,
1501                     ensightGeometryFile
1502                 );
1504                 writeAllFacePrimsBinary
1505                 (
1506                     "nsided",
1507                     polys,
1508                     nfp.nPolys,
1509                     patchFaces,
1510                     patchPointOffsets,
1511                     patchProcessors,
1512                     ensightGeometryFile
1513                 );
1514             }
1515         }
1516     }
1519     if (Pstream::master())
1520     {
1521         delete ensightGeometryFilePtr;
1522     }
1525 // ************************************************************************* //