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-2021 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 scalar sign = reverse_ ? -1 : 1;
39 
40  if (fanCurve_.valid())
41  {
42  // Preferred method
43 
44  scalar volFlowRate = 0;
45 
46  if (phip.internalField().dimensions() == dimFlux)
47  {
48  volFlowRate = gSum(phip);
49  }
50  else if
51  (
52  phip.internalField().dimensions() == dimMassFlux
53  )
54  {
55  const scalarField& rhop =
56  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
57 
58  volFlowRate = gSum(phip/rhop);
59  }
60  else
61  {
63  << "dimensions of phi are not correct"
64  << "\n on patch " << patch().name()
65  << " of field " << internalField().name()
66  << " in file " << internalField().objectPath() << nl
67  << exit(FatalError);
68  }
69 
70  jump_ = sign*max(fanCurve_->value(max(sign*volFlowRate, 0)), 0);
71  }
72  else
73  {
74  // Backwards compatibility fallback
75 
76  scalarField Un(max(sign*phip/patch().magSf(), scalar(0)));
77 
78  if (phip.internalField().dimensions() == dimFlux)
79  {
80  // Do nothing
81  }
82  else if
83  (
84  phip.internalField().dimensions() == dimMassFlux
85  )
86  {
87  const fvPatchField<scalar>& rhop =
88  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
89 
90  Un /= rhop;
91  }
92  else
93  {
95  << "dimensions of phi are not correct"
96  << "\n on patch " << patch().name()
97  << " of field " << internalField().name()
98  << " in file " << internalField().objectPath() << nl
99  << exit(FatalError);
100  }
101 
102  jump_ = sign*max(jumpTable_->value(Un), scalar(0));
103  }
104 }
105 
106 
107 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
108 
110 (
111  const fvPatch& p,
113 )
114 :
115  fixedJumpFvPatchScalarField(p, iF),
116  fanCurve_(),
117  jumpTable_(),
118  reverse_(false),
119  phiName_("phi"),
120  rhoName_("rho")
121 {}
122 
123 
125 (
126  const fvPatch& p,
128  const dictionary& dict
129 )
130 :
131  fixedJumpFvPatchScalarField(p, iF),
132  fanCurve_(),
133  jumpTable_(),
134  reverse_(dict.lookupOrDefault<Switch>("reverse", false)),
135  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
136  rhoName_(dict.lookupOrDefault<word>("rho", "rho"))
137 {
138  if (cyclicPatch().owner())
139  {
140  if (dict.found("fanCurve"))
141  {
142  // Preferred entry name
143  fanCurve_ = Function1<scalar>::New("fanCurve", dict);
144  }
145  else if (dict.found("jumpTable"))
146  {
147  // Backwards compatibility fallback
148  jumpTable_ = Function1<scalar>::New("jumpTable", dict);
149  }
150  else
151  {
152  // Call function construction for error message
153  fanCurve_ = Function1<scalar>::New("fanCurve", dict);
154  }
155  }
156 
157  if (dict.found("value"))
158  {
160  (
161  scalarField("value", dict, p.size())
162  );
163  }
164  else
165  {
167  }
168 }
169 
170 
172 (
174  const fvPatch& p,
176  const fvPatchFieldMapper& mapper
177 )
178 :
179  fixedJumpFvPatchScalarField(ptf, p, iF, mapper),
180  fanCurve_(ptf.fanCurve_, false),
181  jumpTable_(ptf.jumpTable_, false),
182  reverse_(ptf.reverse_),
183  phiName_(ptf.phiName_),
184  rhoName_(ptf.rhoName_)
185 {}
186 
187 
189 (
192 )
193 :
194  fixedJumpFvPatchScalarField(ptf, iF),
195  fanCurve_(ptf.fanCurve_, false),
196  jumpTable_(ptf.jumpTable_, false),
197  reverse_(ptf.reverse_),
198  phiName_(ptf.phiName_),
199  rhoName_(ptf.rhoName_)
200 {}
201 
202 
203 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
204 
206 {
207  if (updated())
208  {
209  return;
210  }
211 
212  if (cyclicPatch().owner())
213  {
214  calcFanJump();
215  }
216 
217  fixedJumpFvPatchScalarField::updateCoeffs();
218 }
219 
220 
222 {
224 
225  if (cyclicPatch().owner())
226  {
227  if (fanCurve_.valid())
228  {
229  writeEntry(os, fanCurve_());
230  }
231  else
232  {
233  writeEntry(os, jumpTable_());
234  }
235  }
236 
237  writeEntryIfDifferent<Switch>(os, "reverse", false, reverse_);
238  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
239  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
240 }
241 
242 
243 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
244 
245 namespace Foam
246 {
248  (
251  );
252 };
253 
254 
255 // ************************************************************************* //
Foam::surfaceFields.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:643
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
void evaluate(GeometricField< Type, PatchField, GeoMesh > &result, const Function1< Type > &func, const GeometricField< Type, PatchField, GeoMesh > &x)
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:323
friend Ostream & operator(Ostream &, const fvPatchField< Type > &)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:62
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.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Type gSum(const FieldField< Field, Type > &f)
virtual Type value(const scalar x) const =0
Return value as a function of scalar x.
A class for handling words, derived from string.
Definition: word.H:59
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Foam::fvPatchFieldMapper.
const dimensionSet dimFlux
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
bool valid() const
Return true if the autoPtr valid (ie, the pointer is set)
Definition: autoPtrI.H:83
virtual label size() const
Return size.
Definition: fvPatch.H:156
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
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 dimMassFlux
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)
fanPressureJumpFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Namespace for OpenFOAM.
static autoPtr< Function1< Type > > New(const word &name, const dictionary &dict)
Selector.
Definition: Function1New.C:32