ensightField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
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.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "ensightField.H"
27 #include "fvMesh.H"
28 #include "volFields.H"
29 #include "OFstream.H"
30 #include "IOmanip.H"
31 #include "itoa.H"
32 #include "volPointInterpolation.H"
33 #include "ensightBinaryStream.H"
34 #include "ensightAsciiStream.H"
35 #include "globalIndex.H"
36 #include "ensightPTraits.H"
38 
39 using namespace Foam;
40 
41 // * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
42 
43 template<class Type>
46 (
47  const fvMeshSubset& meshSubsetter,
49 )
50 {
51  if (meshSubsetter.hasSubMesh())
52  {
54  (
55  meshSubsetter.interpolate(vf)
56  );
57  tfld().checkOut();
58  tfld().rename(vf.name());
59  return tfld;
60  }
61  else
62  {
63  return vf;
64  }
65 }
66 
67 
68 template<class Type>
71 (
72  const fvMeshSubset& meshSubsetter,
73  const typename GeometricField
74  <
75  Type,
77  volMesh
78  >::DimensionedInternalField& df
79 )
80 {
81  // Construct volField (with zeroGradient) from dimensioned field
82 
83  IOobject io(df);
84  io.readOpt() = IOobject::NO_READ;
85 
87  (
89  (
90  io,
91  df.mesh(),
92  df.dimensions(),
94  )
95  );
96  tvf().internalField() = df;
97  tvf().correctBoundaryConditions();
99 
100  if (meshSubsetter.hasSubMesh())
101  {
103  (
104  meshSubsetter.interpolate(vf)
105  );
106  tfld().checkOut();
107  tfld().rename(vf.name());
108  return tfld;
109  }
110  else
111  {
112  return tvf;
113  }
114 }
115 
116 
117 //template<class Container>
118 //void readAndConvertField
119 //(
120 // const fvMeshSubset& meshSubsetter,
121 // const IOobject& io,
122 // const fvMesh& mesh,
123 // const ensightMesh& eMesh,
124 // const fileName& postProcPath,
125 // const word& prepend,
126 // const label timeIndex,
127 // const bool binary,
128 // const bool nodeValues,
129 // Ostream& ensightCaseFile
130 //)
131 //{
132 // Container fld(io, mesh);
133 // ensightField<typename Container::value_type>
134 // (
135 // volField<typename Container::value_type>(meshSubsetter, fld),
136 // eMesh,
137 // postProcPath,
138 // prepend,
139 // timeIndex,
140 // binary,
141 // nodeValues,
142 // ensightCaseFile
143 // );
144 //}
145 
146 
147 template<class Type>
148 Field<Type> map
149 (
150  const Field<Type>& vf,
151  const labelList& map1,
152  const labelList& map2
153 )
154 {
155  Field<Type> mf(map1.size() + map2.size());
156 
157  forAll(map1, i)
158  {
159  mf[i] = vf[map1[i]];
160  }
161 
162  label offset = map1.size();
163 
164  forAll(map2, i)
165  {
166  mf[i + offset] = vf[map2[i]];
167  }
168 
169  return mf;
170 }
171 
172 
173 template<class Type>
174 void writeField
175 (
176  const char* key,
177  const Field<Type>& vf,
179 )
180 {
181  if (returnReduce(vf.size(), sumOp<label>()) > 0)
182  {
183  if (Pstream::master())
184  {
185  ensightFile.write(key);
186 
187  for (direction i=0; i < pTraits<Type>::nComponents; ++i)
188  {
190 
191  ensightFile.write(vf.component(cmpt));
192 
193  for (int slave=1; slave<Pstream::nProcs(); slave++)
194  {
195  IPstream fromSlave(Pstream::scheduled, slave);
196  scalarField slaveData(fromSlave);
197  ensightFile.write(slaveData);
198  }
199  }
200  }
201  else
202  {
203  for (direction i=0; i < pTraits<Type>::nComponents; ++i)
204  {
206 
208  toMaster<< vf.component(cmpt);
209  }
210  }
211  }
212 }
213 
214 
215 template<class Type>
216 bool writePatchField
217 (
218  const Field<Type>& pf,
219  const label patchi,
220  const label ensightPatchI,
221  const faceSets& boundaryFaceSet,
222  const ensightMesh::nFacePrimitives& nfp,
224 )
225 {
226  if (nfp.nTris || nfp.nQuads || nfp.nPolys)
227  {
228  if (Pstream::master())
229  {
230  ensightFile.writePartHeader(ensightPatchI);
231  }
232 
233  writeField
234  (
235  "tria3",
236  Field<Type>(pf, boundaryFaceSet.tris),
238  );
239 
240  writeField
241  (
242  "quad4",
243  Field<Type>(pf, boundaryFaceSet.quads),
245  );
246 
247  writeField
248  (
249  "nsided",
250  Field<Type>(pf, boundaryFaceSet.polys),
252  );
253 
254  return true;
255  }
256  else
257  {
258  return false;
259  }
260 }
261 
262 
263 template<class Type>
264 void writePatchField
265 (
266  const word& fieldName,
267  const Field<Type>& pf,
268  const word& patchName,
269  const ensightMesh& eMesh,
270  const fileName& postProcPath,
271  const word& prepend,
272  const label timeIndex,
273  const bool binary,
274  Ostream& ensightCaseFile
275 )
276 {
277  const Time& runTime = eMesh.mesh().time();
278 
279  const List<faceSets>& boundaryFaceSets = eMesh.boundaryFaceSets();
280  const wordList& allPatchNames = eMesh.allPatchNames();
282  nPatchPrims = eMesh.nPatchPrims();
283 
284  label ensightPatchI = eMesh.patchPartOffset();
285 
286  label patchi = -1;
287 
288  forAll(allPatchNames, i)
289  {
290  if (allPatchNames[i] == patchName)
291  {
292  patchi = i;
293  break;
294  }
295  ensightPatchI++;
296  }
297 
298 
299  word pfName = patchName + '.' + fieldName;
300 
301  word timeFile = prepend + itoa(timeIndex);
302 
303  ensightStream* ensightFilePtr = NULL;
304  if (Pstream::master())
305  {
306  if (timeIndex == 0)
307  {
308  ensightCaseFile.setf(ios_base::left);
309 
310  ensightCaseFile
312  << " per element: 1 "
313  << setw(15) << pfName
314  << (' ' + prepend + "****." + pfName).c_str()
315  << nl;
316  }
317 
318  // set the filename of the ensight file
319  fileName ensightFileName(timeFile + "." + pfName);
320 
321  if (binary)
322  {
323  ensightFilePtr = new ensightBinaryStream
324  (
325  postProcPath/ensightFileName,
326  runTime
327  );
328  }
329  else
330  {
331  ensightFilePtr = new ensightAsciiStream
332  (
333  postProcPath/ensightFileName,
334  runTime
335  );
336  }
337  }
338 
339  ensightStream& ensightFile = *ensightFilePtr;
340 
341  if (Pstream::master())
342  {
344  }
345 
346  if (patchi >= 0)
347  {
349  (
350  pf,
351  patchi,
352  ensightPatchI,
353  boundaryFaceSets[patchi],
354  nPatchPrims.find(patchName)(),
356  );
357  }
358  else
359  {
360  faceSets nullFaceSets;
361 
363  (
364  Field<Type>(),
365  -1,
366  ensightPatchI,
367  nullFaceSets,
368  nPatchPrims.find(patchName)(),
370  );
371  }
372 
373  if (Pstream::master())
374  {
375  delete ensightFilePtr;
376  }
377 }
378 
379 
380 template<class Type>
381 void ensightField
382 (
384  const ensightMesh& eMesh,
385  const fileName& postProcPath,
386  const word& prepend,
387  const label timeIndex,
388  const bool binary,
389  Ostream& ensightCaseFile
390 )
391 {
392  Info<< "Converting field " << vf.name() << endl;
393 
394  word timeFile = prepend + itoa(timeIndex);
395 
396  const fvMesh& mesh = eMesh.mesh();
397  const Time& runTime = mesh.time();
398 
399  const cellSets& meshCellSets = eMesh.meshCellSets();
400  const List<faceSets>& boundaryFaceSets = eMesh.boundaryFaceSets();
401  const wordList& allPatchNames = eMesh.allPatchNames();
402  const wordHashSet& patchNames = eMesh.patchNames();
404  nPatchPrims = eMesh.nPatchPrims();
405  const List<faceSets>& faceZoneFaceSets = eMesh.faceZoneFaceSets();
406  const wordHashSet& faceZoneNames = eMesh.faceZoneNames();
408  nFaceZonePrims = eMesh.nFaceZonePrims();
409 
410  const labelList& tets = meshCellSets.tets;
411  const labelList& pyrs = meshCellSets.pyrs;
412  const labelList& prisms = meshCellSets.prisms;
413  const labelList& wedges = meshCellSets.wedges;
414  const labelList& hexes = meshCellSets.hexes;
415  const labelList& polys = meshCellSets.polys;
416 
417  ensightStream* ensightFilePtr = NULL;
418  if (Pstream::master())
419  {
420  // set the filename of the ensight file
421  fileName ensightFileName(timeFile + "." + vf.name());
422 
423  if (binary)
424  {
425  ensightFilePtr = new ensightBinaryStream
426  (
427  postProcPath/ensightFileName,
428  runTime
429  );
430  }
431  else
432  {
433  ensightFilePtr = new ensightAsciiStream
434  (
435  postProcPath/ensightFileName,
436  runTime
437  );
438  }
439  }
440 
441  ensightStream& ensightFile = *ensightFilePtr;
442 
443  if (Pstream::master())
444  {
445  if (timeIndex == 0)
446  {
447  ensightCaseFile.setf(ios_base::left);
448 
449  ensightCaseFile
451  << " per element: 1 "
452  << setw(15) << vf.name()
453  << (' ' + prepend + "****." + vf.name()).c_str()
454  << nl;
455  }
456 
458  }
459 
460  if (patchNames.empty())
461  {
462  eMesh.barrier();
463 
464  if (Pstream::master())
465  {
466  ensightFile.writePartHeader(1);
467  }
468 
469  writeField
470  (
471  "hexa8",
472  map(vf, hexes, wedges),
474  );
475 
476  writeField
477  (
478  "penta6",
479  Field<Type>(vf, prisms),
481  );
482 
483  writeField
484  (
485  "pyramid5",
486  Field<Type>(vf, pyrs),
488  );
489 
490  writeField
491  (
492  "tetra4",
493  Field<Type>(vf, tets),
495  );
496 
497  writeField
498  (
499  "nfaced",
500  Field<Type>(vf, polys),
502  );
503  }
504 
505  label ensightPatchI = eMesh.patchPartOffset();
506 
507  forAll(allPatchNames, patchi)
508  {
509  const word& patchName = allPatchNames[patchi];
510 
511  eMesh.barrier();
512 
513  if (patchNames.empty() || patchNames.found(patchName))
514  {
515  if
516  (
518  (
519  vf.boundaryField()[patchi],
520  patchi,
521  ensightPatchI,
522  boundaryFaceSets[patchi],
523  nPatchPrims.find(patchName)(),
525  )
526  )
527  {
528  ensightPatchI++;
529  }
530  }
531  }
532 
533  // write faceZones, if requested
534  if (faceZoneNames.size())
535  {
536  // Interpolates cell values to faces - needed only when exporting
537  // faceZones...
539  (
541  );
542 
543  forAllConstIter(wordHashSet, faceZoneNames, iter)
544  {
545  const word& faceZoneName = iter.key();
546 
547  eMesh.barrier();
548 
549  label zoneID = mesh.faceZones().findZoneID(faceZoneName);
550 
551  const faceZone& fz = mesh.faceZones()[zoneID];
552 
553  // Prepare data to write
554  label nIncluded = 0;
555  forAll(fz, i)
556  {
557  if (eMesh.faceToBeIncluded(fz[i]))
558  {
559  ++nIncluded;
560  }
561  }
562 
563  Field<Type> values(nIncluded);
564 
565  // Loop on the faceZone and store the needed field values
566  label j = 0;
567  forAll(fz, i)
568  {
569  label faceI = fz[i];
570  if (mesh.isInternalFace(faceI))
571  {
572  values[j] = sf[faceI];
573  ++j;
574  }
575  else
576  {
577  if (eMesh.faceToBeIncluded(faceI))
578  {
579  label patchI = mesh.boundaryMesh().whichPatch(faceI);
580  const polyPatch& pp = mesh.boundaryMesh()[patchI];
581  label patchFaceI = pp.whichFace(faceI);
582  Type value = sf.boundaryField()[patchI][patchFaceI];
583  values[j] = value;
584  ++j;
585  }
586  }
587  }
588 
589  if
590  (
592  (
593  values,
594  zoneID,
595  ensightPatchI,
596  faceZoneFaceSets[zoneID],
597  nFaceZonePrims.find(faceZoneName)(),
599  )
600  )
601  {
602  ensightPatchI++;
603  }
604  }
605  }
606  if (Pstream::master())
607  {
608  delete ensightFilePtr;
609  }
610 }
611 
612 
613 template<class Type>
614 void ensightPointField
615 (
617  const ensightMesh& eMesh,
618  const fileName& postProcPath,
619  const word& prepend,
620  const label timeIndex,
621  const bool binary,
622  Ostream& ensightCaseFile
623 )
624 {
625  Info<< "Converting field " << pf.name() << endl;
626 
627  word timeFile = prepend + itoa(timeIndex);
628 
629  const fvMesh& mesh = eMesh.mesh();
630  const wordList& allPatchNames = eMesh.allPatchNames();
631  const wordHashSet& patchNames = eMesh.patchNames();
632  const wordHashSet& faceZoneNames = eMesh.faceZoneNames();
633 
634 
635  ensightStream* ensightFilePtr = NULL;
636  if (Pstream::master())
637  {
638  // set the filename of the ensight file
639  fileName ensightFileName(timeFile + "." + pf.name());
640 
641  if (binary)
642  {
643  ensightFilePtr = new ensightBinaryStream
644  (
645  postProcPath/ensightFileName,
646  mesh.time()
647  );
648  }
649  else
650  {
651  ensightFilePtr = new ensightAsciiStream
652  (
653  postProcPath/ensightFileName,
654  mesh.time()
655  );
656  }
657  }
658 
659  ensightStream& ensightFile = *ensightFilePtr;
660 
661  if (Pstream::master())
662  {
663  if (timeIndex == 0)
664  {
665  ensightCaseFile.setf(ios_base::left);
666 
667  ensightCaseFile
669  << " per node: 1 "
670  << setw(15) << pf.name()
671  << (' ' + prepend + "****." + pf.name()).c_str()
672  << nl;
673  }
674 
676  }
677 
678  if (eMesh.patchNames().empty())
679  {
680  eMesh.barrier();
681 
682  if (Pstream::master())
683  {
684  ensightFile.writePartHeader(1);
685  }
686 
687  writeField
688  (
689  "coordinates",
692  );
693  }
694 
695 
696  label ensightPatchI = eMesh.patchPartOffset();
697 
698  forAll(allPatchNames, patchi)
699  {
700  const word& patchName = allPatchNames[patchi];
701 
702  eMesh.barrier();
703 
704  if (patchNames.empty() || patchNames.found(patchName))
705  {
706  const fvPatch& p = mesh.boundary()[patchi];
707  if
708  (
709  returnReduce(p.size(), sumOp<label>())
710  > 0
711  )
712  {
713  // Renumber the patch points/faces into unique points
714  labelList pointToGlobal;
715  labelList uniqueMeshPointLabels;
716  autoPtr<globalIndex> globalPointsPtr =
718  (
719  p.patch().meshPoints(),
720  p.patch().meshPointMap(),
721  pointToGlobal,
722  uniqueMeshPointLabels
723  );
724 
725  if (Pstream::master())
726  {
727  ensightFile.writePartHeader(ensightPatchI);
728  }
729 
730  writeField
731  (
732  "coordinates",
733  Field<Type>(pf.internalField(), uniqueMeshPointLabels),
735  );
736 
737  ensightPatchI++;
738  }
739  }
740  }
741 
742  // write faceZones, if requested
743  if (faceZoneNames.size())
744  {
745  forAllConstIter(wordHashSet, faceZoneNames, iter)
746  {
747  const word& faceZoneName = iter.key();
748 
749  eMesh.barrier();
750 
751  label zoneID = mesh.faceZones().findZoneID(faceZoneName);
752 
753  const faceZone& fz = mesh.faceZones()[zoneID];
754 
755  if (returnReduce(fz().nPoints(), sumOp<label>()) > 0)
756  {
757  // Renumber the faceZone points/faces into unique points
758  labelList pointToGlobal;
759  labelList uniqueMeshPointLabels;
760  autoPtr<globalIndex> globalPointsPtr =
762  (
763  fz().meshPoints(),
764  fz().meshPointMap(),
765  pointToGlobal,
766  uniqueMeshPointLabels
767  );
768 
769  if (Pstream::master())
770  {
771  ensightFile.writePartHeader(ensightPatchI);
772  }
773 
774  writeField
775  (
776  "coordinates",
778  (
779  pf.internalField(),
780  uniqueMeshPointLabels
781  ),
783  );
784 
785  ensightPatchI++;
786  }
787  }
788  }
789 
790  if (Pstream::master())
791  {
792  delete ensightFilePtr;
793  }
794 }
795 
796 
797 template<class Type>
798 void ensightField
799 (
801  const ensightMesh& eMesh,
802  const fileName& postProcPath,
803  const word& prepend,
804  const label timeIndex,
805  const bool binary,
806  const bool nodeValues,
807  Ostream& ensightCaseFile
808 )
809 {
810  if (nodeValues)
811  {
813  (
815  );
816  pfld().rename(vf.name());
817 
818  ensightPointField<Type>
819  (
820  pfld,
821  eMesh,
822  postProcPath,
823  prepend,
824  timeIndex,
825  binary,
826  ensightCaseFile
827  );
828  }
829  else
830  {
831  ensightField<Type>
832  (
833  vf,
834  eMesh,
835  postProcPath,
836  prepend,
837  timeIndex,
838  binary,
839  ensightCaseFile
840  );
841  }
842 }
843 
844 
845 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
volFields.H
Foam::fvMeshSubset::interpolate
static tmp< GeometricField< Type, fvPatchField, volMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &, const fvMesh &sMesh, const labelList &patchMap, const labelList &cellMap, const labelList &faceMap)
Map volume field.
Definition: fvMeshSubsetInterpolate.C:43
Foam::ensightStream
Abstract base class for writing Ensight data.
Definition: ensightStream.H:50
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::ensightMesh
Definition: ensightMesh.H:62
Foam::ensightAsciiStream
Definition: ensightAsciiStream.H:48
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
Foam::ensightMesh::faceZoneNames
const wordHashSet & faceZoneNames() const
Definition: ensightMesh.H:316
p
p
Definition: pEqn.H:62
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::fileName
A class for handling file names.
Definition: fileName.H:69
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:86
Foam::UPstream::scheduled
@ scheduled
Definition: UPstream.H:67
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::fvMeshSubset
Post-processing mesh subset tool. Given the original mesh and the list of selected cells,...
Definition: fvMeshSubset.H:73
Foam::OPstream
Output inter-processor communications stream.
Definition: OPstream.H:50
Foam::volMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: volMesh.H:47
globalIndex.H
Foam::UPstream::nProcs
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:387
Foam::GeometricField::boundaryField
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
Definition: GeometricField.C:735
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::ensightMesh::faceZoneFaceSets
const List< faceSets > & faceZoneFaceSets() const
Definition: ensightMesh.H:311
Foam::ensightMesh::nFacePrimitives::nQuads
label nQuads
Definition: ensightMesh.H:70
Foam::ensightMesh::mesh
const fvMesh & mesh() const
Definition: ensightMesh.H:281
ensightBinaryStream.H
zeroGradientFvPatchField.H
Foam::ensightMesh::faceToBeIncluded
bool faceToBeIncluded(const label faceI) const
When exporting faceZones, check if a given face has to be included.
Foam::HashSet
A HashTable with keys but without contents.
Definition: HashSet.H:59
Foam::MeshObject< fvMesh, UpdateableMeshObject, volPointInterpolation >::New
static const volPointInterpolation & New(const fvMesh &mesh)
Definition: MeshObject.C:44
Foam::zeroGradientFvPatchField
This boundary condition applies a zero-gradient condition from the patch internal field onto the patc...
Definition: zeroGradientFvPatchField.H:63
Foam::ensightMesh::nFaceZonePrims
const HashTable< nFacePrimitives > & nFaceZonePrims() const
Definition: ensightMesh.H:321
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
OFstream.H
Foam::polyMesh::faceZones
const faceZoneMesh & faceZones() const
Return face zone mesh.
Definition: polyMesh.H:463
Foam::faceSets
Definition: faceSets.H:44
Foam::globalMeshData::mergePoints
autoPtr< globalIndex > mergePoints(labelList &pointToGlobal, labelList &uniquePoints) const
Helper for merging (collocated!) mesh point data.
Definition: globalMeshData.C:2374
nPoints
label nPoints
Definition: gmvOutputHeader.H:2
Foam::ensightPTraits
Conversion of OpenFOAM pTraits into the Ensight equivalent.
Definition: ensightPTraits.H:48
Foam::ensightMesh::patchNames
const wordHashSet & patchNames() const
Definition: ensightMesh.H:301
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::Field< Type >
Foam::faceZone
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:64
Foam::HashTable::empty
bool empty() const
Return true if the hash table is empty.
Definition: HashTableI.H:72
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::ensightMesh::nFacePrimitives
Definition: ensightMesh.H:65
Foam::ensightMesh::nFacePrimitives::nPolys
label nPolys
Definition: ensightMesh.H:71
Foam::Info
messageStream Info
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
ensightPTraits.H
Foam::faceSets::quads
labelList quads
Definition: faceSets.H:53
writePatchField
void writePatchField(const Foam::word &fieldName, const Foam::Field< Type > &pf, const Foam::word &patchName, const Foam::ensightMesh &eMesh, const Foam::fileName &postProcPath, const Foam::word &prepend, const Foam::label timeIndex, Foam::Ostream &ensightCaseFile)
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
IOmanip.H
Istream and Ostream manipulators taking arguments.
Foam::ensightFile::write
virtual Ostream & write(const char *buf, std::streamsize count)
Binary write.
Definition: ensightFile.C:136
Foam::GeometricField::internalField
InternalField & internalField()
Return internal field.
Definition: GeometricField.C:724
Foam::polyBoundaryMesh::whichPatch
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
Definition: polyBoundaryMesh.C:703
ensightAsciiStream.H
patchNames
wordList patchNames(nPatches)
Foam::ensightMesh::boundaryFaceSets
const List< faceSets > & boundaryFaceSets() const
Definition: ensightMesh.H:291
Foam::UPstream::masterNo
static int masterNo()
Process index of the master.
Definition: UPstream.H:393
Foam::cellSets
Definition: cellSets.H:44
Foam::cellSets::pyrs
labelList pyrs
Definition: cellSets.H:55
timeIndex
label timeIndex
Definition: getTimeIndex.H:4
Foam::HashTable::size
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
Foam::ensightMesh::patchPartOffset
label patchPartOffset() const
The ensight part id for the first patch.
Definition: ensightMesh.H:327
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::ensightBinaryStream
Definition: ensightBinaryStream.H:47
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
fvMesh.H
Foam::ensightFile
Ensight output with specialized write() for strings, integers and floats. Correctly handles binary wr...
Definition: ensightFile.H:47
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
sf
volScalarField sf(fieldObject, mesh)
Foam::setw
Omanip< int > setw(const int i)
Definition: IOmanip.H:199
Foam::ZoneMesh::findZoneID
label findZoneID(const word &zoneName) const
Find zone index given a name.
Definition: ZoneMesh.C:348
Foam::Field::component
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
Definition: Field.C:644
Foam::cellSets::tets
labelList tets
Definition: cellSets.H:54
Foam::HashTable::find
iterator find(const Key &)
Find and return an iterator set at the hashedEntry.
Foam::volPointInterpolation::interpolate
tmp< GeometricField< Type, pointPatchField, pointMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &) const
Interpolate volField using inverse distance weighting.
Definition: volPointInterpolate.C:502
Foam::IOstream::setf
ios_base::fmtflags setf(const ios_base::fmtflags f)
Set flags of stream.
Definition: IOstream.H:496
Foam::itoa
word itoa(const label n)
Foam::UPstream::master
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:399
Foam::HashTable
An STL-conforming hash table.
Definition: HashTable.H:61
volPointInterpolation.H
Foam::fvMesh::boundary
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:550
Foam::cellSets::wedges
labelList wedges
Definition: cellSets.H:57
Foam::autoPtr
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
Foam::cellSets::polys
labelList polys
Definition: cellSets.H:59
Foam::cellSets::hexes
labelList hexes
Definition: cellSets.H:58
Foam::polyPatch::whichFace
label whichFace(const label l) const
Return label of face in patch from global face label.
Definition: polyPatch.H:389
Foam::ensightMesh::meshCellSets
const cellSets & meshCellSets() const
Definition: ensightMesh.H:286
Foam::linearInterpolate
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > linearInterpolate(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: linear.H:107
Foam::sumOp
Definition: ops.H:162
ensightField
void ensightField(const Foam::GeometricField< Type, Foam::fvPatchField, Foam::volMesh > &vf, const Foam::ensightMesh &eMesh, const Foam::fileName &postProcPath, const Foam::word &prepend, const Foam::label timeIndex, const bool binary, const bool nodeValues, Foam::Ostream &ensightCaseFile)
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::ensightMesh::allPatchNames
const wordList & allPatchNames() const
Definition: ensightMesh.H:296
Foam::primitiveMesh::isInternalFace
bool isInternalFace(const label faceIndex) const
Return true if given face label is internal to the mesh.
Definition: primitiveMeshI.H:70
patchi
label patchi
Definition: getPatchFieldScalar.H:1
Foam::direction
unsigned char direction
Definition: direction.H:43
Foam::ensightMesh::nFacePrimitives::nTris
label nTris
Definition: ensightMesh.H:69
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
ensightField.H
Foam::IPstream
Input inter-processor communications stream.
Definition: IPstream.H:50
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::polyMesh::globalData
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1143
Foam::ensightMesh::nPatchPrims
const HashTable< nFacePrimitives > & nPatchPrims() const
Definition: ensightMesh.H:306
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
volField
Foam::tmp< Foam::GeometricField< Type, Foam::fvPatchField, Foam::volMesh > > volField(const Foam::fvMeshSubset &, const Foam::GeometricField< Type, Foam::fvPatchField, Foam::volMesh > &vf)
Wrapper to get hold of the field or the subsetted field.
Foam::faceSets::tris
labelList tris
Definition: faceSets.H:52
Foam::faceSets::polys
labelList polys
Definition: faceSets.H:54
Foam::fvMeshSubset::hasSubMesh
bool hasSubMesh() const
Have subMesh?
Definition: fvMeshSubset.C:1466
Foam::ensightMesh::barrier
static void barrier()
Helper to cause barrier. Necessary on Quadrics.
Foam::ensightMesh::uniquePointMap
const labelList & uniquePointMap() const
Local points that are unique.
Definition: ensightMesh.H:348
Foam::cellSets::prisms
labelList prisms
Definition: cellSets.H:56