fanPressureFvPatchScalarField.C
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) 2011-2025 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 \*---------------------------------------------------------------------------*/
25 
28 #include "volFields.H"
29 #include "surfaceFields.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 const Foam::NamedEnum
34 <
36  2
38 {
39  "in",
40  "out"
41 };
42 
43 
44 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 
47 (
48  const fvPatch& p,
50  const dictionary& dict
51 )
52 :
54  fanCurve_
55  (
56  Function1<scalar>::New
57  (
58  "fanCurve",
60  iF.dimensions(),
61  dict
62  )
63  ),
64  direction_(fanFlowDirectionNames_.read(dict.lookup("direction")))
65 {}
66 
67 
69 (
70  const fanPressureFvPatchScalarField& pfopsf,
71  const fvPatch& p,
73  const fieldMapper& mapper
74 )
75 :
76  totalPressureFvPatchScalarField(pfopsf, p, iF, mapper),
77  fanCurve_(pfopsf.fanCurve_, false),
78  direction_(pfopsf.direction_)
79 {}
80 
81 
83 (
84  const fanPressureFvPatchScalarField& pfopsf,
86 )
87 :
89  fanCurve_(pfopsf.fanCurve_, false),
90  direction_(pfopsf.direction_)
91 {}
92 
93 
94 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
95 
97 {
98  if (updated())
99  {
100  return;
101  }
102 
103  const fvsPatchField<scalar>& phip =
104  patch().lookupPatchField<surfaceScalarField, scalar>(phiName_);
105 
106  const fvPatchField<vector>& Up =
107  patch().lookupPatchField<volVectorField, vector>(UName_);
108 
109  const int sign = direction_ ? +1 : -1;
110 
111  // Get the volumetric flow rate
112  scalar volFlowRate = 0;
113  if (phip.internalField().dimensions() == dimVolumetricFlux)
114  {
115  volFlowRate = sign*gSum(phip);
116  }
117  else if
118  (
119  phip.internalField().dimensions() == dimMassFlux
120  )
121  {
122  const scalarField& rhop =
123  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
124 
125  volFlowRate = sign*gSum(phip/rhop);
126  }
127  else
128  {
130  << "dimensions of phi are not correct"
131  << "\n on patch " << patch().name()
132  << " of field " << internalField().name()
133  << " in file " << internalField().objectPath() << nl
134  << exit(FatalError);
135  }
136 
137  // Pressure drop for this flow rate
138  const scalar dp0 = fanCurve_->value(max(volFlowRate, scalar(0)));
139 
141  (
142  p0_ - sign*dp0,
143  -0.5*neg(phip)*magSqr(Up)
144  );
145 }
146 
147 
149 {
151  writeEntry(os, fanCurve_());
152  writeEntry(os, "direction", fanFlowDirectionNames_[direction_]);
153 }
154 
155 
156 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
157 
158 namespace Foam
159 {
161  (
164  );
165 };
166 
167 
168 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Run-time selectable general function of one variable.
Definition: Function1.H:125
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:55
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
void updateCoeffs(const scalarField &p0p, const scalarField &K0mKp)
Update the coefficients associated with the patch field.
This boundary condition can be applied to assign either a pressure inlet or outlet total pressure con...
virtual void write(Ostream &) const
Write.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
fanPressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
static const NamedEnum< fanFlowDirection, 2 > fanFlowDirectionNames_
Fan flow directions names.
Abstract base class for field mapping.
Definition: fieldMapper.H:48
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:91
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:86
const DimensionedField< Type, surfaceMesh > & internalField() const
Return dimensioned internal field reference.
Inflow, outflow and entrainment pressure boundary condition based on a constant total pressure assump...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Type gSum(const FieldField< Field, Type > &f)
VolField< vector > volVectorField
Definition: volFieldsFwd.H:63
dimensionedScalar sign(const dimensionedScalar &ds)
const HashTable< dimensionSet > & dimensions()
Get the table of dimension sets.
Definition: dimensionSets.C:96
const dimensionSet dimMassFlux
const dimensionSet dimVolumetricFlux
SurfaceField< scalar > surfaceScalarField
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:62
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar neg(const dimensionedScalar &ds)
error FatalError
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
void magSqr(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
static const char nl
Definition: Ostream.H:267
dictionary dict
volScalarField & p
Foam::surfaceFields.