Newer
Older
mattijs
committed
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2015 OpenCFD Ltd.
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
\*---------------------------------------------------------------------------*/
#include "pressurePIDControlInletVelocityFvPatchVectorField.H"
#include "volFields.H"
#include "addToRunTimeSelectionTable.H"
#include "fvPatchFieldMapper.H"
#include "surfaceFields.H"
#include "linear.H"
#include "steadyStateDdtScheme.H"
mattijs
committed
#include "syncTools.H"
mattijs
committed
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
// * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
const Foam::surfaceScalarField&
Foam::pressurePIDControlInletVelocityFvPatchVectorField::facePressure() const
{
const volScalarField& p(db().lookupObject<volScalarField>(pName_));
const word pfName(pName_ + "f");
if (!db().foundObject<surfaceScalarField>(pfName))
{
surfaceScalarField* pfPtr
(
new surfaceScalarField(pfName, linearInterpolate(p))
);
pfPtr->store();
}
surfaceScalarField& pf
(
const_cast<surfaceScalarField&>
(
db().lookupObject<surfaceScalarField>(pfName)
)
);
if (!pf.upToDate(p))
{
pf = linearInterpolate(p);
}
return pf;
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::pressurePIDControlInletVelocityFvPatchVectorField::
pressurePIDControlInletVelocityFvPatchVectorField
(
const fvPatch& p,
const DimensionedField<vector, volMesh>& iF
)
:
fixedValueFvPatchField<vector>(p, iF),
upstreamName_(word::null),
downstreamName_(word::null),
deltaP_(1),
shapeFactor_(0),
pName_("p"),
phiName_("phi"),
rhoName_("none"),
P_(0),
I_(0),
D_(0),
Q_(- gSum(*this & patch().Sf())),
error_(0),
errorIntegral_(0),
oldQ_(0),
oldError_(0),
oldErrorIntegral_(0),
timeIndex_(db().time().timeIndex())
{}
Foam::pressurePIDControlInletVelocityFvPatchVectorField::
pressurePIDControlInletVelocityFvPatchVectorField
(
const pressurePIDControlInletVelocityFvPatchVectorField& ptf,
const fvPatch& p,
const DimensionedField<vector, volMesh>& iF,
const fvPatchFieldMapper& mapper
)
:
fixedValueFvPatchField<vector>(ptf, p, iF, mapper),
upstreamName_(ptf.upstreamName_),
downstreamName_(ptf.downstreamName_),
deltaP_(ptf.deltaP_),
shapeFactor_(ptf.shapeFactor_),
pName_(ptf.pName_),
phiName_(ptf.phiName_),
rhoName_(ptf.rhoName_),
P_(ptf.P_),
I_(ptf.I_),
D_(ptf.D_),
Q_(ptf.Q_),
error_(ptf.error_),
errorIntegral_(ptf.errorIntegral_),
oldQ_(ptf.oldQ_),
oldError_(ptf.oldError_),
oldErrorIntegral_(ptf.oldErrorIntegral_),
timeIndex_(ptf.timeIndex_)
{}
Foam::pressurePIDControlInletVelocityFvPatchVectorField::
pressurePIDControlInletVelocityFvPatchVectorField
(
const fvPatch& p,
const DimensionedField<vector, volMesh>& iF,
const dictionary& dict
)
:
fixedValueFvPatchField<vector>(p, iF, dict),
upstreamName_(dict.lookup("upstream")),
downstreamName_(dict.lookup("downstream")),
deltaP_(readScalar(dict.lookup("deltaP"))),
shapeFactor_(dict.lookupOrDefault<scalar>("shapeFactor", 0)),
pName_(dict.lookupOrDefault<word>("p", "p")),
phiName_(dict.lookupOrDefault<word>("phi", "phi")),
rhoName_(dict.lookupOrDefault<word>("rho", "none")),
P_(readScalar(dict.lookup("P"))),
I_(readScalar(dict.lookup("I"))),
D_(readScalar(dict.lookup("D"))),
Q_(- gSum(*this & patch().Sf())),
error_(dict.lookupOrDefault<scalar>("error", 0)),
errorIntegral_(dict.lookupOrDefault<scalar>("errorIntegral", 0)),
oldQ_(0),
oldError_(0),
oldErrorIntegral_(0),
timeIndex_(db().time().timeIndex())
{}
Foam::pressurePIDControlInletVelocityFvPatchVectorField::
pressurePIDControlInletVelocityFvPatchVectorField
(
const pressurePIDControlInletVelocityFvPatchVectorField& ptf
)
:
fixedValueFvPatchField<vector>(ptf),
upstreamName_(ptf.upstreamName_),
downstreamName_(ptf.downstreamName_),
deltaP_(ptf.deltaP_),
shapeFactor_(ptf.shapeFactor_),
pName_(ptf.pName_),
phiName_(ptf.phiName_),
rhoName_(ptf.rhoName_),
P_(ptf.P_),
I_(ptf.I_),
D_(ptf.D_),
Q_(ptf.Q_),
error_(ptf.error_),
errorIntegral_(ptf.errorIntegral_),
oldQ_(ptf.oldQ_),
oldError_(ptf.oldError_),
oldErrorIntegral_(ptf.oldErrorIntegral_),
timeIndex_(ptf.timeIndex_)
{}
Foam::pressurePIDControlInletVelocityFvPatchVectorField::
pressurePIDControlInletVelocityFvPatchVectorField
(
const pressurePIDControlInletVelocityFvPatchVectorField& ptf,
const DimensionedField<vector, volMesh>& iF
)
:
fixedValueFvPatchField<vector>(ptf, iF),
upstreamName_(ptf.upstreamName_),
downstreamName_(ptf.downstreamName_),
deltaP_(ptf.deltaP_),
shapeFactor_(ptf.shapeFactor_),
pName_(ptf.pName_),
phiName_(ptf.phiName_),
rhoName_(ptf.rhoName_),
P_(ptf.P_),
I_(ptf.I_),
D_(ptf.D_),
Q_(ptf.Q_),
error_(ptf.error_),
errorIntegral_(ptf.errorIntegral_),
oldQ_(ptf.oldQ_),
oldError_(ptf.oldError_),
oldErrorIntegral_(ptf.oldErrorIntegral_),
timeIndex_(ptf.timeIndex_)
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
void Foam::pressurePIDControlInletVelocityFvPatchVectorField::updateCoeffs()
{
if (updated())
{
return;
}
// Get the mesh
const fvMesh& mesh(patch().boundaryMesh().mesh());
// Get the time step
const scalar deltaT(db().time().deltaTValue());
// Get the flux field
const surfaceScalarField& phi
(
db().lookupObject<surfaceScalarField>(phiName_)
);
// Update the old-time quantities
if (timeIndex_ != db().time().timeIndex())
{
timeIndex_ = db().time().timeIndex();
oldQ_ = Q_;
oldError_ = error_;
oldErrorIntegral_ = errorIntegral_;
}
// Get the density
scalar rho = 1;
if (phi.dimensions() == dimVelocity*dimArea)
{
// do nothing ...
}
else if (phi.dimensions() == dimDensity*dimVelocity*dimArea)
{
const fvPatchField<scalar>& rhoField =
patch().lookupPatchField<volScalarField, scalar>(rhoName_);
rho = gSum(rhoField*patch().magSf())/gSum(patch().magSf());
}
else
{
FatalErrorIn
(
"void Foam::"
"pressurePIDControlInletVelocityFvPatchVectorField::"
"updateCoeffs()"
) << "The dimensions of the field " << phiName_
<< "are not recognised. The dimensions are " << phi.dimensions()
<< ". The dimensions should be either " << dimVelocity*dimArea
<< " for an incompressible case, or "
<< dimDensity*dimVelocity*dimArea << " for a compressible case."
<< exit(FatalError);
}
// Patch properties
const scalar patchA = gSum(patch().magSf());
Q_ = - gSum(*this & patch().Sf());
// Face-zone properties (a is upstream, b is downstream)
scalar Aa, Ab;
vector xa, xb;
faceZoneAverage(upstreamName_, mesh.Cf(), Aa, xa);
faceZoneAverage(downstreamName_, mesh.Cf(), Ab, xb);
const scalar L = mag(xa - xb);
const scalar LbyALinear = L/(Aa - Ab)*log(Aa/Ab);
const scalar LbyAStep = L/2*(1/Aa + 1/Ab);
const scalar LbyA = (1 - shapeFactor_)*LbyALinear + shapeFactor_*LbyAStep;
// Initialise the pressure drop. If the pressure field does not exist, the
// pressure drop is assumed to be that specified. This removes the error,
// so there is no control and the analytic inlet velocity is applied. This
// scenario only ever going to be applicable to potentialFoam.
scalar deltaP = deltaP_;
if (db().foundObject<volScalarField>(pName_))
{
scalar pa, pb;
faceZoneAverage(upstreamName_, facePressure(), Aa, pa);
faceZoneAverage(downstreamName_, facePressure(), Ab, pb);
deltaP = pa - pb;
}
else
{
WarningIn
(
"void Foam::pressurePIDControlInletVelocityFvPatchVectorField::"
"updateCoeffs()"
) << "The pressure field name, \"pName\", is \"" << pName_ << "\", "
<< "but a field of that name was not found. The inlet velocity "
<< "will be set to an analytical value calculated from the "
<< "specified pressure drop. No PID control will be done and "
<< "transient effects will be ignored. This behaviour is designed "
<< "to be appropriate for potentialFoam solutions. If you are "
<< "getting this warning from another solver, you have probably "
<< "specified an incorrect pressure name."
<< endl << endl;
}
// Target and measured flow rates
scalar QTarget, QMeasured;
const scalar a = (1/sqr(Ab) - 1/sqr(Aa))/(2*rho);
if (!mesh.steady() && db().foundObject<volScalarField>(pName_))
{
const scalar b = LbyA/deltaT;
const scalar c = - LbyA/deltaT*oldQ_ /* - deltaP */;
QTarget = (- b + sqrt(sqr(b) - 4*a*(c - deltaP_)))/(2*a);
QMeasured = (- b + sqrt(sqr(b) - 4*a*(c - deltaP)))/(2*a);
}
else
{
QTarget = sqrt(deltaP_/a);
QMeasured = sqrt(deltaP/a);
}
// Errors
error_ = QTarget - QMeasured;
errorIntegral_ = oldErrorIntegral_ + 0.5*(error_ + oldError_);
const scalar errorDifferential = oldError_ - error_;
// Update the field
operator==
(
- patch().nf()
*(
QTarget
+ P_*error_
+ I_*errorIntegral_
+ D_*errorDifferential
)/patchA
);
// Log output
if (debug)
{
const dimensionSet pDimensions(phi.dimensions()*dimVelocity/dimArea);
const scalar error = deltaP/deltaP_ - 1;
const scalar newQ = - gSum(*this & patch().Sf());
Info<< "pressurePIDControlInletVelocityFvPatchField " << patch().name()
<< endl << " "
<< dimensionedScalar("U", dimVelocity, newQ/patchA)
<< endl << " "
<< dimensionedScalar("deltaP", pDimensions, deltaP)
<< " (" << mag(error)*100 << "\% "
<< (error < 0 ? "below" : "above") << " the target)" << endl;
}
fixedValueFvPatchField<vector>::updateCoeffs();
}
void Foam::pressurePIDControlInletVelocityFvPatchVectorField::write
(
Ostream& os
) const
{
fvPatchField<vector>::write(os);
os.writeKeyword("deltaP") << deltaP_ << token::END_STATEMENT << nl;
os.writeKeyword("upstream") << upstreamName_ << token::END_STATEMENT << nl;
os.writeKeyword("downstream")
<< downstreamName_ << token::END_STATEMENT << nl;
os.writeKeyword("shapeFactor") << shapeFactor_
<< token::END_STATEMENT << nl;
writeEntryIfDifferent<word>(os, "p", "p", pName_);
writeEntryIfDifferent<word>(os, "rho", "none", rhoName_);
os.writeKeyword("P") << P_ << token::END_STATEMENT << nl;
os.writeKeyword("I") << I_ << token::END_STATEMENT << nl;
os.writeKeyword("D") << D_ << token::END_STATEMENT << nl;
os.writeKeyword("error") << error_ << token::END_STATEMENT << nl;
os.writeKeyword("errorIntegral")
<< errorIntegral_ << token::END_STATEMENT << nl;
writeEntry("value", os);
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
makePatchTypeField
(
fvPatchVectorField,
pressurePIDControlInletVelocityFvPatchVectorField
);
}
// ************************************************************************* //