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-2023 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 namespace Foam
34 {
35  template<>
36  const char* NamedEnum
37  <
39  2
40  >::names[] =
41  {
42  "in",
43  "out"
44  };
45 }
46 
47 const Foam::NamedEnum
48 <
50  2
52 
53 
54 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
55 
57 (
58  const fvPatch& p,
60  const dictionary& dict
61 )
62 :
64  fanCurve_(Function1<scalar>::New("fanCurve", dict)),
65  direction_(fanFlowDirectionNames_.read(dict.lookup("direction")))
66 {}
67 
68 
70 (
71  const fanPressureFvPatchScalarField& pfopsf,
72  const fvPatch& p,
74  const fvPatchFieldMapper& mapper
75 )
76 :
77  totalPressureFvPatchScalarField(pfopsf, p, iF, mapper),
78  fanCurve_(pfopsf.fanCurve_, false),
79  direction_(pfopsf.direction_)
80 {}
81 
82 
84 (
85  const fanPressureFvPatchScalarField& pfopsf,
87 )
88 :
90  fanCurve_(pfopsf.fanCurve_, false),
91  direction_(pfopsf.direction_)
92 {}
93 
94 
95 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
96 
98 {
99  if (updated())
100  {
101  return;
102  }
103 
104  const fvsPatchField<scalar>& phip =
105  patch().lookupPatchField<surfaceScalarField, scalar>(phiName_);
106 
107  const fvPatchField<vector>& Up =
108  patch().lookupPatchField<volVectorField, vector>(UName_);
109 
110  const int sign = direction_ ? +1 : -1;
111 
112  // Get the volumetric flow rate
113  scalar volFlowRate = 0;
114  if (phip.internalField().dimensions() == dimFlux)
115  {
116  volFlowRate = sign*gSum(phip);
117  }
118  else if
119  (
120  phip.internalField().dimensions() == dimMassFlux
121  )
122  {
123  const scalarField& rhop =
124  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
125 
126  volFlowRate = sign*gSum(phip/rhop);
127  }
128  else
129  {
131  << "dimensions of phi are not correct"
132  << "\n on patch " << patch().name()
133  << " of field " << internalField().name()
134  << " in file " << internalField().objectPath() << nl
135  << exit(FatalError);
136  }
137 
138  // Pressure drop for this flow rate
139  const scalar dp0 = fanCurve_->value(max(volFlowRate, scalar(0)));
140 
142  (
143  p0_ - sign*dp0,
144  -0.5*neg(phip)*magSqr(Up)
145  );
146 }
147 
148 
150 {
152  writeEntry(os, fanCurve_());
153  writeEntry(os, "direction", fanFlowDirectionNames_[direction_]);
154 }
155 
156 
157 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
158 
159 namespace Foam
160 {
162  (
165  );
166 };
167 
168 
169 // ************************************************************************* //
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:64
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:54
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
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.
Foam::fvPatchFieldMapper.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:87
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:82
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:306
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:62
dimensionedScalar sign(const dimensionedScalar &ds)
const dimensionSet dimMassFlux
SurfaceField< scalar > surfaceScalarField
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar neg(const dimensionedScalar &ds)
error FatalError
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
static const char nl
Definition: Ostream.H:260
dimensioned< scalar > magSqr(const dimensioned< Type > &)
const dimensionSet dimFlux
dictionary dict
volScalarField & p
Foam::surfaceFields.