1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | foam-extend: Open Source CFD
4 \\ / O peration | Version: 3.2
5 \\ / A nd | Web: http://www.foam-extend.org
6 \\/ M anipulation | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
9 This file is part of foam-extend.
11 foam-extend is free software: you can redistribute it and/or modify it
12 under the terms of the GNU General Public License as published by the
13 Free Software Foundation, either version 3 of the License, or (at your
14 option) any later version.
16 foam-extend is distributed in the hope that it will be useful, but
17 WITHOUT ANY WARRANTY; without even the implied warranty of
18 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 General Public License for more details.
21 You should have received a copy of the GNU General Public License
22 along with foam-extend. If not, see <http://www.gnu.org/licenses/>.
24 \*---------------------------------------------------------------------------*/
26 #include "ensightField.H"
28 #include "volFields.H"
32 #include "ensightWriteBinary.H"
36 // * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
38 void writeData(const scalarField& sf, OFstream& ensightFile)
42 if (mag( sf[i] ) >= scalar(floatScalarVSMALL))
44 ensightFile << setw(12) << sf[i] << nl;
48 ensightFile << setw(12) << scalar(0) << nl;
57 const Field<Type>& vf,
62 scalarField mf(map.size());
66 mf[i] = component(vf[map[i]], cmpt);
76 const Field<Type>& vf,
77 const labelList& map1,
78 const labelList& map2,
82 scalarField mf(map1.size() + map2.size());
86 mf[i] = component(vf[map1[i]], cmpt);
89 label offset = map1.size();
93 mf[i + offset] = component(vf[map2[i]], cmpt);
104 const Field<Type>& vf,
105 const labelList& prims,
107 OFstream& ensightFile
112 if (Pstream::master())
114 ensightFile << key << nl;
116 for (direction cmpt=0; cmpt<pTraits<Type>::nComponents; cmpt++)
118 writeData(ensMap(vf, prims, cmpt), ensightFile);
120 for (int slave=1; slave<Pstream::nProcs(); slave++)
122 IPstream fromSlave(Pstream::scheduled, slave);
123 scalarField data(fromSlave);
124 writeData(data, ensightFile);
130 for (direction cmpt=0; cmpt<pTraits<Type>::nComponents; cmpt++)
132 OPstream toMaster(Pstream::scheduled, Pstream::masterNo());
133 toMaster<< ensMap(vf, prims, cmpt);
141 void writeAllDataBinary
144 const Field<Type>& vf,
145 const labelList& prims,
147 std::ofstream& ensightFile
152 if (Pstream::master())
154 writeEnsDataBinary(key,ensightFile);
156 for (direction cmpt=0; cmpt<pTraits<Type>::nComponents; cmpt++)
158 writeEnsDataBinary(ensMap(vf, prims, cmpt), ensightFile);
160 for (int slave=1; slave<Pstream::nProcs(); slave++)
162 IPstream fromSlave(Pstream::scheduled, slave);
163 scalarField data(fromSlave);
164 writeEnsDataBinary(data, ensightFile);
170 for (direction cmpt=0; cmpt<pTraits<Type>::nComponents; cmpt++)
172 OPstream toMaster(Pstream::scheduled, Pstream::masterNo());
173 toMaster<< ensMap(vf, prims, cmpt);
181 void writeAllFaceData
184 const labelList& prims,
186 const Field<Type>& pf,
187 const labelList& patchProcessors,
188 OFstream& ensightFile
193 if (Pstream::master())
195 ensightFile << key << nl;
197 for (direction cmpt=0; cmpt<pTraits<Type>::nComponents; cmpt++)
199 writeData(ensMap(pf, prims, cmpt), ensightFile);
201 forAll (patchProcessors, i)
203 if (patchProcessors[i] != 0)
205 label slave = patchProcessors[i];
206 IPstream fromSlave(Pstream::scheduled, slave);
207 scalarField pf(fromSlave);
209 writeData(pf, ensightFile);
216 for (direction cmpt=0; cmpt<pTraits<Type>::nComponents; cmpt++)
218 OPstream toMaster(Pstream::scheduled, Pstream::masterNo());
219 toMaster<< ensMap(pf, prims, cmpt);
227 void writeAllFaceDataBinary
230 const labelList& prims,
232 const Field<Type>& pf,
233 const labelList& patchProcessors,
234 std::ofstream& ensightFile
239 if (Pstream::master())
241 writeEnsDataBinary(key,ensightFile);
243 for (direction cmpt=0; cmpt<pTraits<Type>::nComponents; cmpt++)
245 writeEnsDataBinary(ensMap(pf, prims, cmpt), ensightFile);
247 forAll (patchProcessors, i)
249 if (patchProcessors[i] != 0)
251 label slave = patchProcessors[i];
252 IPstream fromSlave(Pstream::scheduled, slave);
253 scalarField pf(fromSlave);
255 writeEnsDataBinary(pf, ensightFile);
262 for (direction cmpt=0; cmpt<pTraits<Type>::nComponents; cmpt++)
264 OPstream toMaster(Pstream::scheduled, Pstream::masterNo());
265 toMaster<< ensMap(pf, prims, cmpt);
275 const Foam::Field<Type>& pf,
276 const Foam::label patchi,
277 const Foam::label ensightPatchI,
278 const Foam::faceSets& boundaryFaceSet,
279 const Foam::ensightMesh::nFacePrimitives& nfp,
280 const Foam::labelList& patchProcessors,
281 Foam::OFstream& ensightFile
284 if (nfp.nTris || nfp.nQuads || nfp.nPolys)
286 if (Pstream::master())
290 << setw(10) << ensightPatchI << nl;
296 boundaryFaceSet.tris,
306 boundaryFaceSet.quads,
316 boundaryFaceSet.polys,
333 bool writePatchFieldBinary
335 const Foam::Field<Type>& pf,
336 const Foam::label patchi,
337 const Foam::label ensightPatchI,
338 const Foam::faceSets& boundaryFaceSet,
339 const Foam::ensightMesh::nFacePrimitives& nfp,
340 const Foam::labelList& patchProcessors,
341 std::ofstream& ensightFile
344 if (nfp.nTris || nfp.nQuads || nfp.nPolys)
346 if (Pstream::master())
348 writeEnsDataBinary("part",ensightFile);
349 writeEnsDataBinary(ensightPatchI,ensightFile);
352 writeAllFaceDataBinary
355 boundaryFaceSet.tris,
362 writeAllFaceDataBinary
365 boundaryFaceSet.quads,
372 writeAllFaceDataBinary
375 boundaryFaceSet.polys,
394 const Foam::word& fieldName,
395 const Foam::Field<Type>& pf,
396 const Foam::word& patchName,
397 const Foam::ensightMesh& eMesh,
398 const Foam::fileName& postProcPath,
399 const Foam::word& prepend,
400 const Foam::label timeIndex,
401 Foam::Ostream& ensightCaseFile
404 const Time& runTime = eMesh.mesh().time();
406 const List<faceSets>& boundaryFaceSets = eMesh.boundaryFaceSets();
407 const wordList& allPatchNames = eMesh.allPatchNames();
408 const List<labelList>& allPatchProcs = eMesh.allPatchProcs();
409 const HashTable<ensightMesh::nFacePrimitives>&
410 nPatchPrims = eMesh.nPatchPrims();
412 label ensightPatchI = eMesh.patchPartOffset();
416 forAll(allPatchNames, i)
418 if (allPatchNames[i] == patchName)
427 const labelList& patchProcessors = allPatchProcs[patchi];
429 word pfName = patchName + '.' + fieldName;
431 word timeFile = prepend + itoa(timeIndex);
433 OFstream *ensightFilePtr = NULL;
434 if (Pstream::master())
438 ensightCaseFile.setf(ios_base::left);
441 << pTraits<Type>::typeName
442 << " per element: 1 "
443 << setw(15) << pfName
444 << (' ' + prepend + "***." + pfName).c_str()
448 // set the filename of the ensight file
449 fileName ensightFileName(timeFile + "." + pfName);
450 ensightFilePtr = new OFstream
452 postProcPath/ensightFileName,
453 ios_base::out|ios_base::trunc,
454 runTime.writeFormat(),
455 runTime.writeVersion(),
456 runTime.writeCompression()
460 OFstream& ensightFile = *ensightFilePtr;
462 if (Pstream::master())
464 ensightFile << pTraits<Type>::typeName << nl;
474 boundaryFaceSets[patchi],
475 nPatchPrims.find(patchName)(),
482 faceSets nullFaceSets;
490 nPatchPrims.find(patchName)(),
496 if (Pstream::master())
498 delete ensightFilePtr;
504 void ensightFieldAscii
506 const Foam::IOobject& fieldObject,
507 const Foam::ensightMesh& eMesh,
508 const Foam::fileName& postProcPath,
509 const Foam::word& prepend,
510 const Foam::label timeIndex,
511 Foam::Ostream& ensightCaseFile
514 Info<< "Converting field " << fieldObject.name() << endl;
516 word timeFile = prepend + itoa(timeIndex);
518 const fvMesh& mesh = eMesh.mesh();
519 const Time& runTime = mesh.time();
521 const cellSets& meshCellSets = eMesh.meshCellSets();
522 const List<faceSets>& boundaryFaceSets = eMesh.boundaryFaceSets();
523 const wordList& allPatchNames = eMesh.allPatchNames();
524 const List<labelList>& allPatchProcs = eMesh.allPatchProcs();
525 const wordHashSet& patchNames = eMesh.patchNames();
526 const HashTable<ensightMesh::nFacePrimitives>&
527 nPatchPrims = eMesh.nPatchPrims();
529 const labelList& tets = meshCellSets.tets;
530 const labelList& pyrs = meshCellSets.pyrs;
531 const labelList& prisms = meshCellSets.prisms;
532 const labelList& wedges = meshCellSets.wedges;
533 const labelList& hexes = meshCellSets.hexes;
534 const labelList& polys = meshCellSets.polys;
536 OFstream *ensightFilePtr = NULL;
537 if (Pstream::master())
539 // set the filename of the ensight file
540 fileName ensightFileName(timeFile + "." + fieldObject.name());
541 ensightFilePtr = new OFstream
543 postProcPath/ensightFileName,
544 ios_base::out|ios_base::trunc,
545 runTime.writeFormat(),
546 runTime.writeVersion(),
547 runTime.writeCompression()
551 OFstream& ensightFile = *ensightFilePtr;
553 GeometricField<Type, fvPatchField, volMesh> vf(fieldObject, mesh);
555 if (patchNames.empty())
557 if (Pstream::master())
561 ensightCaseFile.setf(ios_base::left);
564 << pTraits<Type>::typeName
565 << " per element: 1 "
566 << setw(15) << vf.name()
567 << (' ' + prepend + "***." + vf.name()).c_str()
572 << pTraits<Type>::typeName << nl
574 << setw(10) << 1 << nl;
576 ensightFile.setf(ios_base::scientific, ios_base::floatfield);
577 ensightFile.precision(5);
580 if (meshCellSets.nHexesWedges)
582 if (Pstream::master())
584 ensightFile << "hexa8" << nl;
586 for (direction cmpt=0; cmpt<pTraits<Type>::nComponents; cmpt++)
590 ensMap(vf, hexes, wedges, cmpt),
594 for (int slave=1; slave<Pstream::nProcs(); slave++)
596 IPstream fromSlave(Pstream::scheduled, slave);
597 scalarField data(fromSlave);
598 writeData(data, ensightFile);
604 for (direction cmpt=0; cmpt<pTraits<Type>::nComponents; cmpt++)
606 OPstream toMaster(Pstream::scheduled, Pstream::masterNo());
607 toMaster<< ensMap(vf, hexes, wedges, cmpt);
612 writeAllData("penta6", vf, prisms, meshCellSets.nPrisms, ensightFile);
613 writeAllData("pyramid5", vf, pyrs, meshCellSets.nPyrs, ensightFile);
614 writeAllData("tetra4", vf, tets, meshCellSets.nTets, ensightFile);
615 writeAllData("nfaced", vf, polys, meshCellSets.nPolys, ensightFile);
618 label ensightPatchI = eMesh.patchPartOffset();
620 forAll(allPatchNames, patchi)
622 const word& patchName = allPatchNames[patchi];
623 const labelList& patchProcessors = allPatchProcs[patchi];
625 if (patchNames.empty() || patchNames.found(patchName))
627 if (mesh.boundary()[patchi].size())
633 vf.boundaryField()[patchi],
636 boundaryFaceSets[patchi],
637 nPatchPrims.find(patchName)(),
647 else if (Pstream::master())
649 faceSets nullFaceSet;
659 nPatchPrims.find(patchName)(),
671 if (Pstream::master())
673 delete ensightFilePtr;
679 void ensightFieldBinary
681 const Foam::IOobject& fieldObject,
682 const Foam::ensightMesh& eMesh,
683 const Foam::fileName& postProcPath,
684 const Foam::word& prepend,
685 const Foam::label timeIndex,
686 Foam::Ostream& ensightCaseFile
689 Info<< "Converting field (binary) " << fieldObject.name() << endl;
691 word timeFile = prepend + itoa(timeIndex);
693 const fvMesh& mesh = eMesh.mesh();
694 //const Time& runTime = mesh.time();
696 const cellSets& meshCellSets = eMesh.meshCellSets();
697 const List<faceSets>& boundaryFaceSets = eMesh.boundaryFaceSets();
698 const wordList& allPatchNames = eMesh.allPatchNames();
699 const List<labelList>& allPatchProcs = eMesh.allPatchProcs();
700 const wordHashSet& patchNames = eMesh.patchNames();
701 const HashTable<ensightMesh::nFacePrimitives>&
702 nPatchPrims = eMesh.nPatchPrims();
704 const labelList& tets = meshCellSets.tets;
705 const labelList& pyrs = meshCellSets.pyrs;
706 const labelList& prisms = meshCellSets.prisms;
707 const labelList& wedges = meshCellSets.wedges;
708 const labelList& hexes = meshCellSets.hexes;
709 const labelList& polys = meshCellSets.polys;
711 std::ofstream *ensightFilePtr = NULL;
712 if (Pstream::master())
714 // set the filename of the ensight file
715 fileName ensightFileName(timeFile + "." + fieldObject.name());
716 ensightFilePtr = new std::ofstream
718 (postProcPath/ensightFileName).c_str(),
719 ios_base::out | ios_base::binary | ios_base::trunc
721 // Check on file opened?
724 std::ofstream& ensightFile = *ensightFilePtr;
726 GeometricField<Type, fvPatchField, volMesh> vf(fieldObject, mesh);
728 if (patchNames.empty())
730 if (Pstream::master())
734 ensightCaseFile.setf(ios_base::left);
737 << pTraits<Type>::typeName
738 << " per element: 1 "
739 << setw(15) << vf.name()
740 << (' ' + prepend + "***." + vf.name()).c_str()
744 writeEnsDataBinary(pTraits<Type>::typeName,ensightFile);
745 writeEnsDataBinary("part",ensightFile);
746 writeEnsDataBinary(1,ensightFile);
749 if (meshCellSets.nHexesWedges)
751 if (Pstream::master())
753 writeEnsDataBinary("hexa8",ensightFile);
755 for (direction cmpt=0; cmpt<pTraits<Type>::nComponents; cmpt++)
759 ensMap(vf, hexes, wedges, cmpt),
763 for (int slave=1; slave<Pstream::nProcs(); slave++)
765 IPstream fromSlave(Pstream::scheduled, slave);
766 scalarField data(fromSlave);
767 writeEnsDataBinary(data, ensightFile);
773 for (direction cmpt=0; cmpt<pTraits<Type>::nComponents; cmpt++)
775 OPstream toMaster(Pstream::scheduled, Pstream::masterNo());
776 toMaster<< ensMap(vf, hexes, wedges, cmpt);
786 meshCellSets.nPrisms,
818 label ensightPatchI = eMesh.patchPartOffset();
820 forAll(allPatchNames, patchi)
822 const word& patchName = allPatchNames[patchi];
823 const labelList& patchProcessors = allPatchProcs[patchi];
825 if (patchNames.empty() || patchNames.found(patchName))
827 if (mesh.boundary()[patchi].size())
831 writePatchFieldBinary
833 vf.boundaryField()[patchi],
836 boundaryFaceSets[patchi],
837 nPatchPrims.find(patchName)(),
847 else if (Pstream::master())
849 faceSets nullFaceSet;
853 writePatchFieldBinary
859 nPatchPrims.find(patchName)(),
871 if (Pstream::master())
881 const Foam::IOobject& fieldObject,
882 const Foam::ensightMesh& eMesh,
883 const Foam::fileName& postProcPath,
884 const Foam::word& prepend,
885 const Foam::label timeIndex,
887 Foam::Ostream& ensightCaseFile
892 ensightFieldBinary<Type>
904 ensightFieldAscii<Type>
917 // ************************************************************************* //