fanPressureJumpFvPatchScalarField.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-2019 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 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
32 
33 void Foam::fanPressureJumpFvPatchScalarField::calcFanJump()
34 {
35  const fvsPatchField<scalar>& phip =
36  patch().lookupPatchField<surfaceScalarField, scalar>(phiName_);
37 
38  const bool massBasedFlux =
40 
41  const scalar dir = reverse_ ? -1 : 1;
42 
43  //- Jump table is defined in backward compatibility mode
44  if(jumpTable_.valid())
45  {
46  //- Patch normal velocity field
47  scalarField Un(max(dir*phip/patch().magSf(), scalar(0)));
48 
49  if (massBasedFlux)
50  {
51  const fvPatchField<scalar>& rhop =
52  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
53 
54  Un /= rhop;
55  }
56 
57  jump_ = dir*max(jumpTable_->value(Un), scalar(0));
58  }
59  else
60  {
61  scalar volFlowRate = 0;
62 
63  if (massBasedFlux)
64  {
65  const scalarField& rhop =
66  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
67 
68  volFlowRate = gSum(phip/rhop);
69  }
70  else
71  {
72  volFlowRate = gSum(phip);
73  }
74 
75 
76  jump_ = dir*max(fanCurve_->value(max(dir*volFlowRate, 0)), 0);
77  }
78 }
79 
80 
81 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
82 
84 (
85  const fvPatch& p,
87 )
88 :
89  fixedJumpFvPatchScalarField(p, iF),
90  fanCurve_(),
91  jumpTable_(),
92  reverse_(false),
93  phiName_("phi"),
94  rhoName_("rho")
95 {}
96 
97 
99 (
100  const fvPatch& p,
102  const dictionary& dict
103 )
104 :
105  fixedJumpFvPatchScalarField(p, iF),
106  fanCurve_(),
107  jumpTable_(),
108  reverse_(dict.lookupOrDefault<Switch>("reverse", false)),
109  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
110  rhoName_(dict.lookupOrDefault<word>("rho", "rho"))
111 {
112  if (cyclicPatch().owner())
113  {
114  if (dict.found("jumpTable"))
115  {
116  //- backward compatibility model
117  jumpTable_ = Function1<scalar>::New("jumpTable", dict);
118  }
119  else
120  {
121  fanCurve_ = Function1<scalar>::New("fanCurve", dict);
122  }
123  }
124 
125  if (dict.found("value"))
126  {
128  (
129  scalarField("value", dict, p.size())
130  );
131  }
132  else
133  {
135  }
136 }
137 
138 
140 (
142  const fvPatch& p,
144  const fvPatchFieldMapper& mapper
145 )
146 :
147  fixedJumpFvPatchScalarField(ptf, p, iF, mapper),
148  fanCurve_(ptf.fanCurve_, false),
149  jumpTable_(ptf.jumpTable_, false),
150  reverse_(ptf.reverse_),
151  phiName_(ptf.phiName_),
152  rhoName_(ptf.rhoName_)
153 {}
154 
155 
157 (
159 )
160 :
161  fixedJumpFvPatchScalarField(ptf),
162  fanCurve_(ptf.fanCurve_, false),
163  jumpTable_(ptf.jumpTable_, false),
164  reverse_(ptf.reverse_),
165  phiName_(ptf.phiName_),
166  rhoName_(ptf.rhoName_)
167 {}
168 
169 
171 (
174 )
175 :
176  fixedJumpFvPatchScalarField(ptf, iF),
177  fanCurve_(ptf.fanCurve_, false),
178  jumpTable_(ptf.jumpTable_, false),
179  reverse_(ptf.reverse_),
180  phiName_(ptf.phiName_),
181  rhoName_(ptf.rhoName_)
182 {}
183 
184 
185 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
186 
188 {
189  if (updated())
190  {
191  return;
192  }
193 
194  if (cyclicPatch().owner())
195  {
196  calcFanJump();
197  }
198 
199  fixedJumpFvPatchScalarField::updateCoeffs();
200 }
201 
202 
204 {
206 
207  if (jumpTable_.valid()) writeEntry(os, jumpTable_());
208  if (fanCurve_.valid()) writeEntry(os, fanCurve_());
209 
210  writeEntryIfDifferent<Switch>(os, "reverse", false, reverse_);
211  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
212  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
213 }
214 
215 
216 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
217 
218 namespace Foam
219 {
221  (
224  );
225 };
226 
227 
228 // ************************************************************************* //
Foam::surfaceFields.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:438
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 > &)
friend Ostream & operator(Ostream &, const fvPatchField< Type > &)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, y/n, t/f, or none/any.
Definition: Switch.H:60
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
Macros for easy insertion into run-time selection tables.
const Internal & internalField() const
Return a const-reference to the dimensioned internal field.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
void write(Ostream &, const label, const dictionary &)
Write with dictionary lookup.
const dimensionSet & dimensions() const
Return dimensions.
Type gSum(const FieldField< Field, Type > &f)
A class for handling words, derived from string.
Definition: word.H:59
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Foam::fvPatchFieldMapper.
virtual label size() const
Return size.
Definition: fvPatch.H:155
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
This boundary condition provides a pressure jump condition, using the cyclic condition as a base...
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
const dimensionSet dimDensity
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
static autoPtr< Function1< Type > > New(const word &entryName, const dictionary &dict)
Selector.
Definition: Function1New.C:32
fanPressureJumpFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
Namespace for OpenFOAM.
const dimensionSet dimVelocity