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-------------------------------------------------------------------------------
10License
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"
35#include "syncTools.H"
36
37// * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
38
40Foam::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
67Foam::pressurePIDControlInletVelocityFvPatchVectorField::
68pressurePIDControlInletVelocityFvPatchVectorField
69(
70 const fvPatch& p,
72)
73:
75 upstreamName_(word::null),
76 downstreamName_(word::null),
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
95Foam::pressurePIDControlInletVelocityFvPatchVectorField::
96pressurePIDControlInletVelocityFvPatchVectorField
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
125Foam::pressurePIDControlInletVelocityFvPatchVectorField::
126pressurePIDControlInletVelocityFvPatchVectorField
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
154Foam::pressurePIDControlInletVelocityFvPatchVectorField::
155pressurePIDControlInletVelocityFvPatchVectorField
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
181Foam::pressurePIDControlInletVelocityFvPatchVectorField::
182pressurePIDControlInletVelocityFvPatchVectorField
183(
186)
187:
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
226 (
227 db().lookupObject<surfaceScalarField>(phiName_)
228 );
229
230 // Update the old-time quantities
231 if (timeIndex_ != db().time().timeIndex())
232 {
233 timeIndex_ = db().time().timeIndex();
234 oldQ_ = Q_;
235 oldError_ = error_;
236 oldErrorIntegral_ = errorIntegral_;
237 }
238
239 // Get the density
240 scalar rho = 1;
242 {
243 // do nothing ...
244 }
246 {
247 const fvPatchField<scalar>& rhoField =
248 patch().lookupPatchField<volScalarField, scalar>(rhoName_);
249
250 rho = gSum(rhoField*patch().magSf())/gSum(patch().magSf());
251 }
252 else
253 {
255 << "The dimensions of the field " << phiName_
256 << "are not recognised. The dimensions are " << phi.dimensions()
257 << ". The dimensions should be either " << dimVelocity*dimArea
258 << " for an incompressible case, or "
259 << dimDensity*dimVelocity*dimArea << " for a compressible case."
260 << exit(FatalError);
261 }
262
263 // Patch properties
264 const scalar patchA = gSum(patch().magSf());
265 Q_ = - gSum(*this & patch().Sf());
266
267 // Face-zone properties (a is upstream, b is downstream)
268 scalar Aa, Ab;
269 vector xa, xb;
270 faceZoneAverage(upstreamName_, mesh.Cf(), Aa, xa);
271 faceZoneAverage(downstreamName_, mesh.Cf(), Ab, xb);
272 const scalar L = mag(xa - xb);
273 const scalar LbyALinear = L/(Aa - Ab)*log(Aa/Ab);
274 const scalar LbyAStep = L/2*(1/Aa + 1/Ab);
275 const scalar LbyA = (1 - shapeFactor_)*LbyALinear + shapeFactor_*LbyAStep;
276
277 // Initialise the pressure drop. If the pressure field does not exist, the
278 // pressure drop is assumed to be that specified. This removes the error,
279 // so there is no control and the analytic inlet velocity is applied. This
280 // scenario only ever going to be applicable to potentialFoam.
281 scalar deltaP = deltaP_;
282 if (db().foundObject<volScalarField>(pName_))
283 {
284 scalar pa, pb;
285 faceZoneAverage(upstreamName_, facePressure(), Aa, pa);
286 faceZoneAverage(downstreamName_, facePressure(), Ab, pb);
287 deltaP = pa - pb;
288 }
289 else
290 {
292 << "The pressure field name, 'p' is \"" << pName_ << "\", "
293 << "but a field of that name was not found. The inlet velocity "
294 << "will be set to an analytical value calculated from the "
295 << "specified pressure drop. No PID control will be done and "
296 << "transient effects will be ignored. This behaviour is designed "
297 << "to be appropriate for potentialFoam solutions. If you are "
298 << "getting this warning from another solver, you have probably "
299 << "specified an incorrect pressure name."
300 << endl << endl;
301 }
302
303 // Target and measured flow rates
304 scalar QTarget, QMeasured;
305 const scalar a = (1/sqr(Ab) - 1/sqr(Aa))/(2*rho);
306 if (!mesh.steady() && db().foundObject<volScalarField>(pName_))
307 {
308 const scalar b = LbyA/deltaT;
309 const scalar c = - LbyA/deltaT*oldQ_ /* - deltaP */;
310 QTarget = (- b + sqrt(sqr(b) - 4*a*(c - deltaP_)))/(2*a);
311 QMeasured = (- b + sqrt(sqr(b) - 4*a*(c - deltaP)))/(2*a);
312 }
313 else
314 {
315 QTarget = sqrt(deltaP_/a);
316 QMeasured = sqrt(deltaP/a);
317 }
318
319 // Errors
320 error_ = QTarget - QMeasured;
321 errorIntegral_ = oldErrorIntegral_ + 0.5*(error_ + oldError_);
322 const scalar errorDifferential = oldError_ - error_;
323
324 // Update the field
325 operator==
326 (
327 - patch().nf()
328 *(
329 QTarget
330 + P_*error_
331 + I_*errorIntegral_
332 + D_*errorDifferential
333 )/patchA
334 );
335
336 // Log output
337 if (debug)
338 {
339 const dimensionSet pDimensions(phi.dimensions()*dimVelocity/dimArea);
340 const scalar error = deltaP/deltaP_ - 1;
341 const scalar newQ = - gSum(*this & patch().Sf());
342 Info<< "pressurePIDControlInletVelocityFvPatchField " << patch().name()
343 << endl << " "
344 << dimensionedScalar("U", dimVelocity, newQ/patchA)
345 << endl << " "
346 << dimensionedScalar("deltaP", pDimensions, deltaP)
347 << " (" << mag(error)*100 << "% "
348 << (error < 0 ? "below" : "above") << " the target)" << endl;
349 }
350
352}
353
354
356(
357 Ostream& os
358) const
359{
361
362 os.writeEntry("deltaP", deltaP_);
363 os.writeEntry("upstream", upstreamName_);
364 os.writeEntry("downstream", downstreamName_);
365 os.writeEntry("shapeFactor", shapeFactor_);
366 os.writeEntryIfDifferent<word>("p", "p", pName_);
367 os.writeEntryIfDifferent<word>("rho", "none", rhoName_);
368 os.writeEntry("P", P_);
369 os.writeEntry("I", I_);
370 os.writeEntry("D", D_);
371 os.writeEntry("error", error_);
372 os.writeEntry("errorIntegral", errorIntegral_);
373
374 writeEntry("value", os);
375}
376
377
378// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
379
380namespace Foam
381{
383 (
386 );
387}
388
389
390// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
surfaceScalarField & phi
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet & dimensions() const
Return dimensions.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:239
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:251
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Definition: boundaryMesh.H:63
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:109
Class to handle errors and exceptions in a simple, consistent stream-based manner.
Definition: error.H:77
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
virtual bool write()
Write the output fields.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
const surfaceVectorField & Cf() const
Return face centres as surfaceVectorField.
A FieldMapper for finite-volume patch fields.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:82
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
This boundary condition tries to generate an inlet velocity that maintains a specified pressure drop ...
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Lookup type of boundary radiation properties.
Definition: lookup.H:66
bool steady() const noexcept
True if default ddtScheme is steady-state.
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Definition: fvPatchField.H:676
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Type gSum(const FieldField< Field, Type > &f)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:59
const dimensionSet dimVelocity
messageStream Info
Information stream (stdout output on master, null elsewhere)
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
dimensionedScalar log(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > linearInterpolate(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: linear.H:112
error FatalError
const dimensionSet dimDensity
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
label timeIndex
Definition: getTimeIndex.H:30
dictionary dict
volScalarField & b
Definition: createFields.H:27
Foam::surfaceFields.
const vector L(dict.get< vector >("L"))