FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

dynamicRefineFvMesh.H

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------*\
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 1991-2010 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software: you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
00023 
00024 Class
00025     Foam::dynamicRefineFvMesh
00026 
00027 Description
00028     A fvMesh with built-in refinement.
00029 
00030     Determines which cells to refine/unrefine and does all in update().
00031 
00032 SourceFiles
00033     dynamicRefineFvMesh.C
00034 
00035 \*---------------------------------------------------------------------------*/
00036 
00037 #ifndef dynamicRefineFvMesh_H
00038 #define dynamicRefineFvMesh_H
00039 
00040 #include <dynamicFvMesh/dynamicFvMesh.H>
00041 #include <dynamicMesh/hexRef8.H>
00042 #include <OpenFOAM/PackedBoolList.H>
00043 #include <OpenFOAM/Switch.H>
00044 
00045 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00046 
00047 namespace Foam
00048 {
00049 
00050 /*---------------------------------------------------------------------------*\
00051                            Class dynamicRefineFvMesh Declaration
00052 \*---------------------------------------------------------------------------*/
00053 
00054 class dynamicRefineFvMesh
00055 :
00056     public dynamicFvMesh
00057 {
00058 protected:
00059 
00060         //- Mesh cutting engine
00061         hexRef8 meshCutter_;
00062 
00063         //- Dump cellLevel for postprocessing
00064         Switch dumpLevel_;
00065 
00066         //- Fluxes to map
00067         List<Pair<word> > correctFluxes_;
00068 
00069         //- Number of refinement/unrefinement steps done so far.
00070         label nRefinementIterations_;
00071 
00072         //- Protected cells (usually since not hexes)
00073         PackedBoolList protectedCell_;
00074 
00075 
00076     // Private Member Functions
00077 
00078         //- Count set/unset elements in packedlist.
00079         static label count(const PackedBoolList&, const unsigned int);
00080 
00081         //- Calculate cells that cannot be refined since would trigger
00082         //  refinement of protectedCell_ (since 2:1 refinement cascade)
00083         void calculateProtectedCells(PackedBoolList& unrefineableCell) const;
00084 
00085         //- Read the projection parameters from dictionary
00086         void readDict();
00087 
00088 
00089         //- Refine cells. Update mesh and fields.
00090         autoPtr<mapPolyMesh> refine(const labelList&);
00091 
00092         //- Unrefine cells. Gets passed in centre points of cells to combine.
00093         autoPtr<mapPolyMesh> unrefine(const labelList&);
00094 
00095 
00096         // Selection of cells to un/refine
00097 
00098             //- Calculates approximate value for refinement level so
00099             //  we don't go above maxCell
00100             scalar getRefineLevel
00101             (
00102                 const label maxCells,
00103                 const label maxRefinement,
00104                 const scalar refineLevel,
00105                 const scalarField&
00106             ) const;
00107 
00108             //- Get per cell max of connected point
00109             scalarField maxPointField(const scalarField&) const;
00110 
00111             //- Get point min of connected cell
00112             scalarField minCellField(const volScalarField&) const;
00113 
00114             scalarField cellToPoint(const scalarField& vFld) const;
00115 
00116             scalarField error
00117             (
00118                 const scalarField& fld,
00119                 const scalar minLevel,
00120                 const scalar maxLevel
00121             ) const;
00122 
00123             //- Select candidate cells for refinement
00124             virtual void selectRefineCandidates
00125             (
00126                 const scalar lowerRefineLevel,
00127                 const scalar upperRefineLevel,
00128                 const scalarField& vFld,
00129                 PackedBoolList& candidateCell
00130             ) const;
00131 
00132             //- Subset candidate cells for refinement
00133             virtual labelList selectRefineCells
00134             (
00135                 const label maxCells,
00136                 const label maxRefinement,
00137                 const PackedBoolList& candidateCell
00138             ) const;
00139 
00140             //- Select points that can be unrefined.
00141             virtual labelList selectUnrefinePoints
00142             (
00143                 const scalar unrefineLevel,
00144                 const PackedBoolList& markedCell,
00145                 const scalarField& pFld
00146             ) const;
00147 
00148             //- Extend markedCell with cell-face-cell.
00149             void extendMarkedCells(PackedBoolList& markedCell) const;
00150 
00151 
00152 private:
00153 
00154         //- Disallow default bitwise copy construct
00155         dynamicRefineFvMesh(const dynamicRefineFvMesh&);
00156 
00157         //- Disallow default bitwise assignment
00158         void operator=(const dynamicRefineFvMesh&);
00159 
00160 public:
00161 
00162     //- Runtime type information
00163     TypeName("dynamicRefineFvMesh");
00164 
00165 
00166     // Constructors
00167 
00168         //- Construct from IOobject
00169         explicit dynamicRefineFvMesh(const IOobject& io);
00170 
00171 
00172     // Destructor
00173 
00174         virtual ~dynamicRefineFvMesh();
00175 
00176 
00177     // Member Functions
00178 
00179         //- Direct access to the refinement engine
00180         const hexRef8& meshCutter() const
00181         {
00182             return meshCutter_;
00183         }
00184 
00185         //- Cells which should not be refined/unrefined
00186         const PackedBoolList& protectedCell() const
00187         {
00188             return protectedCell_;
00189         }
00190 
00191         //- Cells which should not be refined/unrefined
00192         PackedBoolList& protectedCell()
00193         {
00194             return protectedCell_;
00195         }
00196 
00197         //- Update the mesh for both mesh motion and topology change
00198         virtual bool update();
00199 
00200 
00201     // Writing
00202 
00203         //- Write using given format, version and compression
00204         virtual bool writeObject
00205         (
00206             IOstream::streamFormat fmt,
00207             IOstream::versionNumber ver,
00208             IOstream::compressionType cmp
00209         ) const;
00210 
00211 };
00212 
00213 
00214 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00215 
00216 } // End namespace Foam
00217 
00218 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00219 
00220 #endif
00221 
00222 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines