MRFZone.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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.
15 
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.
20 
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/>.
23 
24 Class
25  Foam::MRFZone
26 
27 Description
28  MRF zone definition based on cell zone and parameters
29  obtained from a control dictionary constructed from the given stream.
30 
31  The rotation of the MRF region is defined by an origin and axis of
32  rotation and an angular speed.
33 
34 SourceFiles
35  MRFZone.C
36 
37 \*---------------------------------------------------------------------------*/
38 
39 #ifndef MRFZone_H
40 #define MRFZone_H
41 
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 "surfaceFields.H"
49 #include "fvMatricesFwd.H"
50 #include "mapPolyMesh.H"
51 #include "DataEntry.H"
52 #include "autoPtr.H"
53 
54 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
55 
56 namespace Foam
57 {
58 
59 // Forward declaration of classes
60 class fvMesh;
61 
62 /*---------------------------------------------------------------------------*\
63  Class MRFZone Declaration
64 \*---------------------------------------------------------------------------*/
65 
66 class MRFZone
67 {
68  // Private data
69 
70  //- Reference to the mesh database
71  const fvMesh& mesh_;
72 
73  //- Name of the MRF region
74  const word name_;
75 
76  //- Coefficients dictionary
78 
79  //- MRF region active flag
80  bool active_;
81 
82  //- Name of cell zone
84 
85  //- Cell zone ID
87 
89 
91 
92  //- Internal faces that are part of MRF
94 
95  //- Outside faces (per patch) that move with the MRF
97 
98  //- Excluded faces (per patch) that do not move with the MRF
100 
101  //- Origin of the axis
102  const vector origin_;
103 
104  //- Axis vector
105  vector axis_;
106 
107  //- Angular velocty (rad/sec)
109 
110 
111  // Private Member Functions
112 
113  //- Divide faces in frame according to patch
114  void setMRFFaces();
115 
116  //- Make the given absolute mass/vol flux relative within the MRF region
117  template<class RhoFieldType>
119  (
120  const RhoFieldType& rho,
122  ) const;
123 
124  //- Make the given absolute mass/vol flux relative within the MRF region
125  template<class RhoFieldType>
127  (
128  const RhoFieldType& rho,
130  ) const;
131 
132  //- Make the given relative mass/vol flux absolute within the MRF region
133  template<class RhoFieldType>
135  (
136  const RhoFieldType& rho,
138  ) const;
139 
140  //- Disallow default bitwise copy construct
141  MRFZone(const MRFZone&);
142 
143  //- Disallow default bitwise assignment
144  void operator=(const MRFZone&);
145 
146 
147 public:
148 
149  // Declare name of the class and its debug switch
150  ClassName("MRFZone");
151 
152 
153  // Constructors
154 
155  //- Construct from fvMesh
156  MRFZone
157  (
158  const word& name,
159  const fvMesh& mesh,
160  const dictionary& dict,
161  const word& cellZoneName = word::null
162  );
163 
164  //- Return clone
165  autoPtr<MRFZone> clone() const
166  {
168  return autoPtr<MRFZone>(NULL);
169  }
170 
171 
172  // Member Functions
173 
174  // Access
175 
176  //- Return const access to the MRF region name
177  inline const word& name() const;
178 
179  //- Return const access to the MRF active flag
180  inline bool active() const;
181 
182  //- Return the current Omega vector
183  vector Omega() const;
184 
185 
186  // Evaluation
187 
188  //- Update the mesh corresponding to given map
189  void updateMesh(const mapPolyMesh& mpm)
190  {
191  // Only updates face addressing
192  setMRFFaces();
193  }
194 
195  //- Add the Coriolis force contribution to the acceleration field
196  void addCoriolis
197  (
198  const volVectorField& U,
199  volVectorField& ddtU
200  ) const;
201 
202  //- Add the Coriolis force contribution to the momentum equation
203  // Adds to the lhs of the equation; optionally add to rhs
204  void addCoriolis
205  (
207  const bool rhs = false
208  ) const;
209 
210  //- Add the Coriolis force contribution to the momentum equation
211  // Adds to the lhs of the equation; optionally add to rhs
212  void addCoriolis
213  (
214  const volScalarField& rho,
216  const bool rhs = false
217  ) const;
218 
219  //- Make the given absolute velocity relative within the MRF region
220  void makeRelative(volVectorField& U) const;
221 
222  //- Make the given absolute flux relative within the MRF region
223  void makeRelative(surfaceScalarField& phi) const;
224 
225  //- Make the given absolute boundary flux relative
226  // within the MRF region
228 
229  //- Make the given absolute mass-flux relative within the MRF region
230  void makeRelative
231  (
232  const surfaceScalarField& rho,
234  ) const;
235 
236  //- Make the given relative velocity absolute within the MRF region
237  void makeAbsolute(volVectorField& U) const;
238 
239  //- Make the given relative flux absolute within the MRF region
240  void makeAbsolute(surfaceScalarField& phi) const;
241 
242  //- Make the given relative mass-flux absolute within the MRF region
243  void makeAbsolute
244  (
245  const surfaceScalarField& rho,
247  ) const;
248 
249  //- Correct the boundary velocity for the rotation of the MRF region
251 
252 
253  // I-O
254 
255  //- Write
256  void writeData(Ostream& os) const;
257 
258  //- Read MRF dictionary
259  bool read(const dictionary& dict);
260 };
261 
262 
263 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
264 
265 } // End namespace Foam
266 
267 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
268 
269 #ifdef NoRepository
270 # include "MRFZoneTemplates.C"
271 #endif
272 
273 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
274 
275 #include "MRFZoneI.H"
276 
277 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
278 
279 #endif
280 
281 // ************************************************************************* //
Foam::MRFZone::name_
const word name_
Name of the MRF region.
Definition: MRFZone.H:73
Foam::MRFZone
MRF zone definition based on cell zone and parameters obtained from a control dictionary constructed ...
Definition: MRFZone.H:65
volFieldsFwd.H
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::FieldField
Generic field type.
Definition: FieldField.H:51
Foam::MRFZone::coeffs_
dictionary coeffs_
Coefficients dictionary.
Definition: MRFZone.H:76
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::MRFZone::origin_
const vector origin_
Origin of the axis.
Definition: MRFZone.H:101
Foam::MRFZone::MRFZone
MRFZone(const MRFZone &)
Disallow default bitwise copy construct.
mapPolyMesh.H
Foam::MRFZone::updateMesh
void updateMesh(const mapPolyMesh &mpm)
Update the mesh corresponding to given map.
Definition: MRFZone.H:188
Foam::MRFZone::setMRFFaces
void setMRFFaces()
Divide faces in frame according to patch.
Definition: MRFZone.C:45
Foam::MRFZone::makeAbsoluteRhoFlux
void makeAbsoluteRhoFlux(const RhoFieldType &rho, surfaceScalarField &phi) const
Make the given relative mass/vol flux absolute within the MRF region.
Definition: MRFZoneTemplates.C:102
Foam::MRFZone::internalFaces_
labelList internalFaces_
Internal faces that are part of MRF.
Definition: MRFZone.H:92
Foam::MRFZone::makeRelativeRhoFlux
void makeRelativeRhoFlux(const RhoFieldType &rho, surfaceScalarField &phi) const
Make the given absolute mass/vol flux relative within the MRF region.
Definition: MRFZoneTemplates.C:36
Foam::MRFZone::makeAbsolute
void makeAbsolute(volVectorField &U) const
Make the given relative velocity absolute within the MRF region.
Definition: MRFZone.C:466
fvMatricesFwd.H
Forward declarations of fvMatrix specializations.
surfaceFields.H
Foam::surfaceFields.
wordList.H
U
U
Definition: pEqn.H:46
Foam::MRFZone::axis_
vector axis_
Axis vector.
Definition: MRFZone.H:104
Foam::MRFZone::read
bool read(const dictionary &dict)
Read MRF dictionary.
Definition: MRFZone.C:564
NotImplemented
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:365
DataEntry.H
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::MRFZone::excludedPatchNames_
const wordReList excludedPatchNames_
Definition: MRFZone.H:87
labelList.H
UEqn
tmp< fvVectorMatrix > UEqn(fvm::div(phi, U)+MRF.DDt(U)+turbulence->divDevReff(U)==fvOptions(U))
Foam::MRFZone::active
bool active() const
Return const access to the MRF active flag.
Definition: MRFZoneI.H:32
dict
dictionary dict
Definition: searchingEngine.H:14
dimensionedVector.H
MRFZoneI.H
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::MRFZone::writeData
void writeData(Ostream &os) const
Write.
Definition: MRFZone.C:543
Foam::MRFZone::active_
bool active_
MRF region active flag.
Definition: MRFZone.H:79
Foam::MRFZone::cellZoneID_
label cellZoneID_
Cell zone ID.
Definition: MRFZone.H:85
Foam::MRFZone::Omega
vector Omega() const
Return the current Omega vector.
Definition: MRFZone.C:304
rho
rho
Definition: pEqn.H:3
Foam::MRFZone::correctBoundaryVelocity
void correctBoundaryVelocity(volVectorField &U) const
Correct the boundary velocity for the rotation of the MRF region.
Definition: MRFZone.C:520
Foam::MRFZone::excludedPatchLabels_
labelList excludedPatchLabels_
Definition: MRFZone.H:89
dimensionedScalar.H
Foam::autoPtr
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
Foam::MRFZone::clone
autoPtr< MRFZone > clone() const
Return clone.
Definition: MRFZone.H:164
Foam::MRFZone::makeRelative
void makeRelative(volVectorField &U) const
Make the given absolute velocity relative within the MRF region.
Definition: MRFZone.C:406
Foam::Vector< scalar >
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::MRFZone::ClassName
ClassName("MRFZone")
Foam::MRFZone::excludedFaces_
labelListList excludedFaces_
Excluded faces (per patch) that do not move with the MRF.
Definition: MRFZone.H:98
dictionary.H
Foam::MRFZone::mesh_
const fvMesh & mesh_
Reference to the mesh database.
Definition: MRFZone.H:70
Foam::MRFZone::includedFaces_
labelListList includedFaces_
Outside faces (per patch) that move with the MRF.
Definition: MRFZone.H:95
Foam::word::null
static const word null
An empty word.
Definition: word.H:77
Foam::MRFZone::cellZoneName_
word cellZoneName_
Name of cell zone.
Definition: MRFZone.H:82
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
Foam::MRFZone::operator=
void operator=(const MRFZone &)
Disallow default bitwise assignment.
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::MRFZone::name
const word & name() const
Return const access to the MRF region name.
Definition: MRFZoneI.H:26
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
MRFZoneTemplates.C
Foam::MRFZone::omega_
autoPtr< DataEntry< scalar > > omega_
Angular velocty (rad/sec)
Definition: MRFZone.H:107
Foam::MRFZone::addCoriolis
void addCoriolis(const volVectorField &U, volVectorField &ddtU) const
Add the Coriolis force contribution to the acceleration field.
Definition: MRFZone.C:311
autoPtr.H