43 fixedValueFvPatchScalarField(
p, iF),
57 Ap_(
dict.get<scalar>(
"Ap")),
58 Sp_(
dict.get<scalar>(
"Sp")),
59 VsI_(
dict.get<scalar>(
"VsI")),
60 tas_(
dict.get<scalar>(
"tas")),
61 tae_(
dict.get<scalar>(
"tae")),
62 tds_(
dict.get<scalar>(
"tds")),
63 tde_(
dict.get<scalar>(
"tde")),
64 psI_(
dict.get<scalar>(
"psI")),
65 psi_(
dict.get<scalar>(
"psi")),
66 ams_(
dict.get<scalar>(
"ams")),
68 phiName_(
dict.getOrDefault<
word>(
"phi",
"phi")),
71 scalar ps = (psI_*VsI_ + ams_/psi_)/Vs(db().time().value());
84 fixedValueFvPatchScalarField(sppsf,
p, iF, mapper),
96 phiName_(sppsf.phiName_),
107 fixedValueFvPatchScalarField(sppsf, iF),
119 phiName_(sppsf.phiName_),
129 fixedValueFvPatchScalarField(sppsf),
141 phiName_(sppsf.phiName_),
148Foam::scalar Foam::syringePressureFvPatchScalarField::Vs(
const scalar t)
const
158 + 0.5*Ap_*Sp_*
sqr(t - tas_)/(tae_ - tas_);
164 + 0.5*Ap_*Sp_*(tae_ - tas_)
165 + Ap_*Sp_*(t - tae_);
171 + 0.5*Ap_*Sp_*(tae_ - tas_)
172 + Ap_*Sp_*(tds_ - tae_)
174 - 0.5*Ap_*Sp_*
sqr(t - tds_)/(tde_ - tds_);
180 + 0.5*Ap_*Sp_*(tae_ - tas_)
181 + Ap_*Sp_*(tds_ - tae_)
182 + 0.5*Ap_*Sp_*(tde_ - tds_);
194 if (curTimeIndex_ != db().time().
timeIndex())
197 curTimeIndex_ = db().time().timeIndex();
200 scalar t = db().time().value();
201 scalar deltaT = db().time().deltaTValue();
207 ams_ = ams0_ + deltaT*
sum((*
this*psi_)*phip);
211 ams_ = ams0_ + deltaT*
sum(phip);
216 <<
"dimensions of phi are not correct"
217 <<
"\n on patch " << this->
patch().name()
218 <<
" of field " << this->internalField().name()
219 <<
" in file " << this->internalField().objectPath()
223 scalar ps = (psI_*VsI_ + ams_/psi_)/Vs(t);
227 fixedValueFvPatchScalarField::updateCoeffs();
235 os.writeEntry(
"Ap", Ap_);
236 os.writeEntry(
"Sp", Sp_);
237 os.writeEntry(
"VsI", VsI_);
238 os.writeEntry(
"tas", tas_);
239 os.writeEntry(
"tae", tae_);
240 os.writeEntry(
"tds", tds_);
241 os.writeEntry(
"tde", tde_);
242 os.writeEntry(
"psI", psI_);
243 os.writeEntry(
"psi", psi_);
244 os.writeEntry(
"ams", ams_);
245 os.writeEntryIfDifferent<
word>(
"phi",
"phi", phiName_);
258 syringePressureFvPatchScalarField
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...
A simple container of IOobject preferences. Can also be used for general handling of read/no-read/rea...
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
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 operator=(const UList< Type > &)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
This boundary condition provides a pressure condition, obtained from a zero-D model of the cylinder o...
virtual void write(Ostream &) const
Write.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
syringePressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
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,...
const std::string patch
OpenFOAM patch number as a std::string.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1, const label comm)
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
const dimensionSet dimVolume(pow3(dimLength))
errorManipArg< error, int > exit(error &err, const int errNo=1)
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
fvPatchField< scalar > fvPatchScalarField