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-2024 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 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32 
34 (
35  const fvPatch& p,
37  const dictionary& dict
38 )
39 :
40  fixedJumpFvPatchScalarField(p, iF, dict, true),
41  phiName_
42  (
43  cyclicPatch().owner()
44  ? dict.lookupOrDefault<word>("phi", "phi")
45  : word::null
46  ),
47  rhoName_
48  (
49  cyclicPatch().owner()
50  ? dict.lookupOrDefault<word>("rho", "rho")
51  : word::null
52  ),
53  fanCurve_(),
54  jumpTable_(),
55  reverse_(dict.lookupOrDefault<Switch>("reverse", false))
56 {
57  if (cyclicPatch().owner())
58  {
59  if (!dict.found("fanCurve") && dict.found("jumpTable"))
60  {
61  // Backwards compatibility fallback
62  jumpTable_ =
64  (
65  "jumpTable",
67  iF.dimensions(),
68  dict
69  );
70  }
71  else
72  {
73  // Preferred entry name
74  fanCurve_ =
76  (
77  "fanCurve",
79  iF.dimensions(),
80  dict
81  );
82  }
83  }
84 }
85 
86 
88 (
90  const fvPatch& p,
92  const fieldMapper& mapper
93 )
94 :
95  fixedJumpFvPatchScalarField(ptf, p, iF, mapper),
96  phiName_(ptf.phiName_),
97  rhoName_(ptf.rhoName_),
98  fanCurve_(ptf.fanCurve_, false),
99  jumpTable_(ptf.jumpTable_, false),
100  reverse_(ptf.reverse_)
101 {}
102 
103 
105 (
108 )
109 :
110  fixedJumpFvPatchScalarField(ptf, iF),
111  phiName_(ptf.phiName_),
112  rhoName_(ptf.rhoName_),
113  fanCurve_(ptf.fanCurve_, false),
114  jumpTable_(ptf.jumpTable_, false),
115  reverse_(ptf.reverse_)
116 {}
117 
118 
119 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
120 
122 {
123  if (updated())
124  {
125  return;
126  }
127 
128  if (cyclicPatch().owner())
129  {
130  const fvsPatchField<scalar>& phip =
131  patch().lookupPatchField<surfaceScalarField, scalar>(phiName_);
132 
133  const scalar sign = reverse_ ? -1 : 1;
134 
135  if (fanCurve_.valid())
136  {
137  // Preferred method
138 
139  scalar volFlowRate = 0;
140 
141  if (phip.internalField().dimensions() == dimVolumetricFlux)
142  {
143  volFlowRate = gSum(phip);
144  }
145  else if
146  (
147  phip.internalField().dimensions() == dimMassFlux
148  )
149  {
150  const scalarField& rhop =
151  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
152 
153  volFlowRate = gSum(phip/rhop);
154  }
155  else
156  {
158  << "dimensions of phi are not correct"
159  << "\n on patch " << patch().name()
160  << " of field " << internalField().name()
161  << " in file " << internalField().objectPath() << nl
162  << exit(FatalError);
163  }
164 
165  jumpRef() = sign*max(fanCurve_->value(max(sign*volFlowRate, 0)), 0);
166  }
167  else
168  {
169  // Backwards compatibility fallback
170 
171  scalarField Un(max(sign*phip/patch().magSf(), scalar(0)));
172 
173  if (phip.internalField().dimensions() == dimVolumetricFlux)
174  {
175  // Do nothing
176  }
177  else if
178  (
179  phip.internalField().dimensions() == dimMassFlux
180  )
181  {
182  const fvPatchField<scalar>& rhop =
183  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
184 
185  Un /= rhop;
186  }
187  else
188  {
190  << "dimensions of phi are not correct"
191  << "\n on patch " << patch().name()
192  << " of field " << internalField().name()
193  << " in file " << internalField().objectPath() << nl
194  << exit(FatalError);
195  }
196 
197  jumpRef() = sign*max(jumpTable_->value(Un), scalar(0));
198  }
199  }
200 
201  fixedJumpFvPatchScalarField::updateCoeffs();
202 }
203 
204 
206 {
208 
209  if (cyclicPatch().owner())
210  {
211  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
212  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
213 
214  if (fanCurve_.valid())
215  {
216  writeEntry(os, fanCurve_());
217  }
218  else
219  {
220  writeEntry(os, jumpTable_());
221  }
222 
223  writeEntryIfDifferent<Switch>(os, "reverse", false, reverse_);
224  }
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...
const dimensionSet & dimensions() const
Return dimensions.
static autoPtr< Function1< Type > > New(const word &name, const Function1s::unitConversions &units, const dictionary &dict)
Select from dictionary.
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:162
This boundary condition provides a fan 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.
Abstract base class for field mapping.
Definition: fieldMapper.H:48
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:88
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:83
const DimensionedField< Type, surfaceMesh > & internalField() const
Return dimensioned internal field reference.
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
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
const dimensionSet dimVolumetricFlux
SurfaceField< scalar > surfaceScalarField
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:64
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimVelocity
error FatalError
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
static const char nl
Definition: Ostream.H:266
dictionary dict
volScalarField & p
Foam::surfaceFields.