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 | Copyright (C) 2012-2015 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  PtrList<MRFZone>
59 {
60 private:
61 
62  // Private Member Functions
63 
64  //- Disallow default bitwise copy construct
65  MRFZoneList(const MRFZoneList&);
66 
67  //- Disallow default bitwise assignment
68  void operator=(const MRFZoneList&);
69 
70 
71 protected:
72 
73  // Protected data
74 
75  //- Reference to the mesh database
76  const fvMesh& mesh_;
77 
78 
79 public:
80 
81  //- Constructor
82  MRFZoneList(const fvMesh& mesh, const dictionary& dict);
83 
84  //- Destructor
85  ~MRFZoneList();
86 
87 
88  // Member Functions
89 
90  //- Return active status
91  bool active(const bool warn = false) const;
92 
93  //- Reset the source list
94  void reset(const dictionary& dict);
95 
96  //- Add the frame acceleration
97  void addAcceleration
98  (
99  const volVectorField& U,
100  volVectorField& ddtU
101  ) const;
102 
103  //- Add the frame acceleration contribution to the momentum equation
104  void addAcceleration(fvVectorMatrix& UEqn) const;
105 
106  //- Add the frame acceleration contribution to the momentum equation
107  void addAcceleration
108  (
109  const volScalarField& rho,
111  ) const;
112 
113  //- Return the frame acceleration
115  (
116  const volVectorField& U
117  ) const;
118 
119  //- Return the frame acceleration
121  (
122  const volScalarField& rho,
123  const volVectorField& U
124  ) const;
125 
126  //- Make the given absolute velocity relative within the MRF region
127  void makeRelative(volVectorField& U) const;
128 
129  //- Make the given absolute flux relative within the MRF region
130  void makeRelative(surfaceScalarField& phi) const;
131 
132  //- Return the given absolute flux relative within the MRF region
134  (
136  ) const;
137 
138  //- Return the given absolute boundary flux relative within
139  // the MRF region
141  (
143  ) const;
144 
145  //- Make the given absolute mass-flux relative within the MRF region
146  void makeRelative
147  (
148  const surfaceScalarField& rho,
150  ) const;
151 
152  //- Make the given relative velocity absolute within the MRF region
153  void makeAbsolute(volVectorField& U) const;
154 
155  //- Make the given relative flux absolute within the MRF region
156  void makeAbsolute(surfaceScalarField& phi) const;
157 
158  //- Return the given relative flux absolute within the MRF region
160  (
162  ) const;
163 
164  //- Make the given relative mass-flux absolute within the MRF region
165  void makeAbsolute
166  (
167  const surfaceScalarField& rho,
169  ) const;
170 
171  //- Correct the boundary velocity for the rotation of the MRF region
173 
174  //- Correct the boundary flux for the rotation of the MRF region
176  (
177  const volVectorField& U,
179  ) const;
180 
181 
182  // I-O
183 
184  //- Read dictionary
185  bool read(const dictionary& dict);
186 
187  //- Write data to Ostream
188  bool writeData(Ostream& os) const;
189 
190  //- Ostream operator
191  friend Ostream& operator<<
192  (
193  Ostream& os,
194  const MRFZoneList& models
195  );
196 };
197 
198 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
199 
200 } // End namespace Foam
201 
202 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
203 
204 #endif
205 
206 // ************************************************************************* //
Foam::MRFZoneList::addAcceleration
void addAcceleration(const volVectorField &U, volVectorField &ddtU) const
Add the frame acceleration.
Definition: MRFZoneList.C:130
Foam::FieldField
Generic field type.
Definition: FieldField.H:51
Foam::MRFZoneList::MRFZoneList
MRFZoneList(const MRFZoneList &)
Disallow default bitwise copy construct.
Foam::MRFZoneList::writeData
bool writeData(Ostream &os) const
Write data to Ostream.
Definition: MRFZoneList.C:117
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::MRFZoneList
List container for MRF zomes.
Definition: MRFZoneList.H:55
Foam::MRFZoneList::reset
void reset(const dictionary &dict)
Reset the source list.
Definition: MRFZoneList.C:72
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
MRFZone.H
fvMatricesFwd.H
Forward declarations of fvMatrix specializations.
Foam::MRFZoneList::correctBoundaryVelocity
void correctBoundaryVelocity(volVectorField &U) const
Correct the boundary velocity for the rotation of the MRF region.
Definition: MRFZoneList.C:305
Foam::MRFZoneList::active
bool active(const bool warn=false) const
Return active status.
Definition: MRFZoneList.C:55
U
U
Definition: pEqn.H:46
Foam::MRFZoneList::makeRelative
void makeRelative(volVectorField &U) const
Make the given absolute velocity relative within the MRF region.
Definition: MRFZoneList.C:204
Foam::MRFZoneList::absolute
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &phi) const
Return the given relative flux absolute within the MRF region.
Definition: MRFZoneList.C:282
Foam::MRFZoneList::DDt
tmp< volVectorField > DDt(const volVectorField &U) const
Return the frame acceleration.
Definition: MRFZoneList.C:165
UEqn
tmp< fvVectorMatrix > UEqn(fvm::div(phi, U)+MRF.DDt(U)+turbulence->divDevReff(U)==fvOptions(U))
Foam::PtrList
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:61
Foam::operator<<
Ostream & operator<<(Ostream &, const edgeMesh &)
Definition: edgeMeshIO.C:130
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::MRFZoneList::mesh_
const fvMesh & mesh_
Reference to the mesh database.
Definition: MRFZoneList.H:75
Foam::MRFZoneList::relative
tmp< surfaceScalarField > relative(const tmp< surfaceScalarField > &phi) const
Return the given absolute flux relative within the MRF region.
Definition: MRFZoneList.C:223
rho
rho
Definition: pEqn.H:3
Foam::MRFZoneList::makeAbsolute
void makeAbsolute(volVectorField &U) const
Make the given relative velocity absolute within the MRF region.
Definition: MRFZoneList.C:263
Foam::MRFZoneList::operator=
void operator=(const MRFZoneList &)
Disallow default bitwise assignment.
dictionary.H
Foam::MRFZoneList::correctBoundaryFlux
void correctBoundaryFlux(const volVectorField &U, surfaceScalarField &phi) const
Correct the boundary flux for the rotation of the MRF region.
Definition: MRFZoneList.C:315
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
Foam::MRFZoneList::read
bool read(const dictionary &dict)
Read dictionary.
Definition: MRFZoneList.C:104
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::MRFZoneList::~MRFZoneList
~MRFZoneList()
Destructor.
Definition: MRFZoneList.C:49
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52