Go to the documentation of this file.
69 angle0_(
dict.
get<scalar>(
"angle0")),
70 amplitude_(
dict.
get<scalar>(
"amplitude")),
71 omega_(
dict.
get<scalar>(
"omega"))
84 p0_ =
p.localPoints();
100 origin_(ptf.origin_),
101 angle0_(ptf.angle0_),
102 amplitude_(ptf.amplitude_),
117 origin_(ptf.origin_),
118 angle0_(ptf.angle0_),
119 amplitude_(ptf.amplitude_),
145 refCast<const angularOscillatingVelocityPointPatchVectorField>(ptf);
149 p0_.rmap(aOVptf.p0_, addr);
164 scalar angle = angle0_ + amplitude_*
sin(omega_*t.
value());
172 + p0Rel*(
cos(angle) - 1)
173 + (axisHat ^ p0Rel*
sin(angle))
174 + (axisHat & p0Rel)*(1 -
cos(angle))*axisHat
195 writeEntry(
"value",
os);
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Foam::angularOscillatingVelocityPointPatchVectorField.
static constexpr const zero Zero
Global zero (0)
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
dimensionedScalar sin(const dimensionedScalar &ds)
friend Ostream & operator(Ostream &, const Field< vector > &)
const Type & value() const
Return const reference to value.
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Basic pointPatch represents a set of points from the mesh.
virtual void write(Ostream &) const
Write.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Mesh consisting of general polyhedral cells.
const pointPatch & patch() const
Return patch.
Foam::pointPatchFieldMapper.
scalar deltaTValue() const noexcept
Return time step value.
Field< vector > vectorField
Specialisation of Field<T> for vector.
makePointPatchTypeField(pointPatchVectorField, solidBodyMotionDisplacementPointPatchVectorField)
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
const DimensionedField< vector, pointMesh > & internalField() const
Return dimensioned internal field reference.
virtual void write(Ostream &) const
Write.
ITstream & lookup(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
OBJstream os(runTime.globalPath()/outputName)
Macros for easy insertion into run-time selection tables.
virtual void autoMap(const pointPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
bool updated() const
Return true if the boundary condition has already been updated.
virtual void rmap(const pointPatchField< Type > &, const labelList &)
Reverse map the given PointPatchField onto.
angularOscillatingVelocityPointPatchVectorField(const pointPatch &, const DimensionedField< vector, pointMesh > &)
Construct from patch and internal field.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
virtual void autoMap(const pointPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
virtual void rmap(const pointPatchField< vector > &, const labelList &)
Reverse map the given pointPatchField onto this pointPatchField.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
dimensionedScalar cos(const dimensionedScalar &ds)