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-2018 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  // Private Member Functions
61 
62  //- Disallow default bitwise copy construct
63  MRFZoneList(const MRFZoneList&);
64 
65  //- Disallow default bitwise assignment
66  void operator=(const MRFZoneList&);
67 
68 
69 protected:
70 
71  // Protected data
72 
73  //- Reference to the mesh database
74  const fvMesh& mesh_;
75 
76 
77 public:
78 
79  //- Constructor
80  MRFZoneList(const fvMesh& mesh, const dictionary& dict);
81 
82  //- Destructor
83  ~MRFZoneList();
84 
85 
86  // Member Functions
87 
88  //- Return active status
89  bool active(const bool warn = false) const;
90 
91  //- Reset the source list
92  void reset(const dictionary& dict);
93 
94  //- Add the frame acceleration
95  void addAcceleration
96  (
97  const volVectorField& U,
98  volVectorField& ddtU
99  ) const;
100 
101  //- Add the frame acceleration contribution to the momentum equation
102  void addAcceleration(fvVectorMatrix& UEqn) const;
103 
104  //- Add the frame acceleration contribution to the momentum equation
105  void addAcceleration
106  (
107  const volScalarField& rho,
109  ) const;
110 
111  //- Return the frame acceleration
113  (
114  const volVectorField& U
115  ) const;
116 
117  //- Return the frame acceleration
119  (
120  const volScalarField& rho,
121  const volVectorField& U
122  ) const;
123 
124  //- Make the given absolute velocity relative within the MRF region
125  void makeRelative(volVectorField& U) const;
126 
127  //- Make the given absolute flux relative within the MRF region
128  void makeRelative(surfaceScalarField& phi) const;
129 
130  //- Return the given absolute flux relative within the MRF region
132  (
134  ) const;
135 
136  //- Return the given absolute boundary flux relative within
137  // the MRF region
139  (
141  ) const;
142 
143  //- Return the given absolute patch flux relative within
144  // the MRF region
146  (
147  const tmp<Field<scalar>>& tphi,
148  const label patchi
149  ) const;
150 
151  //- Make the given absolute mass-flux relative within the MRF region
152  void makeRelative
153  (
154  const surfaceScalarField& rho,
156  ) const;
157 
158  //- Make the given relative velocity absolute within the MRF region
159  void makeAbsolute(volVectorField& U) const;
160 
161  //- Make the given relative flux absolute within the MRF region
162  void makeAbsolute(surfaceScalarField& phi) const;
163 
164  //- Return the given relative flux absolute within the MRF region
166  (
168  ) const;
169 
170  //- Make the given relative mass-flux absolute within the MRF region
171  void makeAbsolute
172  (
173  const surfaceScalarField& rho,
175  ) const;
176 
177  //- Correct the boundary velocity for the rotation of the MRF region
179 
180  //- Correct the boundary flux for the rotation of the MRF region
182  (
183  const volVectorField& U,
185  ) const;
186 
187  //- Filter-out the MRF region contribution from the given field
188  // setting the corresponding values to zero
189  template<class Type>
191  (
193  ) const;
194 
195  //- Update MRFZone faces if the mesh topology changes
196  void update();
197 
198 
199  // I-O
200 
201  //- Read dictionary
202  bool read(const dictionary& dict);
203 
204  //- Write data to Ostream
205  bool writeData(Ostream& os) const;
206 
207  //- Ostream operator
208  friend Ostream& operator<<
209  (
210  Ostream& os,
211  const MRFZoneList& models
212  );
213 };
214 
215 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
216 
217 } // End namespace Foam
218 
219 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
220 
221 #ifdef NoRepository
222  #include "MRFZoneListTemplates.C"
223 #endif
224 
225 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
226 
227 #endif
228 
229 // ************************************************************************* //
~MRFZoneList()
Destructor.
Definition: MRFZoneList.C:49
void makeAbsolute(volVectorField &U) const
Make the given relative velocity absolute within the MRF region.
Definition: MRFZoneList.C:319
dictionary dict
void update()
Update MRFZone faces if the mesh topology changes.
Definition: MRFZoneList.C:417
void reset(const dictionary &dict)
Reset the source list.
Definition: MRFZoneList.C:72
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
surfaceScalarField & phi
void addAcceleration(const volVectorField &U, volVectorField &ddtU) const
Add the frame acceleration.
Definition: MRFZoneList.C:130
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
const fvMesh & mesh_
Reference to the mesh database.
Definition: MRFZoneList.H:73
Generic field type.
Definition: FieldField.H:51
bool active(const bool warn=false) const
Return active status.
Definition: MRFZoneList.C:55
void correctBoundaryVelocity(volVectorField &U) const
Correct the boundary velocity for the rotation of the MRF region.
Definition: MRFZoneList.C:381
dynamicFvMesh & mesh
void correctBoundaryFlux(const volVectorField &U, surfaceScalarField &phi) const
Correct the boundary flux for the rotation of the MRF region.
Definition: MRFZoneList.C:391
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:338
tmp< volVectorField > DDt(const volVectorField &U) const
Return the frame acceleration.
Definition: MRFZoneList.C:165
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
bool writeData(Ostream &os) const
Write data to Ostream.
Definition: MRFZoneList.C:117
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:223
Forward declarations of fvMatrix specializations.
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:63
fvVectorMatrix & UEqn
Definition: UEqn.H:13
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Ostream & operator<<(Ostream &, const ensightPart &)
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:204
bool read(const dictionary &dict)
Read dictionary.
Definition: MRFZoneList.C:104
Namespace for OpenFOAM.