All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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-2020 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
98 )
99 :
101  fanCurve_(pfopsf.fanCurve_, false),
102  direction_(pfopsf.direction_)
103 {}
104 
105 
107 (
108  const fanPressureFvPatchScalarField& pfopsf,
110 )
111 :
113  fanCurve_(pfopsf.fanCurve_, false),
114  direction_(pfopsf.direction_)
115 {}
116 
117 
118 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
119 
121 {
122  if (updated())
123  {
124  return;
125  }
126 
127  const fvsPatchField<scalar>& phip =
128  patch().lookupPatchField<surfaceScalarField, scalar>(phiName_);
129  const fvPatchField<vector>& Up =
130  patch().lookupPatchField<volVectorField, vector>(UName_);
131 
132  int sign = direction_ ? +1 : -1;
133 
134  // Get the volumetric flow rate
135  scalar volFlowRate = 0;
136  if (phip.internalField().dimensions() == dimVelocity*dimArea)
137  {
138  volFlowRate = sign*gSum(phip);
139  }
140  else if
141  (
142  phip.internalField().dimensions() == dimVelocity*dimArea*dimDensity
143  )
144  {
145  const scalarField& rhop =
146  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
147 
148  volFlowRate = sign*gSum(phip/rhop);
149  }
150  else
151  {
153  << "dimensions of phi are not correct"
154  << "\n on patch " << patch().name()
155  << " of field " << internalField().name()
156  << " in file " << internalField().objectPath() << nl
157  << exit(FatalError);
158  }
159 
160  // Pressure drop for this flow rate
161  const scalar dp0 = fanCurve_->value(max(volFlowRate, scalar(0)));
162 
164 }
165 
166 
168 {
170  writeEntry(os, fanCurve_());
171  writeEntry(os, "direction", fanFlowDirectionNames_[direction_]);
172 }
173 
174 
175 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
176 
177 namespace Foam
178 {
180  (
183  );
184 };
185 
186 
187 // ************************************************************************* //
Foam::surfaceFields.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
dimensionedScalar sign(const dimensionedScalar &ds)
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: Function1.H:52
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:158
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
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:58
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.
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Type gSum(const FieldField< Field, Type > &f)
Foam::fvPatchFieldMapper.
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 dimDensity
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
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
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:812
const dimensionSet dimVelocity
This boundary condition provides a total pressure condition. Four variants are possible: ...