initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / dynamicFvMesh / dynamicRefineFvMesh / dynamicRefineFvMesh.H
blobf8bf0ec3e9b8245182e4922c91723c67a65ba301
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 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::dynamicRefineFvMesh
28 Description
29     A fvMesh with built-in refinement.
31     Determines which cells to refine/unrefine and does all in update().
33 SourceFiles
34     dynamicRefineFvMesh.C
36 \*---------------------------------------------------------------------------*/
38 #ifndef dynamicRefineFvMesh_H
39 #define dynamicRefineFvMesh_H
41 #include "dynamicFvMesh.H"
42 #include "hexRef8.H"
43 #include "PackedBoolList.H"
44 #include "Switch.H"
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 namespace Foam
51 /*---------------------------------------------------------------------------*\
52                            Class dynamicRefineFvMesh Declaration
53 \*---------------------------------------------------------------------------*/
55 class dynamicRefineFvMesh
57     public dynamicFvMesh
59 protected:
61         //- Mesh cutting engine
62         hexRef8 meshCutter_;
64         //- Dump cellLevel for postprocessing
65         Switch dumpLevel_;
67         //- Fluxes to map
68         List<Pair<word> > correctFluxes_;
70         //- Number of refinement/unrefinement steps done so far.
71         label nRefinementIterations_;
73         //- Protected cells (usually since not hexes)
74         PackedBoolList protectedCell_;
77     // Private Member Functions
79         //- Count set/unset elements in packedlist.
80         static label count(const PackedBoolList&, const unsigned int);
82         //- Calculate cells that cannot be refined since would trigger
83         //  refinement of protectedCell_ (since 2:1 refinement cascade)
84         void calculateProtectedCells(PackedBoolList& unrefineableCell) const;
86         //- Read the projection parameters from dictionary
87         void readDict();
90         //- Refine cells. Update mesh and fields.
91         autoPtr<mapPolyMesh> refine(const labelList&);
93         //- Unrefine cells. Gets passed in centre points of cells to combine.
94         autoPtr<mapPolyMesh> unrefine(const labelList&);
97         // Selection of cells to un/refine
99             //- Calculates approximate value for refinement level so
100             //  we don't go above maxCell
101             scalar getRefineLevel
102             (
103                 const label maxCells,
104                 const label maxRefinement,
105                 const scalar refineLevel,
106                 const scalarField&
107             ) const;
109             //- Get per cell max of connected point
110             scalarField maxPointField(const scalarField&) const;
112             //- Get point min of connected cell
113             scalarField minCellField(const volScalarField&) const;
115             scalarField cellToPoint(const scalarField& vFld) const;
117             scalarField error
118             (
119                 const scalarField& fld,
120                 const scalar minLevel,
121                 const scalar maxLevel
122             ) const;
124             //- Select candidate cells for refinement
125             virtual void selectRefineCandidates
126             (
127                 const scalar lowerRefineLevel,
128                 const scalar upperRefineLevel,
129                 const scalarField& vFld,
130                 PackedBoolList& candidateCell
131             ) const;
133             //- Subset candidate cells for refinement
134             virtual labelList selectRefineCells
135             (
136                 const label maxCells,
137                 const label maxRefinement,
138                 const PackedBoolList& candidateCell
139             ) const;
141             //- Select points that can be unrefined.
142             virtual labelList selectUnrefinePoints
143             (
144                 const scalar unrefineLevel,
145                 const PackedBoolList& markedCell,
146                 const scalarField& pFld
147             ) const;
149             //- Extend markedCell with cell-face-cell.
150             void extendMarkedCells(PackedBoolList& markedCell) const;
153 private:
155         //- Disallow default bitwise copy construct
156         dynamicRefineFvMesh(const dynamicRefineFvMesh&);
158         //- Disallow default bitwise assignment
159         void operator=(const dynamicRefineFvMesh&);
161 public:
163     //- Runtime type information
164     TypeName("dynamicRefineFvMesh");
167     // Constructors
169         //- Construct from IOobject
170         explicit dynamicRefineFvMesh(const IOobject& io);
173     // Destructor
175         virtual ~dynamicRefineFvMesh();
178     // Member Functions
180         //- Direct access to the refinement engine
181         const hexRef8& meshCutter() const
182         {
183             return meshCutter_;
184         }
186         //- Cells which should not be refined/unrefined
187         const PackedBoolList& protectedCell() const
188         {
189             return protectedCell_;
190         }
192         //- Cells which should not be refined/unrefined
193         PackedBoolList& protectedCell()
194         {
195             return protectedCell_;
196         }
198         //- Update the mesh for both mesh motion and topology change
199         virtual bool update();
202     // Writing
204         //- Write using given format, version and compression
205         virtual bool writeObject
206         (
207             IOstream::streamFormat fmt,
208             IOstream::versionNumber ver,
209             IOstream::compressionType cmp
210         ) const;
215 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
217 } // End namespace Foam
219 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
221 #endif
223 // ************************************************************************* //