1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright held by original author
7 -------------------------------------------------------------------------------
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
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 "ensightField.H"
29 #include "volFields.H"
33 #include "ensightWriteBinary.H"
37 // * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
39 void writeData(const scalarField& sf, OFstream& ensightFile)
43 if (mag( sf[i] ) >= scalar(floatScalarVSMALL))
45 ensightFile << setw(12) << sf[i] << nl;
49 ensightFile << setw(12) << scalar(0) << nl;
58 const Field<Type>& vf,
63 scalarField mf(map.size());
67 mf[i] = component(vf[map[i]], cmpt);
77 const Field<Type>& vf,
78 const labelList& map1,
79 const labelList& map2,
83 scalarField mf(map1.size() + map2.size());
87 mf[i] = component(vf[map1[i]], cmpt);
90 label offset = map1.size();
94 mf[i + offset] = component(vf[map2[i]], cmpt);
105 const Field<Type>& vf,
106 const labelList& prims,
108 OFstream& ensightFile
113 if (Pstream::master())
115 ensightFile << key << nl;
117 for (direction cmpt=0; cmpt<pTraits<Type>::nComponents; cmpt++)
119 writeData(map(vf, prims, cmpt), ensightFile);
121 for (int slave=1; slave<Pstream::nProcs(); slave++)
123 IPstream fromSlave(Pstream::scheduled, slave);
124 scalarField data(fromSlave);
125 writeData(data, ensightFile);
131 for (direction cmpt=0; cmpt<pTraits<Type>::nComponents; cmpt++)
133 OPstream toMaster(Pstream::scheduled, Pstream::masterNo());
134 toMaster<< map(vf, prims, cmpt);
142 void writeAllDataBinary
145 const Field<Type>& vf,
146 const labelList& prims,
148 std::ofstream& ensightFile
153 if (Pstream::master())
155 writeEnsDataBinary(key,ensightFile);
157 for (direction cmpt=0; cmpt<pTraits<Type>::nComponents; cmpt++)
159 writeEnsDataBinary(map(vf, prims, cmpt), ensightFile);
161 for (int slave=1; slave<Pstream::nProcs(); slave++)
163 IPstream fromSlave(Pstream::scheduled, slave);
164 scalarField data(fromSlave);
165 writeEnsDataBinary(data, ensightFile);
171 for (direction cmpt=0; cmpt<pTraits<Type>::nComponents; cmpt++)
173 OPstream toMaster(Pstream::scheduled, Pstream::masterNo());
174 toMaster<< map(vf, prims, cmpt);
182 void writeAllFaceData
185 const labelList& prims,
187 const Field<Type>& pf,
188 const labelList& patchProcessors,
189 OFstream& ensightFile
194 if (Pstream::master())
196 ensightFile << key << nl;
198 for (direction cmpt=0; cmpt<pTraits<Type>::nComponents; cmpt++)
200 writeData(map(pf, prims, cmpt), ensightFile);
202 forAll (patchProcessors, i)
204 if (patchProcessors[i] != 0)
206 label slave = patchProcessors[i];
207 IPstream fromSlave(Pstream::scheduled, slave);
208 scalarField pf(fromSlave);
210 writeData(pf, ensightFile);
217 for (direction cmpt=0; cmpt<pTraits<Type>::nComponents; cmpt++)
219 OPstream toMaster(Pstream::scheduled, Pstream::masterNo());
220 toMaster<< map(pf, prims, cmpt);
228 void writeAllFaceDataBinary
231 const labelList& prims,
233 const Field<Type>& pf,
234 const labelList& patchProcessors,
235 std::ofstream& ensightFile
240 if (Pstream::master())
242 writeEnsDataBinary(key,ensightFile);
244 for (direction cmpt=0; cmpt<pTraits<Type>::nComponents; cmpt++)
246 writeEnsDataBinary(map(pf, prims, cmpt), ensightFile);
248 forAll (patchProcessors, i)
250 if (patchProcessors[i] != 0)
252 label slave = patchProcessors[i];
253 IPstream fromSlave(Pstream::scheduled, slave);
254 scalarField pf(fromSlave);
256 writeEnsDataBinary(pf, ensightFile);
263 for (direction cmpt=0; cmpt<pTraits<Type>::nComponents; cmpt++)
265 OPstream toMaster(Pstream::scheduled, Pstream::masterNo());
266 toMaster<< map(pf, prims, cmpt);
276 const Foam::Field<Type>& pf,
277 const Foam::label patchi,
278 const Foam::label ensightPatchI,
279 const Foam::faceSets& boundaryFaceSet,
280 const Foam::ensightMesh::nFacePrimitives& nfp,
281 const Foam::labelList& patchProcessors,
282 Foam::OFstream& ensightFile
285 if (nfp.nTris || nfp.nQuads || nfp.nPolys)
287 if (Pstream::master())
291 << setw(10) << ensightPatchI << nl;
297 boundaryFaceSet.tris,
307 boundaryFaceSet.quads,
317 boundaryFaceSet.polys,
334 bool writePatchFieldBinary
336 const Foam::Field<Type>& pf,
337 const Foam::label patchi,
338 const Foam::label ensightPatchI,
339 const Foam::faceSets& boundaryFaceSet,
340 const Foam::ensightMesh::nFacePrimitives& nfp,
341 const Foam::labelList& patchProcessors,
342 std::ofstream& ensightFile
345 if (nfp.nTris || nfp.nQuads || nfp.nPolys)
347 if (Pstream::master())
349 writeEnsDataBinary("part",ensightFile);
350 writeEnsDataBinary(ensightPatchI,ensightFile);
353 writeAllFaceDataBinary
356 boundaryFaceSet.tris,
363 writeAllFaceDataBinary
366 boundaryFaceSet.quads,
373 writeAllFaceDataBinary
376 boundaryFaceSet.polys,
395 const Foam::word& fieldName,
396 const Foam::Field<Type>& pf,
397 const Foam::word& patchName,
398 const Foam::ensightMesh& eMesh,
399 const Foam::fileName& postProcPath,
400 const Foam::word& prepend,
401 const Foam::label timeIndex,
402 Foam::Ostream& ensightCaseFile
405 const Time& runTime = eMesh.mesh().time();
407 const List<faceSets>& boundaryFaceSets = eMesh.boundaryFaceSets();
408 const wordList& allPatchNames = eMesh.allPatchNames();
409 const List<labelList>& allPatchProcs = eMesh.allPatchProcs();
410 const HashTable<ensightMesh::nFacePrimitives>&
411 nPatchPrims = eMesh.nPatchPrims();
413 label ensightPatchI = eMesh.patchPartOffset();
417 forAll(allPatchNames, i)
419 if (allPatchNames[i] == patchName)
428 const labelList& patchProcessors = allPatchProcs[patchi];
430 word pfName = patchName + '.' + fieldName;
432 word timeFile = prepend + itoa(timeIndex);
434 OFstream *ensightFilePtr = NULL;
435 if (Pstream::master())
439 ensightCaseFile.setf(ios_base::left);
442 << pTraits<Type>::typeName
443 << " per element: 1 "
444 << setw(15) << pfName
445 << (' ' + prepend + "***." + pfName).c_str()
449 // set the filename of the ensight file
450 fileName ensightFileName(timeFile + "." + pfName);
451 ensightFilePtr = new OFstream
453 postProcPath/ensightFileName,
454 ios_base::out|ios_base::trunc,
455 runTime.writeFormat(),
456 runTime.writeVersion(),
457 runTime.writeCompression()
461 OFstream& ensightFile = *ensightFilePtr;
463 if (Pstream::master())
465 ensightFile << pTraits<Type>::typeName << nl;
475 boundaryFaceSets[patchi],
476 nPatchPrims.find(patchName)(),
483 faceSets nullFaceSets;
491 nPatchPrims.find(patchName)(),
497 if (Pstream::master())
499 delete ensightFilePtr;
505 void ensightFieldAscii
507 const Foam::IOobject& fieldObject,
508 const Foam::ensightMesh& eMesh,
509 const Foam::fileName& postProcPath,
510 const Foam::word& prepend,
511 const Foam::label timeIndex,
512 Foam::Ostream& ensightCaseFile
515 Info<< "Converting field " << fieldObject.name() << endl;
517 word timeFile = prepend + itoa(timeIndex);
519 const fvMesh& mesh = eMesh.mesh();
520 const Time& runTime = mesh.time();
522 const cellSets& meshCellSets = eMesh.meshCellSets();
523 const List<faceSets>& boundaryFaceSets = eMesh.boundaryFaceSets();
524 const wordList& allPatchNames = eMesh.allPatchNames();
525 const List<labelList>& allPatchProcs = eMesh.allPatchProcs();
526 const wordHashSet& patchNames = eMesh.patchNames();
527 const HashTable<ensightMesh::nFacePrimitives>&
528 nPatchPrims = eMesh.nPatchPrims();
530 const labelList& tets = meshCellSets.tets;
531 const labelList& pyrs = meshCellSets.pyrs;
532 const labelList& prisms = meshCellSets.prisms;
533 const labelList& wedges = meshCellSets.wedges;
534 const labelList& hexes = meshCellSets.hexes;
535 const labelList& polys = meshCellSets.polys;
537 OFstream *ensightFilePtr = NULL;
538 if (Pstream::master())
540 // set the filename of the ensight file
541 fileName ensightFileName(timeFile + "." + fieldObject.name());
542 ensightFilePtr = new OFstream
544 postProcPath/ensightFileName,
545 ios_base::out|ios_base::trunc,
546 runTime.writeFormat(),
547 runTime.writeVersion(),
548 runTime.writeCompression()
552 OFstream& ensightFile = *ensightFilePtr;
554 GeometricField<Type, fvPatchField, volMesh> vf(fieldObject, mesh);
556 if (patchNames.empty())
558 if (Pstream::master())
562 ensightCaseFile.setf(ios_base::left);
565 << pTraits<Type>::typeName
566 << " per element: 1 "
567 << setw(15) << vf.name()
568 << (' ' + prepend + "***." + vf.name()).c_str()
573 << pTraits<Type>::typeName << nl
575 << setw(10) << 1 << nl;
577 ensightFile.setf(ios_base::scientific, ios_base::floatfield);
578 ensightFile.precision(5);
581 if (meshCellSets.nHexesWedges)
583 if (Pstream::master())
585 ensightFile << "hexa8" << nl;
587 for (direction cmpt=0; cmpt<pTraits<Type>::nComponents; cmpt++)
591 map(vf, hexes, wedges, cmpt),
595 for (int slave=1; slave<Pstream::nProcs(); slave++)
597 IPstream fromSlave(Pstream::scheduled, slave);
598 scalarField data(fromSlave);
599 writeData(data, ensightFile);
605 for (direction cmpt=0; cmpt<pTraits<Type>::nComponents; cmpt++)
607 OPstream toMaster(Pstream::scheduled, Pstream::masterNo());
608 toMaster<< map(vf, hexes, wedges, cmpt);
613 writeAllData("penta6", vf, prisms, meshCellSets.nPrisms, ensightFile);
614 writeAllData("pyramid5", vf, pyrs, meshCellSets.nPyrs, ensightFile);
615 writeAllData("tetra4", vf, tets, meshCellSets.nTets, ensightFile);
616 writeAllData("nfaced", vf, polys, meshCellSets.nPolys, ensightFile);
619 label ensightPatchI = eMesh.patchPartOffset();
621 forAll(allPatchNames, patchi)
623 const word& patchName = allPatchNames[patchi];
624 const labelList& patchProcessors = allPatchProcs[patchi];
626 if (patchNames.empty() || patchNames.found(patchName))
628 if (mesh.boundary()[patchi].size())
634 vf.boundaryField()[patchi],
637 boundaryFaceSets[patchi],
638 nPatchPrims.find(patchName)(),
648 else if (Pstream::master())
650 faceSets nullFaceSet;
660 nPatchPrims.find(patchName)(),
672 if (Pstream::master())
674 delete ensightFilePtr;
680 void ensightFieldBinary
682 const Foam::IOobject& fieldObject,
683 const Foam::ensightMesh& eMesh,
684 const Foam::fileName& postProcPath,
685 const Foam::word& prepend,
686 const Foam::label timeIndex,
687 Foam::Ostream& ensightCaseFile
690 Info<< "Converting field (binary) " << fieldObject.name() << endl;
692 word timeFile = prepend + itoa(timeIndex);
694 const fvMesh& mesh = eMesh.mesh();
695 //const Time& runTime = mesh.time();
697 const cellSets& meshCellSets = eMesh.meshCellSets();
698 const List<faceSets>& boundaryFaceSets = eMesh.boundaryFaceSets();
699 const wordList& allPatchNames = eMesh.allPatchNames();
700 const List<labelList>& allPatchProcs = eMesh.allPatchProcs();
701 const wordHashSet& patchNames = eMesh.patchNames();
702 const HashTable<ensightMesh::nFacePrimitives>&
703 nPatchPrims = eMesh.nPatchPrims();
705 const labelList& tets = meshCellSets.tets;
706 const labelList& pyrs = meshCellSets.pyrs;
707 const labelList& prisms = meshCellSets.prisms;
708 const labelList& wedges = meshCellSets.wedges;
709 const labelList& hexes = meshCellSets.hexes;
710 const labelList& polys = meshCellSets.polys;
712 std::ofstream *ensightFilePtr = NULL;
713 if (Pstream::master())
715 // set the filename of the ensight file
716 fileName ensightFileName(timeFile + "." + fieldObject.name());
717 ensightFilePtr = new std::ofstream
719 (postProcPath/ensightFileName).c_str(),
720 ios_base::out | ios_base::binary | ios_base::trunc
722 // Check on file opened?
725 std::ofstream& ensightFile = *ensightFilePtr;
727 GeometricField<Type, fvPatchField, volMesh> vf(fieldObject, mesh);
729 if (patchNames.empty())
731 if (Pstream::master())
735 ensightCaseFile.setf(ios_base::left);
738 << pTraits<Type>::typeName
739 << " per element: 1 "
740 << setw(15) << vf.name()
741 << (' ' + prepend + "***." + vf.name()).c_str()
745 writeEnsDataBinary(pTraits<Type>::typeName,ensightFile);
746 writeEnsDataBinary("part",ensightFile);
747 writeEnsDataBinary(1,ensightFile);
750 if (meshCellSets.nHexesWedges)
752 if (Pstream::master())
754 writeEnsDataBinary("hexa8",ensightFile);
756 for (direction cmpt=0; cmpt<pTraits<Type>::nComponents; cmpt++)
760 map(vf, hexes, wedges, cmpt),
764 for (int slave=1; slave<Pstream::nProcs(); slave++)
766 IPstream fromSlave(Pstream::scheduled, slave);
767 scalarField data(fromSlave);
768 writeEnsDataBinary(data, ensightFile);
774 for (direction cmpt=0; cmpt<pTraits<Type>::nComponents; cmpt++)
776 OPstream toMaster(Pstream::scheduled, Pstream::masterNo());
777 toMaster<< map(vf, hexes, wedges, cmpt);
787 meshCellSets.nPrisms,
819 label ensightPatchI = eMesh.patchPartOffset();
821 forAll(allPatchNames, patchi)
823 const word& patchName = allPatchNames[patchi];
824 const labelList& patchProcessors = allPatchProcs[patchi];
826 if (patchNames.empty() || patchNames.found(patchName))
828 if (mesh.boundary()[patchi].size())
832 writePatchFieldBinary
834 vf.boundaryField()[patchi],
837 boundaryFaceSets[patchi],
838 nPatchPrims.find(patchName)(),
848 else if (Pstream::master())
850 faceSets nullFaceSet;
854 writePatchFieldBinary
860 nPatchPrims.find(patchName)(),
872 if (Pstream::master())
882 const Foam::IOobject& fieldObject,
883 const Foam::ensightMesh& eMesh,
884 const Foam::fileName& postProcPath,
885 const Foam::word& prepend,
886 const Foam::label timeIndex,
888 Foam::Ostream& ensightCaseFile
893 ensightFieldBinary<Type>
905 ensightFieldAscii<Type>
918 // ************************************************************************* //