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-2021 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 zomes
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
49 class MRFZoneList;
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  //- Return active status
85  bool active(const bool warn = false) const;
86 
87  //- Reset the source list
88  void reset(const dictionary& dict);
89 
90  //- Add the frame acceleration
91  void addAcceleration
92  (
93  const volVectorField& U,
94  volVectorField& ddtU
95  ) const;
96 
97  //- Add the frame acceleration contribution to the momentum equation
98  void addAcceleration(fvVectorMatrix& UEqn) const;
99 
100  //- Add the frame acceleration contribution to the momentum equation
101  void addAcceleration
102  (
103  const volScalarField& rho,
105  ) const;
106 
107  //- Return the frame acceleration
109  (
110  const volVectorField& U
111  ) const;
112 
113  //- Return the frame acceleration
115  (
116  const volScalarField& rho,
117  const volVectorField& U
118  ) const;
119 
120  //- Make the given absolute velocity relative within the MRF region
121  void makeRelative(volVectorField& U) const;
122 
123  //- Make the given absolute flux relative within the MRF region
124  void makeRelative(surfaceScalarField& phi) const;
125 
126  //- Return the given absolute flux relative within the MRF region
128  (
130  ) const;
131 
132  //- Return the given absolute boundary flux relative within
133  // the MRF region
135  (
137  ) const;
138 
139  //- Return the given absolute patch flux relative within
140  // the MRF region
142  (
143  const tmp<Field<scalar>>& tphi,
144  const label patchi
145  ) const;
146 
147  //- Make the given absolute mass-flux relative within the MRF region
148  void makeRelative
149  (
150  const surfaceScalarField& rho,
152  ) const;
153 
154  //- Make the given relative velocity absolute within the MRF region
155  void makeAbsolute(volVectorField& U) const;
156 
157  //- Make the given relative flux absolute within the MRF region
158  void makeAbsolute(surfaceScalarField& phi) const;
159 
160  //- Return the given relative flux absolute within the MRF region
162  (
164  ) const;
165 
166  //- Make the given relative mass-flux absolute within the MRF region
167  void makeAbsolute
168  (
169  const surfaceScalarField& rho,
171  ) const;
172 
173  //- Correct the boundary velocity for the rotation of the MRF region
175 
176  //- Filter-out the MRF region contribution from the given field
177  // setting the corresponding values to zero
178  template<class Type>
180  (
182  ) const;
183 
184  //- Update MRFZone faces if the mesh topology changes
185  void update();
186 
187 
188  // I-O
189 
190  //- Read dictionary
191  bool read(const dictionary& dict);
192 
193  //- Write data to Ostream
194  bool writeData(Ostream& os) const;
195 
196  //- Ostream operator
197  friend Ostream& operator<<
198  (
199  Ostream& os,
200  const MRFZoneList& models
201  );
202 
203 
204  // Member Operators
205 
206  //- Disallow default bitwise assignment
207  void operator=(const MRFZoneList&) = delete;
208 };
209 
210 
211 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
212 
213 } // End namespace Foam
214 
215 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
216 
217 #ifdef NoRepository
218  #include "MRFZoneListTemplates.C"
219 #endif
220 
221 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
222 
223 #endif
224 
225 // ************************************************************************* //
~MRFZoneList()
Destructor.
Definition: MRFZoneList.C:47
void makeAbsolute(volVectorField &U) const
Make the given relative velocity absolute within the MRF region.
Definition: MRFZoneList.C:295
dictionary dict
void update()
Update MRFZone faces if the mesh topology changes.
Definition: MRFZoneList.C:366
void reset(const dictionary &dict)
Reset the source list.
Definition: MRFZoneList.C:53
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
void addAcceleration(const volVectorField &U, volVectorField &ddtU) const
Add the frame acceleration.
Definition: MRFZoneList.C:111
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
void operator=(const MRFZoneList &)=delete
Disallow default bitwise assignment.
const fvMesh & mesh_
Reference to the mesh database.
Definition: MRFZoneList.H:64
Generic field type.
Definition: FieldField.H:51
bool active(const bool warn=false) const
Return active status.
void correctBoundaryVelocity(volVectorField &U) const
Correct the boundary velocity for the rotation of the MRF region.
Definition: MRFZoneList.C:357
dynamicFvMesh & mesh
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
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &phi) const
Return the given relative flux absolute within the MRF region.
Definition: MRFZoneList.C:314
tmp< volVectorField > DDt(const volVectorField &U) const
Return the frame acceleration.
Definition: MRFZoneList.C:146
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
bool writeData(Ostream &os) const
Write data to Ostream.
Definition: MRFZoneList.C:98
List container for MRF zomes.
Definition: MRFZoneList.H:55
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > zeroFilter(const tmp< GeometricField< Type, fvsPatchField, surfaceMesh >> &tphi) const
Filter-out the MRF region contribution from the given field.
tmp< surfaceScalarField > relative(const tmp< surfaceScalarField > &phi) const
Return the given absolute flux relative within the MRF region.
Definition: MRFZoneList.C:199
Forward declarations of fvMatrix specialisations.
label patchi
U
Definition: pEqn.H:72
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Ostream & operator<<(Ostream &, const ensightPart &)
MRFZoneList(const fvMesh &mesh, const dictionary &dict)
Definition: MRFZoneList.C:33
A class for managing temporary objects.
Definition: PtrList.H:53
void makeRelative(volVectorField &U) const
Make the given absolute velocity relative within the MRF region.
Definition: MRFZoneList.C:180
bool read(const dictionary &dict)
Read dictionary.
Definition: MRFZoneList.C:85
Namespace for OpenFOAM.
fvVectorMatrix & UEqn
Definition: UEqn.H:13