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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2018 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Class
27  Foam::MRFZone
28 
29 Description
30  MRF zone definition based on cell zone and parameters
31  obtained from a control dictionary constructed from the given stream.
32 
33  The rotation of the MRF region is defined by an origin and axis of
34  rotation and an angular speed.
35 
36 SourceFiles
37  MRFZone.C
38 
39 \*---------------------------------------------------------------------------*/
40 
41 #ifndef MRFZone_H
42 #define MRFZone_H
43 
44 #include "dictionary.H"
45 #include "wordList.H"
46 #include "labelList.H"
47 #include "dimensionedScalar.H"
48 #include "dimensionedVector.H"
49 #include "volFieldsFwd.H"
50 #include "surfaceFields.H"
51 #include "fvMatricesFwd.H"
52 #include "mapPolyMesh.H"
53 #include "Function1.H"
54 #include "autoPtr.H"
55 
56 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
57 
58 namespace Foam
59 {
60 
61 // Forward declaration of classes
62 class fvMesh;
63 
64 /*---------------------------------------------------------------------------*\
65  Class MRFZone Declaration
66 \*---------------------------------------------------------------------------*/
67 
68 class MRFZone
69 {
70  // Private data
71 
72  //- Reference to the mesh database
73  const fvMesh& mesh_;
74 
75  //- Name of the MRF region
76  const word name_;
77 
78  //- Coefficients dictionary
79  dictionary coeffs_;
80 
81  //- MRF region active flag
82  bool active_;
83 
84  //- Name of cell zone
85  word cellZoneName_;
86 
87  //- Cell zone ID
88  label cellZoneID_;
89 
90  const wordRes excludedPatchNames_;
91 
92  labelList excludedPatchLabels_;
93 
94  //- Internal faces that are part of MRF
95  labelList internalFaces_;
96 
97  //- Outside faces (per patch) that move with the MRF
98  labelListList includedFaces_;
99 
100  //- Excluded faces (per patch) that do not move with the MRF
101  labelListList excludedFaces_;
102 
103  //- Origin of the axis
104  const vector origin_;
105 
106  //- Axis vector
107  vector axis_;
108 
109  //- Angular velocity (rad/sec)
111 
112 
113  // Private Member Functions
114 
115  //- Divide faces in frame according to patch
116  void setMRFFaces();
117 
118  //- Make the given absolute mass/vol flux relative within the MRF region
119  template<class RhoFieldType>
120  void makeRelativeRhoFlux
121  (
122  const RhoFieldType& rho,
124  ) const;
125 
126  //- Make the given absolute mass/vol flux relative within the MRF region
127  template<class RhoFieldType>
128  void makeRelativeRhoFlux
129  (
130  const RhoFieldType& rho,
132  ) const;
133 
134  //- Make the given absolute mass/vol flux relative within the MRF region
135  template<class RhoFieldType>
136  void makeRelativeRhoFlux
137  (
138  const RhoFieldType& rho,
140  const label patchi
141  ) const;
142 
143  //- Make the given relative mass/vol flux absolute within the MRF region
144  template<class RhoFieldType>
145  void makeAbsoluteRhoFlux
146  (
147  const RhoFieldType& rho,
149  ) const;
150 
151  //- No copy construct
152  MRFZone(const MRFZone&) = delete;
153 
154  //- No copy assignment
155  void operator=(const MRFZone&) = delete;
156 
157 
158 public:
159 
160  // Declare name of the class and its debug switch
161  ClassName("MRFZone");
162 
163 
164  // Constructors
165 
166  //- Construct from fvMesh
167  MRFZone
168  (
169  const word& name,
170  const fvMesh& mesh,
171  const dictionary& dict,
172  const word& cellZoneName = word::null
173  );
174 
175  //- Return clone
176  autoPtr<MRFZone> clone() const
177  {
179  return nullptr;
180  }
181 
182 
183  // Member Functions
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  //- Update the mesh corresponding to given map
195  void updateMesh(const mapPolyMesh& mpm)
196  {
197  // Only updates face addressing
198  setMRFFaces();
199  }
200 
201  //- Add the Coriolis force contribution to the acceleration field
202  void addCoriolis
203  (
204  const volVectorField& U,
205  volVectorField& ddtU
206  ) const;
207 
208  //- Add the Coriolis force contribution to the momentum equation
209  // Adds to the lhs of the equation; optionally add to rhs
210  void addCoriolis
211  (
213  const bool rhs = false
214  ) const;
215 
216  //- Add the Coriolis force contribution to the momentum equation
217  // Adds to the lhs of the equation; optionally add to rhs
218  void addCoriolis
219  (
220  const volScalarField& rho,
222  const bool rhs = false
223  ) const;
224 
225  //- Make the given absolute velocity relative within the MRF region
226  void makeRelative(volVectorField& U) const;
227 
228  //- Make the given absolute flux relative within the MRF region
229  void makeRelative(surfaceScalarField& phi) const;
230 
231  //- Make the given absolute boundary flux relative
232  // within the MRF region
234 
235  //- Make the given absolute patch flux relative
236  // within the MRF region
237  void makeRelative(Field<scalar>& phi, const label patchi) const;
238 
239  //- Make the given absolute mass-flux relative within the MRF region
240  void makeRelative
241  (
242  const surfaceScalarField& rho,
244  ) const;
245 
246  //- Make the given relative velocity absolute within the MRF region
247  void makeAbsolute(volVectorField& U) const;
248 
249  //- Make the given relative flux absolute within the MRF region
250  void makeAbsolute(surfaceScalarField& phi) const;
251 
252  //- Make the given relative mass-flux absolute within the MRF region
253  void makeAbsolute
254  (
255  const surfaceScalarField& rho,
257  ) const;
258 
259  //- Correct the boundary velocity for the rotation of the MRF region
261 
262  //- Zero the MRF region of the given field
263  template<class Type>
265 
266  //- Update MRFZone faces if the mesh topology changes
267  void update();
268 
269 
270  // I-O
271 
272  //- Write
273  void writeData(Ostream& os) const;
274 
275  //- Read MRF dictionary
276  bool read(const dictionary& dict);
277 };
278 
279 
280 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
281 
282 } // End namespace Foam
283 
284 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
285 
286 #ifdef NoRepository
287  #include "MRFZoneTemplates.C"
288 #endif
289 
290 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
291 
292 #include "MRFZoneI.H"
293 
294 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
295 
296 #endif
297 
298 // ************************************************************************* //
Foam::MRFZone
MRF zone definition based on cell zone and parameters obtained from a control dictionary constructed ...
Definition: MRFZone.H:67
volFieldsFwd.H
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::MRFZone::zero
void zero(GeometricField< Type, fvsPatchField, surfaceMesh > &phi) const
Zero the MRF region of the given field.
Definition: MRFZoneTemplates.C:214
Foam::FieldField
A field of fields is a PtrList of fields with reference counting.
Definition: FieldField.H:53
Function1.H
mapPolyMesh.H
Foam::MRFZone::updateMesh
void updateMesh(const mapPolyMesh &mpm)
Update the mesh corresponding to given map.
Definition: MRFZone.H:194
Foam::MRFZone::makeAbsolute
void makeAbsolute(volVectorField &U) const
Make the given relative velocity absolute within the MRF region.
Definition: MRFZone.C:481
fvMatricesFwd.H
Forward declarations of fvMatrix specializations.
surfaceFields.H
Foam::surfaceFields.
rho
rho
Definition: readInitialConditions.H:88
wordList.H
Foam::MRFZone::read
bool read(const dictionary &dict)
Read MRF dictionary.
Definition: MRFZone.C:592
NotImplemented
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:445
labelList.H
Foam::Field< scalar >
Foam::MRFZone::active
bool active() const
Return const access to the MRF active flag.
Definition: MRFZoneI.H:34
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
dict
dictionary dict
Definition: searchingEngine.H:14
dimensionedVector.H
MRFZoneI.H
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:83
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::MRFZone::writeData
void writeData(Ostream &os) const
Write.
Definition: MRFZone.C:572
Foam::MRFZone::Omega
vector Omega() const
Return the current Omega vector.
Definition: MRFZone.C:305
Foam::MRFZone::correctBoundaryVelocity
void correctBoundaryVelocity(volVectorField &U) const
Correct the boundary velocity for the rotation of the MRF region.
Definition: MRFZone.C:542
dimensionedScalar.H
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
U
U
Definition: pEqn.H:72
Foam::MRFZone::clone
autoPtr< MRFZone > clone() const
Return clone.
Definition: MRFZone.H:175
Foam::MRFZone::makeRelative
void makeRelative(volVectorField &U) const
Make the given absolute velocity relative within the MRF region.
Definition: MRFZone.C:407
Foam::Vector< scalar >
Foam::List< label >
Foam::MRFZone::ClassName
ClassName("MRFZone")
dictionary.H
Foam::wordRes
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:51
Foam::word::null
static const word null
An empty word.
Definition: word.H:77
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:76
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:161
UEqn
fvVectorMatrix & UEqn
Definition: UEqn.H:13
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::MRFZone::name
const word & name() const
Return const access to the MRF region name.
Definition: MRFZoneI.H:28
Foam::GeometricField< scalar, fvsPatchField, surfaceMesh >
MRFZoneTemplates.C
Foam::MRFZone::update
void update()
Update MRFZone faces if the mesh topology changes.
Definition: MRFZone.C:604
Foam::MRFZone::addCoriolis
void addCoriolis(const volVectorField &U, volVectorField &ddtU) const
Add the Coriolis force contribution to the acceleration field.
Definition: MRFZone.C:312
autoPtr.H