/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | \\ / A nd | Copyright held by original author \\/ M anipulation | ------------------------------------------------------------------------------- License This file is part of OpenFOAM. OpenFOAM is free software; you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation; either version 2 of the License, or (at your option) any later version. OpenFOAM is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details. You should have received a copy of the GNU General Public License along with OpenFOAM; if not, write to the Free Software Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA Class Foam::vtkPV3Foam Description Provides a reader interface for OpenFOAM to VTK interaction. SourceFiles vtkPV3Foam.C vtkPV3Foam.H vtkPV3FoamI.H vtkPV3FoamFields.C vtkPV3FoamMesh.C vtkPV3FoamMeshLagrangian.C vtkPV3FoamMeshPatch.C vtkPV3FoamMeshSet.C vtkPV3FoamMeshVolume.C vtkPV3FoamMeshZone.C vtkPV3FoamFaceField.H vtkPV3FoamLagrangianFields.H vtkPV3FoamPatchField.H vtkPV3FoamPointFields.H vtkPV3FoamPoints.H vtkPV3FoamUpdateInfo.C vtkPV3FoamUpdateInfoFields.H vtkPV3FoamUtilities.C vtkPV3FoamVolFields.H vtkPV3FoamAddToSelection.H // Needed by VTK: vtkDataArrayTemplateImplicit.txx \*---------------------------------------------------------------------------*/ #ifndef vtkPV3Foam_H #define vtkPV3Foam_H // do not include legacy strstream headers #ifndef VTK_EXCLUDE_STRSTREAM_HEADERS # define VTK_EXCLUDE_STRSTREAM_HEADERS #endif #include "className.H" #include "fileName.H" #include "stringList.H" #include "wordList.H" #include "primitivePatch.H" #include "PrimitivePatchInterpolation.H" #include "volPointInterpolation.H" #undef VTKPV3FOAM_DUALPORT // * * * * * * * * * * * * * Forward Declarations * * * * * * * * * * * * * // class vtkDataArraySelection; class vtkDataSet; class vtkPoints; class vtkPV3FoamReader; class vtkRenderer; class vtkTextActor; class vtkMultiBlockDataSet; class vtkPolyData; class vtkUnstructuredGrid; class vtkIndent; // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // namespace Foam { // Foam class forward declarations class argList; class Time; class fvMesh; class IOobjectList; class polyPatch; class faceSet; class pointSet; template class IOField; template class List; /*---------------------------------------------------------------------------*\ Class vtkPV3Foam Declaration \*---------------------------------------------------------------------------*/ class vtkPV3Foam { // Private classes //- Bookkeeping for GUI checklists and the multi-block organization class partInfo { const char *name_; int block_; int start_; int size_; public: partInfo(const char *name, const int blockNo=0) : name_(name), block_(blockNo), start_(-1), size_(0) {} //- Return the block holding these datasets int block() const { return block_; } //- Assign block number, return previous value int block(int blockNo) { int prev = block_; block_ = blockNo; return prev; } const char* name() const { return name_; } int start() const { return start_; } int end() const { return start_ + size_; } int size() const { return size_; } bool empty() const { return !size_; } void reset() { start_ = -1; size_ = 0; } //- Assign new start and reset the size void operator=(const int i) { start_ = i; size_ = 0; } //- Increment the size void operator+=(const int n) { size_ += n; } }; //- bookkeeping for polyhedral cell decomposition // hide in extra pointMap (cellSet/cellZone) for now class polyDecomp { labelList superCells_; labelList addPointCellLabels_; labelList pointMap_; public: polyDecomp() {} //- Label of original cell for decomposed cells labelList& superCells() { return superCells_; } //- Label of original cell for decomposed cells const labelList& superCells() const { return superCells_; } //- Cell-centre labels for additional points of decomposed cells labelList& addPointCellLabels() { return addPointCellLabels_; } //- Cell-centre labels for additional points of decomposed cells const labelList& addPointCellLabels() const { return addPointCellLabels_; } //- Point labels for subsetted meshes labelList& pointMap() { return pointMap_; } //- Point labels for subsetted meshes const labelList& pointMap() const { return pointMap_; } //- Clear void clear() { superCells_.clear(); addPointCellLabels_.clear(); pointMap_.clear(); } }; // Private Data //- Access to the controlling vtkPV3FoamReader vtkPV3FoamReader* reader_; //- Foam time control autoPtr