Merge branch 'master' of ssh://opencfd@repo.or.cz/srv/git/OpenFOAM-1.5.x
[OpenFOAM-1.5.x.git] / src / sampling / cuttingPlane / cuttingPlane.H
blob3e33be9da769c887a842b0a7c7b61084cc0bb671
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2008 OpenCFD Ltd.
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 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
19     for more details.
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 Class
26     Foam::cuttingPlane
28 Description
29     Constructs plane through mesh.
31     No attempt at resolving degenerate cases.
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"
49 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
51 namespace Foam
54 class primitiveMesh;
56 /*---------------------------------------------------------------------------*\
57                            Class cuttingPlane Declaration
58 \*---------------------------------------------------------------------------*/
60 class cuttingPlane
62     public plane
64     // Private data
66         //- List of cells cut by the plane
67         labelList cutCells_;
69         //- Intersection points
70         pointField cutPoints_;
72         //- Cut faces in terms of intersection points
73         faceList cutFaces_;
76     // Private Member Functions
78         //- Determine cut cells, possibly restricted to a list of cells
79         void calcCutCells
80         (
81             const primitiveMesh&,
82             const scalarField& dotProducts,
83             const labelList& cellIdLabels = labelList::null()
84         );
86         //- Determine intersection points (cutPoints_).
87         labelList intersectEdges
88         (
89             const primitiveMesh&,
90             const scalarField& dotProducts
91         );
93         //- Walk around circumference of cell starting from startEdgeI crossing
94         //  only cut edges. Record cutPoint labels in faceVerts.
95         static bool walkCell
96         (
97             const primitiveMesh&,
98             const labelList& edgePoint,
99             const label cellI,
100             const label startEdgeI,
101             DynamicList<label>& faceVerts
102         );
104         //- Determine cuts for all cut cells.
105         void walkCellCuts
106         (
107             const primitiveMesh& mesh,
108             const labelList& edgePoint
109         );
112 protected:
114     // Constructors
116         //- Construct plane description without cutting
117         cuttingPlane(const plane&);
120     // Protected Member Functions
122         //- recut mesh with existing planeDesc
123         void reCut(const primitiveMesh&);
125         //- recut mesh with existing planeDesc, restricted to a list of cells
126         void reCut
127         (
128             const primitiveMesh&,
129             const labelList& cellIdLabels
130         );
133 public:
135     // Constructors
137         //- Construct from components: Mesh reference and plane
138         cuttingPlane(const primitiveMesh&, const plane&);
140         //- Construct from mesh reference and plane,
141         //  restricted to a list of cells
142         cuttingPlane
143         (
144             const primitiveMesh&,
145             const plane&,
146             const labelList& cellIdLabels
147         );
150     // Member Functions
152         //- Return plane used
153         const plane& planeDesc() const;
155         //- Return pointField of cutting points
156         const pointField& points() const;
158         //- Return faceList of points in cells
159         const faceList& faces() const;
161         //- Return labelList of cut cells
162         const labelList& cells() const;
164         //- Return true or false to question: have any cells been cut?
165         bool cut();
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 "cuttingPlaneSample.C"
189 #endif
191 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
193 #endif
195 // ************************************************************************* //