ENH: autoLayerDriver: better layering information message
[OpenFOAM-2.0.x.git] / src / finiteVolume / cfdTools / general / MRF / MRFZone.H
blobc110eb44008d200d8109ac4a885b2ab58ff1aea6
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::MRFZone
27 Description
28     MRF zone definition based on cell zone and parameters
29     obtained from a control dictionary constructed from the given stream.
31     The rotation of the MRF region is defined by an origin and axis of
32     rotation and an angular speed.
34 SourceFiles
35     MRFZone.C
37 \*---------------------------------------------------------------------------*/
39 #ifndef MRFZone_H
40 #define MRFZone_H
42 #include "dictionary.H"
43 #include "wordList.H"
44 #include "labelList.H"
45 #include "dimensionedScalar.H"
46 #include "dimensionedVector.H"
47 #include "volFieldsFwd.H"
48 #include "surfaceFieldsFwd.H"
49 #include "fvMatricesFwd.H"
50 #include "fvMatrices.H"
52 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
54 namespace Foam
57 // Forward declaration of classes
58 class fvMesh;
60 /*---------------------------------------------------------------------------*\
61                            Class MRFZone Declaration
62 \*---------------------------------------------------------------------------*/
64 class MRFZone
66     // Private data
68         const fvMesh& mesh_;
70         const word name_;
72         const dictionary dict_;
74         label cellZoneID_;
76         const wordList excludedPatchNames_;
77         labelList excludedPatchLabels_;
79         //- Internal faces that are part of MRF
80         labelList internalFaces_;
82         //- Outside faces (per patch) that move with the MRF
83         labelListList includedFaces_;
85         //- Excluded faces (per patch) that do not move with the MRF
86         labelListList excludedFaces_;
88         const dimensionedVector origin_;
89         dimensionedVector axis_;
90         const dimensionedScalar omega_;
91         dimensionedVector Omega_;
94     // Private Member Functions
96         //- Divide faces in frame according to patch
97         void setMRFFaces();
99         //- Make the given absolute mass/vol flux relative within the MRF region
100         template<class RhoFieldType>
101         void relativeRhoFlux
102         (
103             const RhoFieldType& rho,
104             surfaceScalarField& phi
105         ) const;
107         //- Make the given relative mass/vol flux absolute within the MRF region
108         template<class RhoFieldType>
109         void absoluteRhoFlux
110         (
111             const RhoFieldType& rho,
112             surfaceScalarField& phi
113         ) const;
115         //- Disallow default bitwise copy construct
116         MRFZone(const MRFZone&);
118         //- Disallow default bitwise assignment
119         void operator=(const MRFZone&);
122 public:
124     // Declare name of the class and its debug switch
125     ClassName("MRFZone");
128     // Constructors
130         //- Construct from fvMesh and Istream
131         MRFZone(const fvMesh& mesh, Istream& is);
133         //- Return clone
134         autoPtr<MRFZone> clone() const
135         {
136             notImplemented("autoPtr<MRFZone> clone() const");
137             return autoPtr<MRFZone>(NULL);
138         }
140         //- Return a pointer to a new MRFZone created on freestore
141         //  from Istream
142         class iNew
143         {
144             const fvMesh& mesh_;
146         public:
148             iNew(const fvMesh& mesh)
149             :
150                 mesh_(mesh)
151             {}
153             autoPtr<MRFZone> operator()(Istream& is) const
154             {
155                 return autoPtr<MRFZone>(new MRFZone(mesh_, is));
156             }
157         };
160     // Member Functions
162         //- Update the mesh corresponding to given map
163         void updateMesh(const mapPolyMesh& mpm)
164         {
165             // Only updates face addressing
166             setMRFFaces();
167         }
169         //- Add the Coriolis force contribution to the momentum equation
170         void addCoriolis(fvVectorMatrix& UEqn) const;
172         //- Add the Coriolis force contribution to the momentum equation
173         void addCoriolis(const volScalarField& rho, fvVectorMatrix& UEqn) const;
175         //- Make the given absolute velocity relative within the MRF region
176         void relativeVelocity(volVectorField& U) const;
178         //- Make the given relative velocity absolute within the MRF region
179         void absoluteVelocity(volVectorField& U) const;
181         //- Make the given absolute flux relative within the MRF region
182         void relativeFlux(surfaceScalarField& phi) const;
184         //- Make the given absolute mass-flux relative within the MRF region
185         void relativeFlux
186         (
187             const surfaceScalarField& rho,
188             surfaceScalarField& phi
189         ) const;
191         //- Make the given relative flux absolute within the MRF region
192         void absoluteFlux(surfaceScalarField& phi) const;
194         //- Make the given relative mass-flux absolute within the MRF region
195         void absoluteFlux
196         (
197             const surfaceScalarField& rho,
198             surfaceScalarField& phi
199         ) const;
201         //- Correct the boundary velocity for the roation of the MRF region
202         void correctBoundaryVelocity(volVectorField& U) const;
205     // Ostream Operator
207         friend Ostream& operator<<(Ostream& os, const MRFZone&)
208         {
209             notImplemented("Ostream& operator<<(Ostream& os, const MRFZone&)");
210             return os;
211         }
215 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
217 } // End namespace Foam
219 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
221 #ifdef NoRepository
222 #   include "MRFZoneTemplates.C"
223 #endif
225 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
227 #endif
229 // ************************************************************************* //