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