fanMomentumSource.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) 2022 Louis Vittoz, SimScale GmbH
9  Copyright (C) 2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 
27 Class
28  Foam::fv::fanMomentumSource
29 
30 Group
31  grpFvOptionsSources
32 
33 Description
34  This source term models the action of a fan on the flow. It calculates
35  flow rate through a set of upstream faces of encompassing the cell zone. The
36  set of upstream faces is automatically calculated based on the flow
37  direction and the surrounding face zone. Based on the flow rate, the
38  pressure gradient is calculated based on the fan pressure curve and the
39  thickness of the fan (measurement section). The pressure gradient is then
40  added to the momentum equation.
41 
42 Usage
43  Minimal example by using \c constant/fvOptions:
44  \verbatim
45  <name>
46  {
47  // Mandatory entries
48  type fanMomentumSource;
49  fanCurve <Function1<scalar>>;
50  flowDir <vector>;
51  thickness <scalar>;
52  cellZone <word>;
53  faceZone <word>;
54 
55  // Optional entries
56  gradient <scalar>;
57  rho <scalar>;
58  U <word>;
59 
60  // Inherited entries
61  ...
62  }
63  \endverbatim
64 
65  where the entries mean:
66  \table
67  Property | Description | Type | Reqd | Deflt
68  type | Type name: fanMomentumSource | word | yes | -
69  fanCurve | Pressure drop vs flow-rate | Function1<scalar>| yes | -
70  flowDir | Direction of the flow through the fan | vector | yes | -
71  cellZone | Cell zone representing the fan | word | yes | -
72  faceZone | Face zone around the cell zone | word | yes | -
73  thickness | Thickness of the fan [m] | scalar | yes | -
74  gradient | Initial pressure gradient | scalar | no | -
75  rho | Reference density for incompressible flow | scalar | no | -
76  U | Name of velocity field | word | no | U
77  \endtable
78 
79  The inherited entries are elaborated in:
80  - \link fvOption.H \endlink
81  - \link cellSetOption.H \endlink
82  - \link Function1.H \endlink
83 
84 Note
85  - The fan curve needs to be provded for a pressure drop in [Pa] and
86  is specified as a function of the volumetric flow rate in [m^3/s].
87 
88 SourceFiles
89  fanMomentumSource.C
90 
91 \*---------------------------------------------------------------------------*/
92 
93 #ifndef Foam_fv_fanMomentumSource_H
94 #define Foam_fv_fanMomentumSource_H
95 
96 #include "cellSetOption.H"
97 
98 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
99 
100 namespace Foam
101 {
102 namespace fv
103 {
104 
105 /*---------------------------------------------------------------------------*\
106  Class fanMomentumSource Declaration
107 \*---------------------------------------------------------------------------*/
108 
109 class fanMomentumSource
110 :
111  public fv::cellSetOption
112 {
113  // Private Data
114 
115  //- Local list of pairs of upstream patch IDs and upstream face IDs
116  // For internal faces, the upstream patch ID is -1
117  List<labelPair> upstreamPatchFaceInfo_;
118 
119  //- Cells in zone
120  labelHashSet cellsInZones_;
121 
122  //- Volumetric flow rate [m3] vs. pressure drop [Pa] table
123  autoPtr<Function1<scalar>> fanCurvePtr_;
124 
125  //- Name of the velocity field this function object operates on
126  word UName_;
127 
128  //- Direction of flow through the fan
129  vector flowDir_;
130 
131  //- Thickness of the fan, in [m]
132  scalar thickness_;
133 
134  //- Pressure drop
135  scalar gradPFan_;
136 
137  //- Id for the surrounding face zone
138  label surroundingFaceZoneID_;
139 
140  //- Reference density for incompressible cases
141  scalar rho_;
142 
143  //- Whether to use reference density
144  bool useRefRho_;
145 
146 
147  // Private Member Functions
148 
149  //- Write the pressure gradient to file (for restarts etc)
150  void writeProps(const scalar gradP, const scalar flowRate) const;
151 
152  //- Collect all faces upstream of
153  //- the centre of gravity of the cell zone
154  void initUpstreamFaces();
155 
156  //- Return the volumetric flow rate given flux and face-interpolated
157  // density. For incompressible cases, use surfaceScalarField::null()
158  // for rhof
159  scalar calcFlowRate
160  (
161  const surfaceScalarField& phi,
162  const surfaceScalarField& rhof
163  ) const;
165 
166 public:
167 
168  //- Runtime type information
169  TypeName("fanMomentumSource");
170 
171 
172  // Constructors
173 
174  //- Construct from explicit source name and mesh
176  (
177  const word& sourceName,
178  const word& modelType,
179  const dictionary& dict,
180  const fvMesh& mesh
181  );
182 
183  //- No copy construct
184  fanMomentumSource(const fanMomentumSource&) = delete;
185 
186  //- No copy assignment
187  void operator=(const fanMomentumSource&) = delete;
188 
189 
190  //- Destructor
191  virtual ~fanMomentumSource() = default;
192 
193 
194  // Member Functions
195 
196  // Evaluation
197 
198  //- Add explicit contribution to momentum equation
199  virtual void addSup
200  (
201  fvMatrix<vector>& eqn,
202  const label fieldi
203  );
204 
205  //- Add explicit contribution to compressible momentum equation
206  virtual void addSup
207  (
208  const volScalarField& rho,
209  fvMatrix<vector>& eqn,
210  const label fieldi
211  );
212 
213 
214  // IO
215 
216  //- Read source dictionary
217  virtual bool read(const dictionary& dict);
218 };
219 
220 
221 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
222 
223 } // End namespace fv
224 } // End namespace Foam
225 
226 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
227 
228 #endif
229 
230 // ************************************************************************* //
dictionary dict
const fvMesh & mesh() const noexcept
Return const access to the mesh database.
Definition: fvOptionI.H:30
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
void operator=(const fanMomentumSource &)=delete
No copy assignment.
volVectorField gradP(fvc::grad(p))
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
fanMomentumSource(const word &sourceName, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from explicit source name and mesh.
A class for handling words, derived from Foam::string.
Definition: word.H:63
labelList fv(nPoints)
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:64
Vector< scalar > vector
Definition: vector.H:57
virtual ~fanMomentumSource()=default
Destructor.
TypeName("fanMomentumSource")
Runtime type information.
virtual void addSup(fvMatrix< vector > &eqn, const label fieldi)
Add explicit contribution to momentum equation.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
This source term models the action of a fan on the flow. It calculates flow rate through a set of ups...
surfaceScalarField rhof(fvc::interpolate(rho, "div(phi,rho)"))
Namespace for OpenFOAM.
virtual bool read(const dictionary &dict)
Read source dictionary.