fanFvPatchFields.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2015 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 
26 #include "fanFvPatchFields.H"
28 #include "volFields.H"
29 #include "surfaceFields.H"
30 #include "Tuple2.H"
31 #include "PolynomialEntry.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
38  (
40  fanFvPatchScalarField
41  );
42 }
43 
44 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
45 
46 template<>
48 {
49  if (this->cyclicPatch().owner())
50  {
51  const surfaceScalarField& phi =
52  db().lookupObject<surfaceScalarField>(phiName_);
53 
54  const fvsPatchField<scalar>& phip =
55  patch().patchField<surfaceScalarField, scalar>(phi);
56 
57  int dir = reversed_ ? -1 : 1;
58 
59  // Average volumetric flow rate
60  scalar volFlowRate = 0;
61 
62  if (phi.dimensions() == dimVelocity*dimArea)
63  {
64  volFlowRate = dir*gSum(phip);
65  }
66  else if (phi.dimensions() == dimVelocity*dimArea*dimDensity)
67  {
68  const scalarField& rhop =
69  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
70  volFlowRate = dir*gSum(phip/rhop);
71  }
72  else
73  {
74  FatalErrorIn("fanFvPatchField<Foam::scalar>::calcFanJump()")
75  << "dimensions of phi are not correct"
76  << "\n on patch " << patch().name()
77  << " of field " << dimensionedInternalField().name()
78  << " in file " << dimensionedInternalField().objectPath() << nl
79  << exit(FatalError);
80  }
81 
82  // Pressure drop for this flow rate
83  this->jump_ = this->jumpTable_->value(max(volFlowRate, scalar(0))) * dir;
84  }
85 }
86 
87 
88 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
89 
90 template<>
92 (
93  const fvPatch& p,
95  const dictionary& dict
96 )
97 :
99  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
100  rhoName_(dict.lookupOrDefault<word>("rho", "none")),
101  reversed_(false)
102 {
103  if (this->cyclicPatch().owner())
104  {
105  if (dict.found("f"))
106  {
107  // Backwards compatibility
108  Istream& is = dict.lookup("f");
109  is.format(IOstream::ASCII);
110  scalarList f(is);
111 
112  label nPows = 0;
113  forAll(f, powI)
114  {
115  if (mag(f[powI]) > VSMALL)
116  {
117  nPows++;
118  }
119  }
120  List<Tuple2<scalar, scalar> > coeffs(nPows);
121  nPows = 0;
122  forAll(f, powI)
123  {
124  if (mag(f[powI]) > VSMALL)
125  {
126  coeffs[nPows++] = Tuple2<scalar, scalar>(f[powI], powI);
127  }
128  }
129 
130  this->jumpTable_.reset
131  (
132  new PolynomialEntry<scalar>("jumpTable", coeffs)
133  );
134  }
135  else
136  {
137  // Generic input constructed from dictionary
138  this->jumpTable_ = DataEntry<scalar>::New("jumpTable", dict);
139  }
140 
141  reversed_.readIfPresent("reversed", dict);
142  }
143 
144  if (dict.found("value"))
145  {
146  fvPatchScalarField::operator=
147  (
148  scalarField("value", dict, p.size())
149  );
150  }
151  else
152  {
153  this->evaluate(Pstream::blocking);
154  }
155 }
156 
157 
158 // ************************************************************************* //
Foam::fvPatchScalarField
fvPatchField< scalar > fvPatchScalarField
Definition: fvPatchFieldsFwd.H:38
volFields.H
Foam::IOstream::format
streamFormat format() const
Return current stream format.
Definition: IOstream.H:377
Foam::PolynomialEntry
PolynomialEntry container data entry for scalars. Items are stored in a list of Tuple2's....
Definition: PolynomialEntry.H:60
p
p
Definition: pEqn.H:62
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::compressible::New
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Definition: turbulentFluidThermoModel.C:36
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Tuple2.H
Foam::dimVelocity
const dimensionSet dimVelocity
Foam::dimDensity
const dimensionSet dimDensity
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
dimensionedInternalField
rDeltaT dimensionedInternalField()
Foam::fanFvPatchField::calcFanJump
void calcFanJump()
Calculate the fan pressure jump.
Definition: fanFvPatchField.C:31
Foam::fvsPatchField
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:65
surfaceFields.H
Foam::surfaceFields.
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:564
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::uniformJumpFvPatchField
This boundary condition provides a jump condition, using the cyclic condition as a base....
Definition: uniformJumpFvPatchField.H:100
Foam::makeTemplatePatchTypeField
makeTemplatePatchTypeField(fvPatchScalarField, fanFvPatchScalarField)
Foam::dimArea
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:57
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
fanFvPatchFields.H
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
PolynomialEntry.H
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
f
labelList f(nPoints)
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:52
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
scalarField
volScalarField scalarField(fieldObject, mesh)
FatalErrorIn
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:313
Foam::Tuple2< scalar, scalar >
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::fanFvPatchField::fanFvPatchField
fanFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
Definition: fanFvPatchField.C:44
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:51