43 fixedValueFvPatchScalarField(
p, iF),
56 fixedValueFvPatchScalarField(
p, iF,
dict, false),
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();
211 ams_ = ams0_ + deltaT*
sum((*
this*psi_)*phip);
215 ams_ = ams0_ + deltaT*
sum(phip);
220 <<
"dimensions of phi are not correct"
221 <<
"\n on patch " << this->patch().name()
222 <<
" of field " << this->internalField().name()
223 <<
" in file " << this->internalField().objectPath()
227 scalar ps = (psI_*VsI_ + ams_/psi_)/Vs(t);
231 fixedValueFvPatchScalarField::updateCoeffs();
251 writeEntry(
"value",
os);
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...
const dimensionSet & dimensions() const
Return dimensions.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
virtual bool write()
Write the output fields.
A FieldMapper for finite-volume patch fields.
virtual void operator=(const UList< Type > &)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
This boundary condition provides a pressure condition, obtained from a zero-D model of the cylinder o...
virtual void updateCoeffs()
Update the coefficients associated with the patch 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)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
const dimensionSet dimArea(sqr(dimLength))
const dimensionSet dimVelocity
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
const dimensionSet dimDensity
errorManipArg< error, int > exit(error &err, const int errNo=1)