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-2022 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 )
61 :
63  fanCurve_(),
64  direction_(ffdOut)
65 {}
66 
67 
69 (
70  const fvPatch& p,
72  const dictionary& dict
73 )
74 :
76  fanCurve_(Function1<scalar>::New("fanCurve", dict)),
77  direction_(fanFlowDirectionNames_.read(dict.lookup("direction")))
78 {}
79 
80 
82 (
83  const fanPressureFvPatchScalarField& pfopsf,
84  const fvPatch& p,
86  const fvPatchFieldMapper& mapper
87 )
88 :
89  totalPressureFvPatchScalarField(pfopsf, p, iF, mapper),
90  fanCurve_(pfopsf.fanCurve_, false),
91  direction_(pfopsf.direction_)
92 {}
93 
94 
96 (
97  const fanPressureFvPatchScalarField& pfopsf,
99 )
100 :
102  fanCurve_(pfopsf.fanCurve_, false),
103  direction_(pfopsf.direction_)
104 {}
105 
106 
107 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
108 
110 {
111  if (updated())
112  {
113  return;
114  }
115 
116  const fvsPatchField<scalar>& phip =
117  patch().lookupPatchField<surfaceScalarField, scalar>(phiName_);
118 
119  const fvPatchField<vector>& Up =
120  patch().lookupPatchField<volVectorField, vector>(UName_);
121 
122  const int sign = direction_ ? +1 : -1;
123 
124  // Get the volumetric flow rate
125  scalar volFlowRate = 0;
126  if (phip.internalField().dimensions() == dimFlux)
127  {
128  volFlowRate = sign*gSum(phip);
129  }
130  else if
131  (
132  phip.internalField().dimensions() == dimMassFlux
133  )
134  {
135  const scalarField& rhop =
136  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
137 
138  volFlowRate = sign*gSum(phip/rhop);
139  }
140  else
141  {
143  << "dimensions of phi are not correct"
144  << "\n on patch " << patch().name()
145  << " of field " << internalField().name()
146  << " in file " << internalField().objectPath() << nl
147  << exit(FatalError);
148  }
149 
150  // Pressure drop for this flow rate
151  const scalar dp0 = fanCurve_->value(max(volFlowRate, scalar(0)));
152 
154  (
155  p0_ - sign*dp0,
156  -0.5*neg(phip)*magSqr(Up)
157  );
158 }
159 
160 
162 {
164  writeEntry(os, fanCurve_());
165  writeEntry(os, "direction", fanFlowDirectionNames_[direction_]);
166 }
167 
168 
169 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
170 
171 namespace Foam
172 {
174  (
177  );
178 };
179 
180 
181 // ************************************************************************* //
Foam::surfaceFields.
dimensionedScalar sign(const dimensionedScalar &ds)
Run-time selectable general function of one variable.
Definition: Function1.H:52
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
This boundary condition can be applied to assign either a pressure inlet or outlet total pressure con...
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
virtual void write(Ostream &) const
Write.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:59
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:51
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
dimensionedScalar neg(const dimensionedScalar &ds)
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
Type gSum(const FieldField< Field, Type > &f)
Foam::fvPatchFieldMapper.
const dimensionSet dimFlux
void updateCoeffs(const scalarField &p0p, const scalarField &K0mKp)
Update the coefficients associated with the patch field.
dimensioned< scalar > magSqr(const dimensioned< Type > &)
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
static const char nl
Definition: Ostream.H:260
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
const dimensionSet dimMassFlux
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const DimensionedField< Type, surfaceMesh > & internalField() const
Return dimensioned internal field reference.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
static const NamedEnum< fanFlowDirection, 2 > fanFlowDirectionNames_
Fan flow directions names.
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:65
Namespace for OpenFOAM.
fanPressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
virtual void write(Ostream &) const
Write.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:864
Inflow, outflow and entrainment pressure boundary condition based on a constant total pressure assump...