40Foam::pressurePIDControlInletVelocityFvPatchVectorField::facePressure()
const
42 const word pfName(pName_ +
"f");
85 Q_(-
gSum(*this & patch().Sf())),
105 upstreamName_(ptf.upstreamName_),
106 downstreamName_(ptf.downstreamName_),
107 deltaP_(ptf.deltaP_),
108 shapeFactor_(ptf.shapeFactor_),
110 phiName_(ptf.phiName_),
111 rhoName_(ptf.rhoName_),
117 errorIntegral_(ptf.errorIntegral_),
119 oldError_(ptf.oldError_),
120 oldErrorIntegral_(ptf.oldErrorIntegral_),
121 timeIndex_(ptf.timeIndex_)
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)),
149 oldErrorIntegral_(0),
161 upstreamName_(ptf.upstreamName_),
162 downstreamName_(ptf.downstreamName_),
163 deltaP_(ptf.deltaP_),
164 shapeFactor_(ptf.shapeFactor_),
166 phiName_(ptf.phiName_),
167 rhoName_(ptf.rhoName_),
173 errorIntegral_(ptf.errorIntegral_),
175 oldError_(ptf.oldError_),
176 oldErrorIntegral_(ptf.oldErrorIntegral_),
177 timeIndex_(ptf.timeIndex_)
189 upstreamName_(ptf.upstreamName_),
190 downstreamName_(ptf.downstreamName_),
191 deltaP_(ptf.deltaP_),
192 shapeFactor_(ptf.shapeFactor_),
194 phiName_(ptf.phiName_),
195 rhoName_(ptf.rhoName_),
201 errorIntegral_(ptf.errorIntegral_),
203 oldError_(ptf.oldError_),
204 oldErrorIntegral_(ptf.oldErrorIntegral_),
205 timeIndex_(ptf.timeIndex_)
219 const fvMesh&
mesh(patch().boundaryMesh().
mesh());
222 const scalar deltaT(db().time().deltaTValue());
228 if (timeIndex_ != db().time().
timeIndex())
230 timeIndex_ = db().time().timeIndex();
233 oldErrorIntegral_ = errorIntegral_;
251 <<
"The dimensions of the field " << phiName_
252 <<
"are not recognised. The dimensions are " <<
phi.dimensions()
254 <<
" for an incompressible case, or "
260 const scalar patchA =
gSum(
patch().magSf());
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;
277 scalar deltaP = deltaP_;
278 if (db().foundObject<volScalarField>(pName_))
281 faceZoneAverage(upstreamName_, facePressure(), Aa, pa);
282 faceZoneAverage(downstreamName_, facePressure(), Ab, pb);
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."
300 scalar QTarget, QMeasured;
301 const scalar a = (1/
sqr(Ab) - 1/
sqr(Aa))/(2*
rho);
302 if (!
mesh.
steady() && db().foundObject<volScalarField>(pName_))
304 const scalar
b = LbyA/deltaT;
305 const scalar
c = - LbyA/deltaT*oldQ_ ;
306 QTarget = (-
b +
sqrt(
sqr(
b) - 4*a*(c - deltaP_)))/(2*a);
307 QMeasured = (-
b +
sqrt(
sqr(
b) - 4*a*(c - deltaP)))/(2*a);
311 QTarget =
sqrt(deltaP_/a);
312 QMeasured =
sqrt(deltaP/a);
316 error_ = QTarget - QMeasured;
317 errorIntegral_ = oldErrorIntegral_ + 0.5*(error_ + oldError_);
318 const scalar errorDifferential = oldError_ - error_;
328 + D_*errorDifferential
336 const scalar
error = deltaP/deltaP_ - 1;
337 const scalar newQ = -
gSum(*
this &
patch().Sf());
338 Info<<
"pressurePIDControlInletVelocityFvPatchField " <<
patch().name()
344 << (
error < 0 ?
"below" :
"above") <<
" the target)" <<
endl;
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_);
381 pressurePIDControlInletVelocityFvPatchVectorField
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...
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Class to handle errors and exceptions in a simple, consistent stream-based manner.
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
Mesh data needed to do the Finite Volume discretisation.
const surfaceVectorField & Cf() const
Return face centres as surfaceVectorField.
A FieldMapper for finite-volume patch fields.
virtual void write(Ostream &) const
Write.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
This boundary condition tries to generate an inlet velocity that maintains a specified pressure drop ...
virtual void write(Ostream &) const
Write.
pressurePIDControlInletVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Lookup type of boundary radiation properties.
bool steady() const noexcept
True if default ddt scheme is steady-state.
A class for handling words, derived from Foam::string.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete fvPatchField type and add to run-time tables Example, (fvPatchScalarField,...
#define WarningInFunction
Report a warning using Foam::Warning.
const dimensionedScalar c
Speed of light in a vacuum.
Namespace for handling debugging switches.
const std::string patch
OpenFOAM patch number as a std::string.
Type gSum(const FieldField< Field, Type > &f)
Type gWeightedAverage(const UList< scalar > &weights, const UList< Type > &fld, const label comm)
The global weighted average of a field, using the mag() of the weights.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const dimensionSet dimArea(sqr(dimLength))
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.
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)
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
const dimensionSet dimVolume(pow3(dimLength))
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
errorManipArg< error, int > exit(error &err, const int errNo=1)
fvPatchField< vector > fvPatchVectorField
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
const vector L(dict.get< vector >("L"))