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 | Copyright (C) 2015 OpenCFD Ltd.
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 
27 #include "volFields.H"
29 #include "fvPatchFieldMapper.H"
30 #include "surfaceFields.H"
31 #include "linear.H"
32 #include "steadyStateDdtScheme.H"
33 #include "syncTools.H"
34 
35 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
36 
38 Foam::pressurePIDControlInletVelocityFvPatchVectorField::facePressure() const
39 {
40  const volScalarField& p(db().lookupObject<volScalarField>(pName_));
41 
42  const word pfName(pName_ + "f");
43 
44  if (!db().foundObject<surfaceScalarField>(pfName))
45  {
46  surfaceScalarField* pfPtr
47  (
49  );
50 
51  pfPtr->store();
52  }
53 
55  (
56  const_cast<surfaceScalarField&>
57  (
58  db().lookupObject<surfaceScalarField>(pfName)
59  )
60  );
61 
62  if (!pf.upToDate(p))
63  {
64  pf = linearInterpolate(p);
65  }
66 
67  return pf;
68 }
69 
70 
71 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
72 
73 Foam::pressurePIDControlInletVelocityFvPatchVectorField::
74 pressurePIDControlInletVelocityFvPatchVectorField
75 (
76  const fvPatch& p,
77  const DimensionedField<vector, volMesh>& iF
78 )
79 :
80  fixedValueFvPatchField<vector>(p, iF),
81  upstreamName_(word::null),
82  downstreamName_(word::null),
83  deltaP_(1),
84  shapeFactor_(0),
85  pName_("p"),
86  phiName_("phi"),
87  rhoName_("none"),
88  P_(0),
89  I_(0),
90  D_(0),
91  Q_(- gSum(*this & patch().Sf())),
92  error_(0),
93  errorIntegral_(0),
94  oldQ_(0),
95  oldError_(0),
96  oldErrorIntegral_(0),
97  timeIndex_(db().time().timeIndex())
98 {}
99 
100 
101 Foam::pressurePIDControlInletVelocityFvPatchVectorField::
102 pressurePIDControlInletVelocityFvPatchVectorField
103 (
104  const pressurePIDControlInletVelocityFvPatchVectorField& ptf,
105  const fvPatch& p,
106  const DimensionedField<vector, volMesh>& iF,
107  const fvPatchFieldMapper& mapper
108 )
109 :
110  fixedValueFvPatchField<vector>(ptf, p, iF, mapper),
111  upstreamName_(ptf.upstreamName_),
112  downstreamName_(ptf.downstreamName_),
113  deltaP_(ptf.deltaP_),
114  shapeFactor_(ptf.shapeFactor_),
115  pName_(ptf.pName_),
116  phiName_(ptf.phiName_),
117  rhoName_(ptf.rhoName_),
118  P_(ptf.P_),
119  I_(ptf.I_),
120  D_(ptf.D_),
121  Q_(ptf.Q_),
122  error_(ptf.error_),
123  errorIntegral_(ptf.errorIntegral_),
124  oldQ_(ptf.oldQ_),
125  oldError_(ptf.oldError_),
126  oldErrorIntegral_(ptf.oldErrorIntegral_),
127  timeIndex_(ptf.timeIndex_)
128 {}
129 
130 
131 Foam::pressurePIDControlInletVelocityFvPatchVectorField::
132 pressurePIDControlInletVelocityFvPatchVectorField
133 (
134  const fvPatch& p,
135  const DimensionedField<vector, volMesh>& iF,
136  const dictionary& dict
137 )
138 :
139  fixedValueFvPatchField<vector>(p, iF, dict),
140  upstreamName_(dict.lookup("upstream")),
141  downstreamName_(dict.lookup("downstream")),
142  deltaP_(readScalar(dict.lookup("deltaP"))),
143  shapeFactor_(dict.lookupOrDefault<scalar>("shapeFactor", 0)),
144  pName_(dict.lookupOrDefault<word>("p", "p")),
145  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
146  rhoName_(dict.lookupOrDefault<word>("rho", "none")),
147  P_(readScalar(dict.lookup("P"))),
148  I_(readScalar(dict.lookup("I"))),
149  D_(readScalar(dict.lookup("D"))),
150  Q_(- gSum(*this & patch().Sf())),
151  error_(dict.lookupOrDefault<scalar>("error", 0)),
152  errorIntegral_(dict.lookupOrDefault<scalar>("errorIntegral", 0)),
153  oldQ_(0),
154  oldError_(0),
155  oldErrorIntegral_(0),
156  timeIndex_(db().time().timeIndex())
157 {}
158 
159 
160 Foam::pressurePIDControlInletVelocityFvPatchVectorField::
161 pressurePIDControlInletVelocityFvPatchVectorField
162 (
163  const pressurePIDControlInletVelocityFvPatchVectorField& ptf
164 )
165 :
166  fixedValueFvPatchField<vector>(ptf),
167  upstreamName_(ptf.upstreamName_),
168  downstreamName_(ptf.downstreamName_),
169  deltaP_(ptf.deltaP_),
170  shapeFactor_(ptf.shapeFactor_),
171  pName_(ptf.pName_),
172  phiName_(ptf.phiName_),
173  rhoName_(ptf.rhoName_),
174  P_(ptf.P_),
175  I_(ptf.I_),
176  D_(ptf.D_),
177  Q_(ptf.Q_),
178  error_(ptf.error_),
179  errorIntegral_(ptf.errorIntegral_),
180  oldQ_(ptf.oldQ_),
181  oldError_(ptf.oldError_),
182  oldErrorIntegral_(ptf.oldErrorIntegral_),
183  timeIndex_(ptf.timeIndex_)
184 {}
185 
186 
187 Foam::pressurePIDControlInletVelocityFvPatchVectorField::
188 pressurePIDControlInletVelocityFvPatchVectorField
189 (
190  const pressurePIDControlInletVelocityFvPatchVectorField& ptf,
191  const DimensionedField<vector, volMesh>& iF
192 )
193 :
194  fixedValueFvPatchField<vector>(ptf, iF),
195  upstreamName_(ptf.upstreamName_),
196  downstreamName_(ptf.downstreamName_),
197  deltaP_(ptf.deltaP_),
198  shapeFactor_(ptf.shapeFactor_),
199  pName_(ptf.pName_),
200  phiName_(ptf.phiName_),
201  rhoName_(ptf.rhoName_),
202  P_(ptf.P_),
203  I_(ptf.I_),
204  D_(ptf.D_),
205  Q_(ptf.Q_),
206  error_(ptf.error_),
207  errorIntegral_(ptf.errorIntegral_),
208  oldQ_(ptf.oldQ_),
209  oldError_(ptf.oldError_),
210  oldErrorIntegral_(ptf.oldErrorIntegral_),
211  timeIndex_(ptf.timeIndex_)
212 {}
213 
214 
215 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
216 
217 void Foam::pressurePIDControlInletVelocityFvPatchVectorField::updateCoeffs()
218 {
219  if (updated())
220  {
221  return;
222  }
223 
224  // Get the mesh
225  const fvMesh& mesh(patch().boundaryMesh().mesh());
226 
227  // Get the time step
228  const scalar deltaT(db().time().deltaTValue());
229 
230  // Get the flux field
231  const surfaceScalarField& phi
232  (
233  db().lookupObject<surfaceScalarField>(phiName_)
234  );
235 
236  // Update the old-time quantities
237  if (timeIndex_ != db().time().timeIndex())
238  {
239  timeIndex_ = db().time().timeIndex();
240  oldQ_ = Q_;
241  oldError_ = error_;
242  oldErrorIntegral_ = errorIntegral_;
243  }
244 
245  // Get the density
246  scalar rho = 1;
247  if (phi.dimensions() == dimVelocity*dimArea)
248  {
249  // do nothing ...
250  }
251  else if (phi.dimensions() == dimDensity*dimVelocity*dimArea)
252  {
253  const fvPatchField<scalar>& rhoField =
254  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
255 
256  rho = gSum(rhoField*patch().magSf())/gSum(patch().magSf());
257  }
258  else
259  {
261  << "The dimensions of the field " << phiName_
262  << "are not recognised. The dimensions are " << phi.dimensions()
263  << ". The dimensions should be either " << dimVelocity*dimArea
264  << " for an incompressible case, or "
265  << dimDensity*dimVelocity*dimArea << " for a compressible case."
266  << exit(FatalError);
267  }
268 
269  // Patch properties
270  const scalar patchA = gSum(patch().magSf());
271  Q_ = - gSum(*this & patch().Sf());
272 
273  // Face-zone properties (a is upstream, b is downstream)
274  scalar Aa, Ab;
275  vector xa, xb;
276  faceZoneAverage(upstreamName_, mesh.Cf(), Aa, xa);
277  faceZoneAverage(downstreamName_, mesh.Cf(), Ab, xb);
278  const scalar L = mag(xa - xb);
279  const scalar LbyALinear = L/(Aa - Ab)*log(Aa/Ab);
280  const scalar LbyAStep = L/2*(1/Aa + 1/Ab);
281  const scalar LbyA = (1 - shapeFactor_)*LbyALinear + shapeFactor_*LbyAStep;
282 
283  // Initialise the pressure drop. If the pressure field does not exist, the
284  // pressure drop is assumed to be that specified. This removes the error,
285  // so there is no control and the analytic inlet velocity is applied. This
286  // scenario only ever going to be applicable to potentialFoam.
287  scalar deltaP = deltaP_;
288  if (db().foundObject<volScalarField>(pName_))
289  {
290  scalar pa, pb;
291  faceZoneAverage(upstreamName_, facePressure(), Aa, pa);
292  faceZoneAverage(downstreamName_, facePressure(), Ab, pb);
293  deltaP = pa - pb;
294  }
295  else
296  {
298  << "The pressure field name, \"pName\", is \"" << pName_ << "\", "
299  << "but a field of that name was not found. The inlet velocity "
300  << "will be set to an analytical value calculated from the "
301  << "specified pressure drop. No PID control will be done and "
302  << "transient effects will be ignored. This behaviour is designed "
303  << "to be appropriate for potentialFoam solutions. If you are "
304  << "getting this warning from another solver, you have probably "
305  << "specified an incorrect pressure name."
306  << endl << endl;
307  }
308 
309  // Target and measured flow rates
310  scalar QTarget, QMeasured;
311  const scalar a = (1/sqr(Ab) - 1/sqr(Aa))/(2*rho);
312  if (!mesh.steady() && db().foundObject<volScalarField>(pName_))
313  {
314  const scalar b = LbyA/deltaT;
315  const scalar c = - LbyA/deltaT*oldQ_ /* - deltaP */;
316  QTarget = (- b + sqrt(sqr(b) - 4*a*(c - deltaP_)))/(2*a);
317  QMeasured = (- b + sqrt(sqr(b) - 4*a*(c - deltaP)))/(2*a);
318  }
319  else
320  {
321  QTarget = sqrt(deltaP_/a);
322  QMeasured = sqrt(deltaP/a);
323  }
324 
325  // Errors
326  error_ = QTarget - QMeasured;
327  errorIntegral_ = oldErrorIntegral_ + 0.5*(error_ + oldError_);
328  const scalar errorDifferential = oldError_ - error_;
329 
330  // Update the field
331  operator==
332  (
333  - patch().nf()
334  *(
335  QTarget
336  + P_*error_
337  + I_*errorIntegral_
338  + D_*errorDifferential
339  )/patchA
340  );
341 
342  // Log output
343  if (debug)
344  {
345  const dimensionSet pDimensions(phi.dimensions()*dimVelocity/dimArea);
346  const scalar error = deltaP/deltaP_ - 1;
347  const scalar newQ = - gSum(*this & patch().Sf());
348  Info<< "pressurePIDControlInletVelocityFvPatchField " << patch().name()
349  << endl << " "
350  << dimensionedScalar("U", dimVelocity, newQ/patchA)
351  << endl << " "
352  << dimensionedScalar("deltaP", pDimensions, deltaP)
353  << " (" << mag(error)*100 << "\% "
354  << (error < 0 ? "below" : "above") << " the target)" << endl;
355  }
356 
358 }
359 
360 
362 (
363  Ostream& os
364 ) const
365 {
367 
368  os.writeKeyword("deltaP") << deltaP_ << token::END_STATEMENT << nl;
369  os.writeKeyword("upstream") << upstreamName_ << token::END_STATEMENT << nl;
370  os.writeKeyword("downstream")
371  << downstreamName_ << token::END_STATEMENT << nl;
372  os.writeKeyword("shapeFactor") << shapeFactor_
373  << token::END_STATEMENT << nl;
374  writeEntryIfDifferent<word>(os, "p", "p", pName_);
375  writeEntryIfDifferent<word>(os, "rho", "none", rhoName_);
376  os.writeKeyword("P") << P_ << token::END_STATEMENT << nl;
377  os.writeKeyword("I") << I_ << token::END_STATEMENT << nl;
378  os.writeKeyword("D") << D_ << token::END_STATEMENT << nl;
379  os.writeKeyword("error") << error_ << token::END_STATEMENT << nl;
380  os.writeKeyword("errorIntegral")
381  << errorIntegral_ << token::END_STATEMENT << nl;
382 
383  writeEntry("value", os);
384 }
385 
386 
387 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
388 
389 namespace Foam
390 {
392  (
394  pressurePIDControlInletVelocityFvPatchVectorField
395  );
396 }
397 
398 
399 // ************************************************************************* //
volFields.H
Foam::fvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:349
Foam::fvMesh::Cf
const surfaceVectorField & Cf() const
Return face centres as surfaceVectorField.
Definition: fvMeshGeometry.C:380
Foam::token::END_STATEMENT
@ END_STATEMENT
Definition: token.H:99
steadyStateDdtScheme.H
p
p
Definition: pEqn.H:62
Foam::dimVelocity
const dimensionSet dimVelocity
Foam::dimDensity
const dimensionSet dimDensity
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
surfaceFields.H
Foam::surfaceFields.
fvPatchFieldMapper.H
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:564
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
syncTools.H
Foam::dimArea
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
Foam::constant::physicoChemical::b
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:28
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
timeIndex
label timeIndex
Definition: getTimeIndex.H:4
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::log
dimensionedScalar log(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:253
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
pressurePIDControlInletVelocityFvPatchVectorField.H
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
Foam::makePatchTypeField
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
Foam::fvPatchField< Type >::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:296
Foam::fvSchemes::steady
bool steady() const
Return true if the default ddtScheme is steadyState.
Definition: fvSchemes.H:137
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
rho
rho
Definition: pEqn.H:3
Foam::fvPatchVectorField
fvPatchField< vector > fvPatchVectorField
Definition: fvPatchFieldsFwd.H:41
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
Foam::linearInterpolate
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > linearInterpolate(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: linear.H:107
Foam::readScalar
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:52
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:259
write
Tcoeff write()
lookup
stressControl lookup("compactNormalStress") >> compactNormalStress