pressurePIDControlInletVelocityFvPatchVectorField.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2015-2020 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
29 #include "volFields.H"
31 #include "fvPatchFieldMapper.H"
32 #include "surfaceFields.H"
33 #include "linear.H"
34 #include "steadyStateDdtScheme.H"
35 #include "syncTools.H"
36 
37 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
38 
40 Foam::pressurePIDControlInletVelocityFvPatchVectorField::facePressure() const
41 {
42  const word pfName(pName_ + "f");
43 
44  const volScalarField& p = db().lookupObject<volScalarField>(pName_);
45 
46  surfaceScalarField* pfPtr = db().getObjectPtr<surfaceScalarField>(pfName);
47 
48  if (!pfPtr)
49  {
50  pfPtr = new surfaceScalarField(pfName, linearInterpolate(p));
51  pfPtr->store();
52  }
53 
54  surfaceScalarField& pf = *pfPtr;
55 
56  if (!pf.upToDate(p))
57  {
58  pf = linearInterpolate(p);
59  }
60 
61  return pf;
62 }
63 
64 
65 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
66 
69 (
70  const fvPatch& p,
72 )
73 :
75  upstreamName_(),
76  downstreamName_(),
77  deltaP_(1),
78  shapeFactor_(0),
79  pName_("p"),
80  phiName_("phi"),
81  rhoName_("none"),
82  P_(0),
83  I_(0),
84  D_(0),
85  Q_(- gSum(*this & patch().Sf())),
86  error_(0),
87  errorIntegral_(0),
88  oldQ_(0),
89  oldError_(0),
90  oldErrorIntegral_(0),
91  timeIndex_(db().time().timeIndex())
92 {}
93 
94 
97 (
99  const fvPatch& p,
101  const fvPatchFieldMapper& mapper
102 )
103 :
104  fixedValueFvPatchField<vector>(ptf, p, iF, mapper),
105  upstreamName_(ptf.upstreamName_),
106  downstreamName_(ptf.downstreamName_),
107  deltaP_(ptf.deltaP_),
108  shapeFactor_(ptf.shapeFactor_),
109  pName_(ptf.pName_),
110  phiName_(ptf.phiName_),
111  rhoName_(ptf.rhoName_),
112  P_(ptf.P_),
113  I_(ptf.I_),
114  D_(ptf.D_),
115  Q_(ptf.Q_),
116  error_(ptf.error_),
117  errorIntegral_(ptf.errorIntegral_),
118  oldQ_(ptf.oldQ_),
119  oldError_(ptf.oldError_),
120  oldErrorIntegral_(ptf.oldErrorIntegral_),
121  timeIndex_(ptf.timeIndex_)
122 {}
123 
124 
127 (
128  const fvPatch& p,
130  const dictionary& dict
131 )
132 :
134  upstreamName_(dict.lookup("upstream")),
135  downstreamName_(dict.lookup("downstream")),
136  deltaP_(dict.get<scalar>("deltaP")),
137  shapeFactor_(dict.getOrDefault<scalar>("shapeFactor", 0)),
138  pName_(dict.getOrDefault<word>("p", "p")),
139  phiName_(dict.getOrDefault<word>("phi", "phi")),
140  rhoName_(dict.getOrDefault<word>("rho", "none")),
141  P_(dict.get<scalar>("P")),
142  I_(dict.get<scalar>("I")),
143  D_(dict.get<scalar>("D")),
144  Q_(- gSum(*this & patch().Sf())),
145  error_(dict.getOrDefault<scalar>("error", 0)),
146  errorIntegral_(dict.getOrDefault<scalar>("errorIntegral", 0)),
147  oldQ_(0),
148  oldError_(0),
149  oldErrorIntegral_(0),
150  timeIndex_(db().time().timeIndex())
151 {}
152 
153 
156 (
158 )
159 :
161  upstreamName_(ptf.upstreamName_),
162  downstreamName_(ptf.downstreamName_),
163  deltaP_(ptf.deltaP_),
164  shapeFactor_(ptf.shapeFactor_),
165  pName_(ptf.pName_),
166  phiName_(ptf.phiName_),
167  rhoName_(ptf.rhoName_),
168  P_(ptf.P_),
169  I_(ptf.I_),
170  D_(ptf.D_),
171  Q_(ptf.Q_),
172  error_(ptf.error_),
173  errorIntegral_(ptf.errorIntegral_),
174  oldQ_(ptf.oldQ_),
175  oldError_(ptf.oldError_),
176  oldErrorIntegral_(ptf.oldErrorIntegral_),
177  timeIndex_(ptf.timeIndex_)
178 {}
179 
180 
183 (
186 )
187 :
188  fixedValueFvPatchField<vector>(ptf, iF),
189  upstreamName_(ptf.upstreamName_),
190  downstreamName_(ptf.downstreamName_),
191  deltaP_(ptf.deltaP_),
192  shapeFactor_(ptf.shapeFactor_),
193  pName_(ptf.pName_),
194  phiName_(ptf.phiName_),
195  rhoName_(ptf.rhoName_),
196  P_(ptf.P_),
197  I_(ptf.I_),
198  D_(ptf.D_),
199  Q_(ptf.Q_),
200  error_(ptf.error_),
201  errorIntegral_(ptf.errorIntegral_),
202  oldQ_(ptf.oldQ_),
203  oldError_(ptf.oldError_),
204  oldErrorIntegral_(ptf.oldErrorIntegral_),
205  timeIndex_(ptf.timeIndex_)
206 {}
207 
208 
209 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
210 
212 {
213  if (updated())
214  {
215  return;
216  }
217 
218  // Get the mesh
219  const fvMesh& mesh(patch().boundaryMesh().mesh());
220 
221  // Get the time step
222  const scalar deltaT(db().time().deltaTValue());
223 
224  // Get the flux field
225  const auto& phi = db().lookupObject<surfaceScalarField>(phiName_);
226 
227  // Update the old-time quantities
228  if (timeIndex_ != db().time().timeIndex())
229  {
230  timeIndex_ = db().time().timeIndex();
231  oldQ_ = Q_;
232  oldError_ = error_;
233  oldErrorIntegral_ = errorIntegral_;
234  }
235 
236  // Get the density
237  scalar rho = 1;
238  if (phi.dimensions() == dimVolume/dimTime)
239  {
240  // do nothing ...
241  }
242  else if (phi.dimensions() == dimMass/dimTime)
243  {
244  const auto& rhop = patch().lookupPatchField<volScalarField>(rhoName_);
245 
246  rho = gSum(rhop*patch().magSf())/gSum(patch().magSf());
247  }
248  else
249  {
251  << "The dimensions of the field " << phiName_
252  << "are not recognised. The dimensions are " << phi.dimensions()
253  << ". The dimensions should be either " << dimVolume/dimTime
254  << " for an incompressible case, or "
255  << dimMass/dimTime << " for a compressible case."
256  << exit(FatalError);
257  }
258 
259  // Patch properties
260  const scalar patchA = gSum(patch().magSf());
261  Q_ = - gSum(*this & patch().Sf());
262 
263  // Face-zone properties (a is upstream, b is downstream)
264  scalar Aa, Ab;
265  vector xa, xb;
266  faceZoneAverage(upstreamName_, mesh.Cf(), Aa, xa);
267  faceZoneAverage(downstreamName_, mesh.Cf(), Ab, xb);
268  const scalar L = mag(xa - xb);
269  const scalar LbyALinear = L/(Aa - Ab)*log(Aa/Ab);
270  const scalar LbyAStep = L/2*(1/Aa + 1/Ab);
271  const scalar LbyA = (1 - shapeFactor_)*LbyALinear + shapeFactor_*LbyAStep;
272 
273  // Initialise the pressure drop. If the pressure field does not exist, the
274  // pressure drop is assumed to be that specified. This removes the error,
275  // so there is no control and the analytic inlet velocity is applied. This
276  // scenario only ever going to be applicable to potentialFoam.
277  scalar deltaP = deltaP_;
278  if (db().foundObject<volScalarField>(pName_))
279  {
280  scalar pa, pb;
281  faceZoneAverage(upstreamName_, facePressure(), Aa, pa);
282  faceZoneAverage(downstreamName_, facePressure(), Ab, pb);
283  deltaP = pa - pb;
284  }
285  else
286  {
288  << "The pressure field name, 'p' is \"" << pName_ << "\", "
289  << "but a field of that name was not found. The inlet velocity "
290  << "will be set to an analytical value calculated from the "
291  << "specified pressure drop. No PID control will be done and "
292  << "transient effects will be ignored. This behaviour is designed "
293  << "to be appropriate for potentialFoam solutions. If you are "
294  << "getting this warning from another solver, you have probably "
295  << "specified an incorrect pressure name."
296  << endl << endl;
297  }
298 
299  // Target and measured flow rates
300  scalar QTarget, QMeasured;
301  const scalar a = (1/sqr(Ab) - 1/sqr(Aa))/(2*rho);
302  if (!mesh.steady() && db().foundObject<volScalarField>(pName_))
303  {
304  const scalar b = LbyA/deltaT;
305  const scalar c = - LbyA/deltaT*oldQ_ /* - deltaP */;
306  QTarget = (- b + sqrt(sqr(b) - 4*a*(c - deltaP_)))/(2*a);
307  QMeasured = (- b + sqrt(sqr(b) - 4*a*(c - deltaP)))/(2*a);
308  }
309  else
310  {
311  QTarget = sqrt(deltaP_/a);
312  QMeasured = sqrt(deltaP/a);
313  }
314 
315  // Errors
316  error_ = QTarget - QMeasured;
317  errorIntegral_ = oldErrorIntegral_ + 0.5*(error_ + oldError_);
318  const scalar errorDifferential = oldError_ - error_;
319 
320  // Update the field
321  operator==
322  (
323  - patch().nf()
324  *(
325  QTarget
326  + P_*error_
327  + I_*errorIntegral_
328  + D_*errorDifferential
329  )/patchA
330  );
331 
332  // Log output
333  if (debug)
334  {
335  const dimensionSet pDimensions(phi.dimensions()*dimVelocity/dimArea);
336  const scalar error = deltaP/deltaP_ - 1;
337  const scalar newQ = - gSum(*this & patch().Sf());
338  Info<< "pressurePIDControlInletVelocityFvPatchField " << patch().name()
339  << endl << " "
340  << dimensionedScalar("U", dimVelocity, newQ/patchA)
341  << endl << " "
342  << dimensionedScalar("deltaP", pDimensions, deltaP)
343  << " (" << mag(error)*100 << "% "
344  << (error < 0 ? "below" : "above") << " the target)" << endl;
345  }
346 
348 }
349 
350 
352 (
353  Ostream& os
354 ) const
355 {
357 
358  os.writeEntry("deltaP", deltaP_);
359  os.writeEntry("upstream", upstreamName_);
360  os.writeEntry("downstream", downstreamName_);
361  os.writeEntry("shapeFactor", shapeFactor_);
362  os.writeEntryIfDifferent<word>("p", "p", pName_);
363  os.writeEntryIfDifferent<word>("rho", "none", rhoName_);
364  os.writeEntry("P", P_);
365  os.writeEntry("I", I_);
366  os.writeEntry("D", D_);
367  os.writeEntry("error", error_);
368  os.writeEntry("errorIntegral", errorIntegral_);
369 
371 }
372 
373 
374 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
375 
376 namespace Foam
377 {
379  (
381  pressurePIDControlInletVelocityFvPatchVectorField
382  );
383 }
384 
385 
386 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
Foam::surfaceFields.
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
dictionary dict
fvPatchField< vector > fvPatchVectorField
dimensionedScalar log(const dimensionedScalar &ds)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
const vector L(dict.get< vector >("L"))
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
const surfaceVectorField & Cf() const
Return face centres as surfaceVectorField.
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > linearInterpolate(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: linear.H:120
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:321
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
Definition: fvPatchField.H:375
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:372
Lookup type of boundary radiation properties.
Definition: lookup.H:57
pressurePIDControlInletVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
dynamicFvMesh & mesh
Type gSum(const FieldField< Field, Type > &f)
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
A class for handling words, derived from Foam::string.
Definition: word.H:63
A FieldMapper for finite-volume patch fields.
label timeIndex() const noexcept
Return the current time index.
Definition: TimeStateI.H:43
Vector< scalar > vector
Definition: vector.H:57
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
Definition: Ostream.H:336
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition: IOobject.C:456
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
int debug
Static debugging option.
OBJstream os(runTime.globalPath()/outputName)
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
#define WarningInFunction
Report a warning using Foam::Warning.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:309
const dimensionedScalar c
Speed of light in a vacuum.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
const std::string patch
OpenFOAM patch number as a std::string.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
messageStream Info
Information stream (stdout output on master, null elsewhere)
bool steady() const noexcept
True if default ddt scheme is steady-state.
This boundary condition tries to generate an inlet velocity that maintains a specified pressure drop ...
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
volScalarField & p
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
Namespace for OpenFOAM.
label timeIndex
Definition: getTimeIndex.H:24
const dimensionSet dimVelocity