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-2023 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  const dictionary& dict
114 )
115 :
116  fixedJumpFvPatchScalarField(p, iF),
117  fanCurve_(),
118  jumpTable_(),
119  reverse_(dict.lookupOrDefault<Switch>("reverse", false)),
120  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
121  rhoName_(dict.lookupOrDefault<word>("rho", "rho"))
122 {
123  if (cyclicPatch().owner())
124  {
125  if (dict.found("fanCurve"))
126  {
127  // Preferred entry name
128  fanCurve_ = Function1<scalar>::New("fanCurve", dict);
129  }
130  else if (dict.found("jumpTable"))
131  {
132  // Backwards compatibility fallback
133  jumpTable_ = Function1<scalar>::New("jumpTable", dict);
134  }
135  else
136  {
137  // Call function construction for error message
138  fanCurve_ = Function1<scalar>::New("fanCurve", dict);
139  }
140  }
141 
142  if (dict.found("value"))
143  {
145  (
146  scalarField("value", dict, p.size())
147  );
148  }
149  else
150  {
152  }
153 }
154 
155 
157 (
159  const fvPatch& p,
161  const fvPatchFieldMapper& mapper
162 )
163 :
164  fixedJumpFvPatchScalarField(ptf, p, iF, mapper),
165  fanCurve_(ptf.fanCurve_, false),
166  jumpTable_(ptf.jumpTable_, false),
167  reverse_(ptf.reverse_),
168  phiName_(ptf.phiName_),
169  rhoName_(ptf.rhoName_)
170 {}
171 
172 
174 (
177 )
178 :
179  fixedJumpFvPatchScalarField(ptf, iF),
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 
188 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
189 
191 {
192  if (updated())
193  {
194  return;
195  }
196 
197  if (cyclicPatch().owner())
198  {
199  calcFanJump();
200  }
201 
202  fixedJumpFvPatchScalarField::updateCoeffs();
203 }
204 
205 
207 {
209 
210  if (cyclicPatch().owner())
211  {
212  if (fanCurve_.valid())
213  {
214  writeEntry(os, fanCurve_());
215  }
216  else
217  {
218  writeEntry(os, jumpTable_());
219  }
220  }
221 
222  writeEntryIfDifferent<Switch>(os, "reverse", false, reverse_);
223  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
224  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
225 }
226 
227 
228 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
229 
230 namespace Foam
231 {
233  (
236  );
237 };
238 
239 
240 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
static autoPtr< Function1< Type > > New(const word &name, const dictionary &dict)
Selector.
Definition: Function1New.C:32
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:61
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
This boundary condition provides a pressure jump condition, using the cyclic condition as a base.
fanPressureJumpFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Foam::fvPatchFieldMapper.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:87
friend Ostream & operator(Ostream &, const fvPatchField< Type > &)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Type gSum(const FieldField< Field, Type > &f)
dimensionedScalar sign(const dimensionedScalar &ds)
const dimensionSet dimMassFlux
SurfaceField< scalar > surfaceScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
void evaluate(GeometricField< Type, PatchField, GeoMesh > &result, const Function1< Type > &func, const GeometricField< Type, PatchField, GeoMesh > &x)
error FatalError
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
static const char nl
Definition: Ostream.H:260
const dimensionSet dimFlux
dictionary dict
volScalarField & p
Foam::surfaceFields.