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-2016 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 "Function1.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
77  dictionary coeffs_;
78 
79  //- MRF region active flag
80  bool active_;
81 
82  //- Name of cell zone
83  word cellZoneName_;
84 
85  //- Cell zone ID
86  label cellZoneID_;
87 
88  const wordReList excludedPatchNames_;
89 
90  labelList excludedPatchLabels_;
91 
92  //- Internal faces that are part of MRF
93  labelList internalFaces_;
94 
95  //- Outside faces (per patch) that move with the MRF
96  labelListList includedFaces_;
97 
98  //- Excluded faces (per patch) that do not move with the MRF
99  labelListList excludedFaces_;
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>
118  void makeRelativeRhoFlux
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>
126  void makeRelativeRhoFlux
127  (
128  const RhoFieldType& rho,
130  ) const;
131 
132  //- Make the given absolute mass/vol flux relative within the MRF region
133  template<class RhoFieldType>
134  void makeRelativeRhoFlux
135  (
136  const RhoFieldType& rho,
137  Field<scalar>& phi,
138  const label patchi
139  ) const;
140 
141  //- Make the given relative mass/vol flux absolute within the MRF region
142  template<class RhoFieldType>
143  void makeAbsoluteRhoFlux
144  (
145  const RhoFieldType& rho,
146  surfaceScalarField& phi
147  ) const;
148 
149  //- Disallow default bitwise copy construct
150  MRFZone(const MRFZone&);
151 
152  //- Disallow default bitwise assignment
153  void operator=(const MRFZone&);
154 
155 
156 public:
157 
158  // Declare name of the class and its debug switch
159  ClassName("MRFZone");
160 
161 
162  // Constructors
163 
164  //- Construct from fvMesh
165  MRFZone
166  (
167  const word& name,
168  const fvMesh& mesh,
169  const dictionary& dict,
170  const word& cellZoneName = word::null
171  );
172 
173  //- Return clone
174  autoPtr<MRFZone> clone() const
175  {
177  return autoPtr<MRFZone>(nullptr);
178  }
179 
180 
181  // Member Functions
182 
183  // Access
184 
185  //- Return const access to the MRF region name
186  inline const word& name() const;
187 
188  //- Return const access to the MRF active flag
189  inline bool active() const;
190 
191  //- Return the current Omega vector
192  vector Omega() const;
193 
194 
195  // Evaluation
196 
197  //- Update the mesh corresponding to given map
198  void updateMesh(const mapPolyMesh& mpm)
199  {
200  // Only updates face addressing
201  setMRFFaces();
202  }
203 
204  //- Add the Coriolis force contribution to the acceleration field
205  void addCoriolis
206  (
207  const volVectorField& U,
208  volVectorField& ddtU
209  ) const;
210 
211  //- Add the Coriolis force contribution to the momentum equation
212  // Adds to the lhs of the equation; optionally add to rhs
213  void addCoriolis
214  (
216  const bool rhs = false
217  ) const;
218 
219  //- Add the Coriolis force contribution to the momentum equation
220  // Adds to the lhs of the equation; optionally add to rhs
221  void addCoriolis
222  (
223  const volScalarField& rho,
224  fvVectorMatrix& UEqn,
225  const bool rhs = false
226  ) const;
227 
228  //- Make the given absolute velocity relative within the MRF region
229  void makeRelative(volVectorField& U) const;
230 
231  //- Make the given absolute flux relative within the MRF region
232  void makeRelative(surfaceScalarField& phi) const;
233 
234  //- Make the given absolute boundary flux relative
235  // within the MRF region
237 
238  //- Make the given absolute patch flux relative
239  // within the MRF region
240  void makeRelative(Field<scalar>& phi, const label patchi) const;
241 
242  //- Make the given absolute mass-flux relative within the MRF region
243  void makeRelative
244  (
245  const surfaceScalarField& rho,
246  surfaceScalarField& phi
247  ) const;
248 
249  //- Make the given relative velocity absolute within the MRF region
250  void makeAbsolute(volVectorField& U) const;
251 
252  //- Make the given relative flux absolute within the MRF region
253  void makeAbsolute(surfaceScalarField& phi) const;
254 
255  //- Make the given relative mass-flux absolute within the MRF region
256  void makeAbsolute
257  (
258  const surfaceScalarField& rho,
259  surfaceScalarField& phi
260  ) const;
261 
262  //- Correct the boundary velocity for the rotation of the MRF region
264 
265 
266  // I-O
267 
268  //- Write
269  void writeData(Ostream& os) const;
270 
271  //- Read MRF dictionary
272  bool read(const dictionary& dict);
273 };
274 
275 
276 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
277 
278 } // End namespace Foam
279 
280 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
281 
282 #ifdef NoRepository
283  #include "MRFZoneTemplates.C"
284 #endif
285 
286 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
287 
288 #include "MRFZoneI.H"
289 
290 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
291 
292 #endif
293 
294 // ************************************************************************* //
bool read(const dictionary &dict)
Read MRF dictionary.
Definition: MRFZone.C:592
Foam::surfaceFields.
bool active() const
Return const access to the MRF active flag.
Definition: MRFZoneI.H:32
surfaceScalarField & phi
void addCoriolis(const volVectorField &U, volVectorField &ddtU) const
Add the Coriolis force contribution to the acceleration field.
Definition: MRFZone.C:311
void makeRelative(volVectorField &U) const
Make the given absolute velocity relative within the MRF region.
Definition: MRFZone.C:406
dictionary dict
U
Definition: pEqn.H:83
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
autoPtr< MRFZone > clone() const
Return clone.
Definition: MRFZone.H:173
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
void updateMesh(const mapPolyMesh &mpm)
Update the mesh corresponding to given map.
Definition: MRFZone.H:197
MRF zone definition based on cell zone and parameters obtained from a control dictionary constructed ...
Definition: MRFZone.H:65
void makeAbsolute(volVectorField &U) const
Make the given relative velocity absolute within the MRF region.
Definition: MRFZone.C:480
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
Generic field type.
Definition: FieldField.H:51
dynamicFvMesh & mesh
A class for handling words, derived from string.
Definition: word.H:59
const word & name() const
Return const access to the MRF region name.
Definition: MRFZoneI.H:26
static const word null
An empty word.
Definition: word.H:77
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:72
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
void correctBoundaryVelocity(volVectorField &U) const
Correct the boundary velocity for the rotation of the MRF region.
Definition: MRFZone.C:541
vector Omega() const
Return the current Omega vector.
Definition: MRFZone.C:304
Forward declarations of fvMatrix specializations.
label patchi
fvVectorMatrix & UEqn
Definition: UEqn.H:13
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
void writeData(Ostream &os) const
Write.
Definition: MRFZone.C:571
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:366
ClassName("MRFZone")
Namespace for OpenFOAM.