31#include "surfaceInterpolate.H"
44Foam::layeredEngineMesh::layeredEngineMesh(
const IOobject&
io)
63 scalar deltaZ = engineDB_.pistonDisplacement().value();
64 Info<<
"deltaZ = " << deltaZ <<
endl;
67 scalar pistonPlusLayers = pistonPosition_.value() + pistonLayers_.value();
73 point&
p = newPoints[pointi];
75 if (
p.z() < pistonPlusLayers)
79 else if (
p.z() < deckHeight_.value())
83 *(deckHeight_.value() -
p.z())
84 /(deckHeight_.value() - pistonPlusLayers);
98 const bool absolutePhi = moving();
105 movePoints(newPoints);
115 movePoints(newPoints);
118 pistonPosition_.value() += deltaZ;
119 scalar pistonSpeed = deltaZ/engineDB_.deltaTValue();
121 Info<<
"clearance: " << deckHeight_.value() - pistonPosition_.value() <<
nl
122 <<
"Piston speed = " << pistonSpeed <<
" m/s" <<
endl;
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
propsDict readIfPresent("fields", acceptFields)
Defines the attributes of an object for which implicit objectRegistry management is supported,...
dimensionedScalar deckHeight_
const engineTime & engineDB_
dimensionedScalar pistonPosition_
const surfaceScalarField & phi() const
Return cell face motion fluxes.
virtual void movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
~layeredEngineMesh()
Destructor.
bool moving() const noexcept
Is mesh moving.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Calculate the mesh motion flux and convert fluxes from absolute to relative and back.
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< surfaceScalarField > meshPhi(const volVectorField &U)
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Ostream & endl(Ostream &os)
Add newline and flush stream.
vector point
Point is a vector.
static constexpr const zero Zero
Global zero (0).
vectorField pointField
pointField is a vectorField.
constexpr char nl
The newline '\n' character (0x0a).
#define forAll(list, i)
Loop across all elements in list.