MRFZoneList.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) 2012-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::MRFZoneList
29
30Description
31 List container for MRF zomes
32
33SourceFiles
34 MRFZoneList.C
35
36\*---------------------------------------------------------------------------*/
37
38#ifndef MRFZoneList_H
39#define MRFZoneList_H
40
41#include "fvMesh.H"
42#include "dictionary.H"
43#include "fvMatricesFwd.H"
44#include "MRFZone.H"
45
46// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47
48namespace Foam
49{
50
51// Forward declaration of friend functions and operators
52class MRFZoneList;
53Ostream& operator<<(Ostream& os, const MRFZoneList& models);
54
55/*---------------------------------------------------------------------------*\
56 Class MRFZoneList Declaration
57\*---------------------------------------------------------------------------*/
59class MRFZoneList
60:
61 public PtrList<MRFZone>
62{
63 // Private Member Functions
64
65 //- No copy construct
66 MRFZoneList(const MRFZoneList&) = delete;
67
68 //- No copy assignment
69 void operator=(const MRFZoneList&) = delete;
70
71
72protected:
73
74 // Protected data
75
76 //- Reference to the mesh database
77 const fvMesh& mesh_;
78
79
80public:
81
82 //- Constructor
83 MRFZoneList(const fvMesh& mesh, const dictionary& dict);
84
85 //- Destructor
86 ~MRFZoneList() = default;
87
88
89 // Member Functions
90
91 //- Return active status
92 bool active(const bool warn = false) const;
93
94 //- Reset the source list
95 void reset(const dictionary& dict);
96
97 //- Return the MRF with a given name
98 const MRFZone& getFromName(const word& name) const;
99
100 //- Add the frame acceleration
101 void addAcceleration
102 (
103 const volVectorField& U,
104 volVectorField& ddtU
105 ) const;
106
107 //- Add the frame acceleration contribution to the momentum equation
109
110 //- Add the frame acceleration contribution to the momentum equation
111 void addAcceleration
112 (
113 const volScalarField& rho,
115 ) const;
116
117 //- Return the frame acceleration
119 (
120 const volVectorField& U
121 ) const;
122
123 //- Return the frame acceleration
125 (
126 const volScalarField& rho,
127 const volVectorField& U
128 ) const;
129
130 //- Make the given absolute velocity relative within the MRF region
131 void makeRelative(volVectorField& U) const;
132
133 //- Make the given absolute flux relative within the MRF region
135
136 //- Return the given absolute flux relative within the MRF region
138 (
140 ) const;
141
142 //- Return the given absolute boundary flux relative within
143 // the MRF region
145 (
147 ) const;
148
149 //- Return the given absolute patch flux relative within
150 // the MRF region
152 (
153 const tmp<Field<scalar>>& tphi,
154 const label patchi
155 ) const;
156
157 //- Make the given absolute mass-flux relative within the MRF region
158 void makeRelative
159 (
160 const surfaceScalarField& rho,
162 ) const;
163
164 //- Make the given relative velocity absolute within the MRF region
165 void makeAbsolute(volVectorField& U) const;
166
167 //- Make the given relative flux absolute within the MRF region
169
170 //- Return the given relative flux absolute within the MRF region
172 (
174 ) const;
175
176 //- Make the given relative mass-flux absolute within the MRF region
177 void makeAbsolute
178 (
179 const surfaceScalarField& rho,
181 ) const;
182
183 //- Correct the boundary velocity for the rotation of the MRF region
185
186 //- Correct the boundary flux for the rotation of the MRF region
188 (
189 const volVectorField& U,
191 ) const;
192
193 //- Filter-out the MRF region contribution from the given field
194 // setting the corresponding values to zero
195 template<class Type>
197 (
199 ) const;
200
201 //- Update MRFZone faces if the mesh topology changes
202 void update();
203
204
205 // I-O
206
207 //- Read dictionary
208 bool read(const dictionary& dict);
209
210 //- Write data to Ostream
211 bool writeData(Ostream& os) const;
212
213 //- Ostream operator
214 friend Ostream& operator<<
215 (
216 Ostream& os,
217 const MRFZoneList& models
218 );
219};
220
221// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
222
223} // End namespace Foam
224
225// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
226
227#ifdef NoRepository
228 #include "MRFZoneListTemplates.C"
229#endif
230
231// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
232
233#endif
234
235// ************************************************************************* //
surfaceScalarField & phi
A field of fields is a PtrList of fields with reference counting.
Definition: FieldField.H:80
Generic templated field type.
Definition: Field.H:82
List container for MRF zomes.
Definition: MRFZoneList.H:61
tmp< surfaceScalarField > relative(const tmp< surfaceScalarField > &phi) const
Return the given absolute flux relative within the MRF region.
Definition: MRFZoneList.C:243
void reset(const dictionary &dict)
Reset the source list.
Definition: MRFZoneList.C:69
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > zeroFilter(const tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > &tphi) const
Filter-out the MRF region contribution from the given field.
const fvMesh & mesh_
Reference to the mesh database.
Definition: MRFZoneList.H:76
~MRFZoneList()=default
Destructor.
bool writeData(Ostream &os) const
Write data to Ostream.
Definition: MRFZoneList.C:139
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &phi) const
Return the given relative flux absolute within the MRF region.
Definition: MRFZoneList.C:358
bool read(const dictionary &dict)
Read dictionary.
Definition: MRFZoneList.C:127
void addAcceleration(const volVectorField &U, volVectorField &ddtU) const
Add the frame acceleration.
Definition: MRFZoneList.C:152
void makeAbsolute(volVectorField &U) const
Make the given relative velocity absolute within the MRF region.
Definition: MRFZoneList.C:339
void correctBoundaryVelocity(volVectorField &U) const
Correct the boundary velocity for the rotation of the MRF region.
Definition: MRFZoneList.C:401
const MRFZone & getFromName(const word &name) const
Return the MRF with a given name.
Definition: MRFZoneList.C:103
void makeRelative(volVectorField &U) const
Make the given absolute velocity relative within the MRF region.
Definition: MRFZoneList.C:224
void update()
Update MRFZone faces if the mesh topology changes.
Definition: MRFZoneList.C:433
tmp< volVectorField > DDt(const volVectorField &U) const
Return the frame acceleration.
Definition: MRFZoneList.C:187
bool active(const bool warn=false) const
Return active status.
Definition: MRFZoneList.C:52
void correctBoundaryFlux(const volVectorField &U, surfaceScalarField &phi) const
Correct the boundary flux for the rotation of the MRF region.
Definition: MRFZoneList.C:411
MRF zone definition based on cell zone and parameters obtained from a control dictionary constructed ...
Definition: MRFZone.H:69
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:121
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
U
Definition: pEqn.H:72
fvVectorMatrix & UEqn
Definition: UEqn.H:13
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
Forward declarations of fvMatrix specializations.
Namespace for OpenFOAM.
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:83
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
dictionary dict