MRFZone.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2020 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  //- Name of cell zone
80  word cellZoneName_;
81 
82  //- Cell zone ID
83  label cellZoneID_;
84 
85  const wordReList excludedPatchNames_;
86 
87  labelList excludedPatchLabels_;
88 
89  //- Internal faces that are part of MRF
90  labelList internalFaces_;
91 
92  //- Outside faces (per patch) that move with the MRF
93  labelListList includedFaces_;
94 
95  //- Excluded faces (per patch) that do not move with the MRF
96  labelListList excludedFaces_;
97 
98  //- Origin of the axis
99  const vector origin_;
100 
101  //- Axis vector
102  vector axis_;
103 
104  //- Angular velocity (rad/sec)
106 
107 
108  // Private Member Functions
109 
110  //- Divide faces in frame according to patch
111  void setMRFFaces();
112 
113  //- Make the given absolute mass/vol flux relative within the MRF region
114  template<class RhoFieldType>
115  void makeRelativeRhoFlux
116  (
117  const RhoFieldType& rho,
119  ) const;
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,
134  Field<scalar>& phi,
135  const label patchi
136  ) const;
137 
138  //- Make the given relative mass/vol flux absolute within the MRF region
139  template<class RhoFieldType>
140  void makeAbsoluteRhoFlux
141  (
142  const RhoFieldType& rho,
143  surfaceScalarField& phi
144  ) const;
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  //- Disallow default bitwise copy construction
165  MRFZone(const MRFZone&) = delete;
166 
167  //- Return clone
168  autoPtr<MRFZone> clone() const
169  {
171  return autoPtr<MRFZone>(nullptr);
172  }
173 
174 
175  // Member Functions
176 
177  //- Return const access to the MRF region name
178  inline const word& name() const;
179 
180  //- Return the current Omega vector
181  vector Omega() const;
182 
183  //- Update the mesh corresponding to given map
184  void updateMesh(const mapPolyMesh& mpm)
185  {
186  // Only updates face addressing
187  setMRFFaces();
188  }
189 
190  //- Add the Coriolis force contribution to the acceleration field
191  void addCoriolis
192  (
193  const volVectorField& U,
194  volVectorField& ddtU
195  ) const;
196 
197  //- Add the Coriolis force contribution to the momentum equation
198  // Adds to the lhs of the equation; optionally add to rhs
199  void addCoriolis
200  (
202  const bool rhs = false
203  ) const;
204 
205  //- Add the Coriolis force contribution to the momentum equation
206  // Adds to the lhs of the equation; optionally add to rhs
207  void addCoriolis
208  (
209  const volScalarField& rho,
210  fvVectorMatrix& UEqn,
211  const bool rhs = false
212  ) const;
213 
214  //- Make the given absolute velocity relative within the MRF region
215  void makeRelative(volVectorField& U) const;
216 
217  //- Make the given absolute flux relative within the MRF region
218  void makeRelative(surfaceScalarField& phi) const;
219 
220  //- Make the given absolute boundary flux relative
221  // within the MRF region
223 
224  //- Make the given absolute patch flux relative
225  // within the MRF region
226  void makeRelative(Field<scalar>& phi, const label patchi) const;
227 
228  //- Make the given absolute mass-flux relative within the MRF region
229  void makeRelative
230  (
231  const surfaceScalarField& rho,
232  surfaceScalarField& phi
233  ) const;
234 
235  //- Make the given relative velocity absolute within the MRF region
236  void makeAbsolute(volVectorField& U) const;
237 
238  //- Make the given relative flux absolute within the MRF region
239  void makeAbsolute(surfaceScalarField& phi) const;
240 
241  //- Make the given relative mass-flux absolute within the MRF region
242  void makeAbsolute
243  (
244  const surfaceScalarField& rho,
245  surfaceScalarField& phi
246  ) const;
247 
248  //- Correct the boundary velocity for the rotation of the MRF region
250 
251  //- Zero the MRF region of the given field
252  template<class Type>
254 
255  //- Update MRFZone faces if the mesh topology changes
256  void update();
257 
258 
259  // I-O
260 
261  //- Write
262  void writeData(Ostream& os) const;
263 
264  //- Read MRF dictionary
265  bool read(const dictionary& dict);
266 
267 
268  // Member Operators
269 
270  //- Disallow default bitwise assignment
271  void operator=(const MRFZone&) = delete;
272 };
273 
274 
275 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
276 
277 } // End namespace Foam
278 
279 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
280 
281 #ifdef NoRepository
282  #include "MRFZoneTemplates.C"
283 #endif
284 
285 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
286 
287 #include "MRFZoneI.H"
288 
289 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
290 
291 #endif
292 
293 // ************************************************************************* //
bool read(const dictionary &dict)
Read MRF dictionary.
Definition: MRFZone.C:577
Foam::surfaceFields.
void addCoriolis(const volVectorField &U, volVectorField &ddtU) const
Add the Coriolis force contribution to the acceleration field.
Definition: MRFZone.C:303
void makeRelative(volVectorField &U) const
Make the given absolute velocity relative within the MRF region.
Definition: MRFZone.C:398
dictionary dict
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:167
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
MRFZone(const word &name, const fvMesh &mesh, const dictionary &dict, const word &cellZoneName=word::null)
Construct from fvMesh.
Definition: MRFZone.C:237
void updateMesh(const mapPolyMesh &mpm)
Update the mesh corresponding to given map.
Definition: MRFZone.H:183
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:472
void update()
Update MRFZone faces if the mesh topology changes.
Definition: MRFZone.C:588
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
void zero(GeometricField< Type, fvsPatchField, surfaceMesh > &phi) const
Zero the MRF region of the given field.
static const word null
An empty word.
Definition: word.H:77
phi
Definition: correctPhi.H:3
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
void operator=(const MRFZone &)=delete
Disallow default bitwise assignment.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
void correctBoundaryVelocity(volVectorField &U) const
Correct the boundary velocity for the rotation of the MRF region.
Definition: MRFZone.C:533
vector Omega() const
Return the current Omega vector.
Definition: MRFZone.C:296
Forward declarations of fvMatrix specialisations.
label patchi
U
Definition: pEqn.H:72
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:558
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:370
ClassName("MRFZone")
Namespace for OpenFOAM.
fvVectorMatrix & UEqn
Definition: UEqn.H:13