All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
MRFZoneList.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) 2012-2023 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::MRFZoneList
26 
27 Description
28  List container for MRF zones
29 
30 SourceFiles
31  MRFZoneList.C
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef MRFZoneList_H
36 #define MRFZoneList_H
37 
38 #include "fvMesh.H"
39 #include "dictionary.H"
40 #include "fvMatricesFwd.H"
41 #include "MRFZone.H"
42 
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 
45 namespace Foam
46 {
47 
48 // Forward declaration of friend functions and operators
50 Ostream& operator<<(Ostream& os, const MRFZoneList& models);
51 
52 /*---------------------------------------------------------------------------*\
53  Class MRFZoneList Declaration
54 \*---------------------------------------------------------------------------*/
55 
56 class MRFZoneList
57 :
58  public PtrList<MRFZone>
59 {
60 protected:
61 
62  // Protected data
63 
64  //- Reference to the mesh database
65  const fvMesh& mesh_;
66 
67 
68 public:
69 
70  // Constructors
71 
72  MRFZoneList(const fvMesh& mesh, const dictionary& dict);
73 
74  //- Disallow default bitwise copy construction
75  MRFZoneList(const MRFZoneList&) = delete;
76 
77 
78  //- Destructor
79  ~MRFZoneList();
80 
81 
82  // Member Functions
83 
84  //- Reset the source list
85  void reset(const dictionary& dict);
86 
87  //- Return the Coriolis acceleration
89  (
90  const volVectorField& U
91  ) const;
92 
93  //- Return the Coriolis acceleration
95  (
96  const volScalarField& rho,
97  const volVectorField& U
98  ) const;
99 
100  //- Return the centrifugal acceleration
102 
103  //- Make the given absolute velocity relative within the MRF region
104  void makeRelative(volVectorField& U) const;
105 
106  //- Make the given absolute flux relative within the MRF region
107  void makeRelative(surfaceScalarField& phi) const;
108 
109  //- Return the given absolute flux relative within the MRF region
111  (
112  const tmp<surfaceScalarField>& phi
113  ) const;
114 
115  //- Return the given absolute boundary flux relative within
116  // the MRF region
118  (
120  ) const;
121 
122  //- Return the given absolute patch flux relative within
123  // the MRF region
125  (
126  const tmp<Field<scalar>>& tphi,
127  const label patchi
128  ) const;
129 
130  //- Make the given absolute mass-flux relative within the MRF region
131  void makeRelative
132  (
133  const surfaceScalarField& rho,
134  surfaceScalarField& phi
135  ) const;
136 
137  //- Make the given relative velocity absolute within the MRF region
138  void makeAbsolute(volVectorField& U) const;
139 
140  //- Make the given relative flux absolute within the MRF region
141  void makeAbsolute(surfaceScalarField& phi) const;
142 
143  //- Return the given relative flux absolute within the MRF region
145  (
146  const tmp<surfaceScalarField>& phi
147  ) const;
148 
149  //- Make the given relative mass-flux absolute within the MRF region
150  void makeAbsolute
151  (
152  const surfaceScalarField& rho,
153  surfaceScalarField& phi
154  ) const;
155 
156  //- Return the given relative flux absolute within the MRF region
158  (
159  const tmp<surfaceScalarField>& phi,
160  const volScalarField& rho
161  ) const;
162 
163  void makeAbsolute(Field<vector>& Up, const label patchi) const;
164 
165  //- Update MRFZone faces if the mesh topology changes
166  void update();
167 
168 
169  // I-O
170 
171  //- Read dictionary
172  bool read(const dictionary& dict);
173 
174 
175  // Member Operators
176 
177  //- Disallow default bitwise assignment
178  void operator=(const MRFZoneList&) = delete;
179 };
180 
181 
182 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
183 
184 } // End namespace Foam
185 
186 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
187 
188 #endif
189 
190 // ************************************************************************* //
Generic field type.
Definition: FieldField.H:77
Generic GeometricField class.
List container for MRF zones.
Definition: MRFZoneList.H:58
tmp< surfaceScalarField > relative(const tmp< surfaceScalarField > &phi) const
Return the given absolute flux relative within the MRF region.
Definition: MRFZoneList.C:176
void reset(const dictionary &dict)
Reset the source list.
Definition: MRFZoneList.C:53
const fvMesh & mesh_
Reference to the mesh database.
Definition: MRFZoneList.H:64
tmp< volVectorField > centrifugalAcceleration() const
Return the centrifugal acceleration.
Definition: MRFZoneList.C:134
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &phi) const
Return the given relative flux absolute within the MRF region.
Definition: MRFZoneList.C:291
bool read(const dictionary &dict)
Read dictionary.
Definition: MRFZoneList.C:85
void makeAbsolute(volVectorField &U) const
Make the given relative velocity absolute within the MRF region.
Definition: MRFZoneList.C:272
void makeRelative(volVectorField &U) const
Make the given absolute velocity relative within the MRF region.
Definition: MRFZoneList.C:157
MRFZoneList(const fvMesh &mesh, const dictionary &dict)
Definition: MRFZoneList.C:33
void update()
Update MRFZone faces if the mesh topology changes.
Definition: MRFZoneList.C:366
void operator=(const MRFZoneList &)=delete
Disallow default bitwise assignment.
~MRFZoneList()
Destructor.
Definition: MRFZoneList.C:47
tmp< volVectorField > DDt(const volVectorField &U) const
Return the Coriolis acceleration.
Definition: MRFZoneList.C:99
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
A class for managing temporary objects.
Definition: tmp.H:55
Forward declarations of fvMatrix specialisations.
label patchi
U
Definition: pEqn.H:72
Namespace for OpenFOAM.
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
Ostream & operator<<(Ostream &os, const fvConstraints &constraints)
dictionary dict