ENH: autoLayerDriver: better layering information message
[OpenFOAM-2.0.x.git] / src / dynamicFvMesh / dynamicRefineFvMesh / dynamicRefineFvMesh.H
blob2dfe09d07554706ef6246eca8d3a904f82fc0fdd
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::dynamicRefineFvMesh
27 Description
28     A fvMesh with built-in refinement.
30     Determines which cells to refine/unrefine and does all in update().
32 SourceFiles
33     dynamicRefineFvMesh.C
35 \*---------------------------------------------------------------------------*/
37 #ifndef dynamicRefineFvMesh_H
38 #define dynamicRefineFvMesh_H
40 #include "dynamicFvMesh.H"
41 #include "hexRef8.H"
42 #include "PackedBoolList.H"
43 #include "Switch.H"
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 namespace Foam
50 /*---------------------------------------------------------------------------*\
51                            Class dynamicRefineFvMesh Declaration
52 \*---------------------------------------------------------------------------*/
54 class dynamicRefineFvMesh
56     public dynamicFvMesh
58 protected:
60         //- Mesh cutting engine
61         hexRef8 meshCutter_;
63         //- Dump cellLevel for postprocessing
64         Switch dumpLevel_;
66         //- Fluxes to map
67         HashTable<word> correctFluxes_;
69         //- Number of refinement/unrefinement steps done so far.
70         label nRefinementIterations_;
72         //- Protected cells (usually since not hexes)
73         PackedBoolList protectedCell_;
76     // Private Member Functions
78         //- Count set/unset elements in packedlist.
79         static label count(const PackedBoolList&, const unsigned int);
81         //- Calculate cells that cannot be refined since would trigger
82         //  refinement of protectedCell_ (since 2:1 refinement cascade)
83         void calculateProtectedCells(PackedBoolList& unrefineableCell) const;
85         //- Read the projection parameters from dictionary
86         void readDict();
89         //- Refine cells. Update mesh and fields.
90         autoPtr<mapPolyMesh> refine(const labelList&);
92         //- Unrefine cells. Gets passed in centre points of cells to combine.
93         autoPtr<mapPolyMesh> unrefine(const labelList&);
96         // Selection of cells to un/refine
98             //- Calculates approximate value for refinement level so
99             //  we don't go above maxCell
100             scalar getRefineLevel
101             (
102                 const label maxCells,
103                 const label maxRefinement,
104                 const scalar refineLevel,
105                 const scalarField&
106             ) const;
108             //- Get per cell max of connected point
109             scalarField maxPointField(const scalarField&) const;
111             //- Get point min of connected cell
112             scalarField minCellField(const volScalarField&) const;
114             scalarField cellToPoint(const scalarField& vFld) const;
116             scalarField error
117             (
118                 const scalarField& fld,
119                 const scalar minLevel,
120                 const scalar maxLevel
121             ) const;
123             //- Select candidate cells for refinement
124             virtual void selectRefineCandidates
125             (
126                 const scalar lowerRefineLevel,
127                 const scalar upperRefineLevel,
128                 const scalarField& vFld,
129                 PackedBoolList& candidateCell
130             ) const;
132             //- Subset candidate cells for refinement
133             virtual labelList selectRefineCells
134             (
135                 const label maxCells,
136                 const label maxRefinement,
137                 const PackedBoolList& candidateCell
138             ) const;
140             //- Select points that can be unrefined.
141             virtual labelList selectUnrefinePoints
142             (
143                 const scalar unrefineLevel,
144                 const PackedBoolList& markedCell,
145                 const scalarField& pFld
146             ) const;
148             //- Extend markedCell with cell-face-cell.
149             void extendMarkedCells(PackedBoolList& markedCell) const;
152 private:
154         //- Disallow default bitwise copy construct
155         dynamicRefineFvMesh(const dynamicRefineFvMesh&);
157         //- Disallow default bitwise assignment
158         void operator=(const dynamicRefineFvMesh&);
160 public:
162     //- Runtime type information
163     TypeName("dynamicRefineFvMesh");
166     // Constructors
168         //- Construct from IOobject
169         explicit dynamicRefineFvMesh(const IOobject& io);
172     //- Destructor
173     virtual ~dynamicRefineFvMesh();
176     // Member Functions
178         //- Direct access to the refinement engine
179         const hexRef8& meshCutter() const
180         {
181             return meshCutter_;
182         }
184         //- Cells which should not be refined/unrefined
185         const PackedBoolList& protectedCell() const
186         {
187             return protectedCell_;
188         }
190         //- Cells which should not be refined/unrefined
191         PackedBoolList& protectedCell()
192         {
193             return protectedCell_;
194         }
196         //- Update the mesh for both mesh motion and topology change
197         virtual bool update();
200     // Writing
202         //- Write using given format, version and compression
203         virtual bool writeObject
204         (
205             IOstream::streamFormat fmt,
206             IOstream::versionNumber ver,
207             IOstream::compressionType cmp
208         ) const;
213 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
215 } // End namespace Foam
217 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
219 #endif
221 // ************************************************************************* //