BUG: UListIO: byteSize overflowing on really big faceLists
[OpenFOAM-2.0.x.git] / src / sampling / cuttingPlane / cuttingPlane.H
blobf23a50f504858260ba2924937819eec3c70c19b0
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
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
13     the Free Software Foundation, either version 3 of the License, or
14     (at your 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, see <http://www.gnu.org/licenses/>.
24 Class
25     Foam::cuttingPlane
27 Description
28     Constructs plane through mesh.
30     No attempt at resolving degenerate cases. Since the cut faces are
31     usually quite ugly, they will always be triangulated.
33 Note
34     When the cutting plane coincides with a mesh face, the cell edge on the
35     positive side of the plane is taken.
37 SourceFiles
38     cuttingPlane.C
40 \*---------------------------------------------------------------------------*/
42 #ifndef cuttingPlane_H
43 #define cuttingPlane_H
45 #include "plane.H"
46 #include "pointField.H"
47 #include "faceList.H"
48 #include "MeshedSurface.H"
50 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52 namespace Foam
55 class primitiveMesh;
57 /*---------------------------------------------------------------------------*\
58                        Class cuttingPlane Declaration
59 \*---------------------------------------------------------------------------*/
61 class cuttingPlane
63     public plane,
64     public MeshedSurface<face>
66     //- Private typedefs for convenience
67         typedef MeshedSurface<face> MeshStorage;
69     // Private data
71         //- List of cells cut by the plane
72         labelList cutCells_;
74     // Private Member Functions
76         //- Determine cut cells, possibly restricted to a list of cells
77         void calcCutCells
78         (
79             const primitiveMesh&,
80             const scalarField& dotProducts,
81             const labelUList& cellIdLabels = labelUList::null()
82         );
84         //- Determine intersection points (cutPoints).
85         void intersectEdges
86         (
87             const primitiveMesh&,
88             const scalarField& dotProducts,
89             List<label>& edgePoint
90         );
92         //- Walk circumference of cell, starting from startEdgeI crossing
93         //  only cut edges. Record cutPoint labels in faceVerts.
94         static bool walkCell
95         (
96             const primitiveMesh&,
97             const labelUList& edgePoint,
98             const label cellI,
99             const label startEdgeI,
100             DynamicList<label>& faceVerts
101         );
103         //- Determine cuts for all cut cells.
104         void walkCellCuts
105         (
106             const primitiveMesh& mesh,
107             const bool triangulate,
108             const labelUList& edgePoint
109         );
112 protected:
114     // Constructors
116         //- Construct plane description without cutting
117         cuttingPlane(const plane&);
119     // Protected Member Functions
121         //- recut mesh with existing planeDesc, restricted to a list of cells
122         void reCut
123         (
124             const primitiveMesh&,
125             const bool triangulate,
126             const labelUList& cellIdLabels = labelUList::null()
127         );
129         //- remap action on triangulation or cleanup
130         virtual void remapFaces(const labelUList& faceMap);
132 public:
134     // Constructors
136         //- Construct from plane and mesh reference,
137         //  possibly restricted to a list of cells
138         cuttingPlane
139         (
140             const plane&,
141             const primitiveMesh&,
142             const bool triangulate,
143             const labelUList& cellIdLabels = labelUList::null()
144         );
147     // Member Functions
149         //- Return plane used
150         const plane& planeDesc() const
151         {
152             return static_cast<const plane&>(*this);
153         }
155         //- Return List of cells cut by the plane
156         const labelList& cutCells() const
157         {
158             return cutCells_;
159         }
161         //- Return true or false to question: have any cells been cut?
162         bool cut() const
163         {
164             return cutCells_.size();
165         }
167         //- Sample the cell field
168         template<class Type>
169         tmp<Field<Type> > sample(const Field<Type>&) const;
171         template<class Type>
172         tmp<Field<Type> > sample(const tmp<Field<Type> >&) const;
175     // Member Operators
177         void operator=(const cuttingPlane&);
181 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
183 } // End namespace Foam
185 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
187 #ifdef NoRepository
188 #   include "cuttingPlaneTemplates.C"
189 #endif
191 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
193 #endif
195 // ************************************************************************* //