ThirdParty: correction of spec files for installation of qt and scotch
[OpenFOAM-1.6-ext.git] / src / sampling / cuttingPlane / cuttingPlane.H
blob87774ab1526858ac768cf05b5545de0e082af414
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright held by original author
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. Since the cut faces are
32     usually quite ugly, they will always be triangulated.
34 Note
35     When the cutting plane coincides with a mesh face, the cell edge on the
36     positive side of the plane is taken.
38 SourceFiles
39     cuttingPlane.C
41 \*---------------------------------------------------------------------------*/
43 #ifndef cuttingPlane_H
44 #define cuttingPlane_H
46 #include "plane.H"
47 #include "pointField.H"
48 #include "faceList.H"
49 #include "MeshedSurface.H"
51 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53 namespace Foam
56 class primitiveMesh;
58 /*---------------------------------------------------------------------------*\
59                        Class cuttingPlane Declaration
60 \*---------------------------------------------------------------------------*/
62 class cuttingPlane
64     public plane,
65     public MeshedSurface<face>
67     //- Private typedefs for convenience
68         typedef MeshedSurface<face> MeshStorage;
70     // Private data
72         //- List of cells cut by the plane
73         labelList cutCells_;
75     // Private Member Functions
77         //- Determine cut cells, possibly restricted to a list of cells
78         void calcCutCells
79         (
80             const primitiveMesh&,
81             const scalarField& dotProducts,
82             const UList<label>& cellIdLabels = UList<label>::null()
83         );
85         //- Determine intersection points (cutPoints).
86         void intersectEdges
87         (
88             const primitiveMesh&,
89             const scalarField& dotProducts,
90             List<label>& edgePoint
91         );
93         //- Walk 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 UList<label>& 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 UList<label>& 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 UList<label>& cellIdLabels = UList<label>::null()
126         );
128         //- remap action on triangulation or cleanup
129         virtual void remapFaces(const UList<label>& faceMap);
131 public:
133     // Constructors
135         //- Construct from plane and mesh reference,
136         //  possibly restricted to a list of cells
137         cuttingPlane
138         (
139             const plane&,
140             const primitiveMesh&,
141             const UList<label>& cellIdLabels = UList<label>::null()
142         );
145     // Member Functions
147         //- Return plane used
148         const plane& planeDesc() const
149         {
150             return static_cast<const plane&>(*this);
151         }
153         //- Return List of cells cut by the plane
154         const labelList& cutCells() const
155         {
156             return cutCells_;
157         }
159         //- Return true or false to question: have any cells been cut?
160         bool cut() const
161         {
162             return cutCells_.size();
163         }
165         //- Sample the cell field
166         template<class Type>
167         tmp<Field<Type> > sample(const Field<Type>&) const;
169         template<class Type>
170         tmp<Field<Type> > sample(const tmp<Field<Type> >&) const;
173     // Member Operators
175         void operator=(const cuttingPlane&);
179 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
181 } // End namespace Foam
183 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
185 #ifdef NoRepository
186 #   include "cuttingPlaneTemplates.C"
187 #endif
189 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
191 #endif
193 // ************************************************************************* //