54 initialPoints_(
p.localPoints()),
74 rhoName_(
dict.getOrDefault<
word>(
"rho",
"rho")),
79 if (rhoName_ ==
"rhoInf")
81 dict.readEntry(
"rhoInf", rhoInf_);
84 if (
dict.readIfPresent(
"g", g_))
94 if (
dict.found(
"initialPoints"))
96 initialPoints_ = vectorField(
"initialPoints", dict , p.size());
100 initialPoints_ = p.localPoints();
115 motion_(ptf.motion_),
116 initialPoints_(ptf.initialPoints_, mapper),
117 rhoInf_(ptf.rhoInf_),
118 rhoName_(ptf.rhoName_),
119 lookupGravity_(ptf.lookupGravity_),
133 motion_(ptf.motion_),
134 initialPoints_(ptf.initialPoints_),
135 rhoInf_(ptf.rhoInf_),
136 rhoName_(ptf.rhoName_),
137 lookupGravity_(ptf.lookupGravity_),
152 initialPoints_.autoMap(m);
167 initialPoints_.rmap(sDoFptf.initialPoints_, addr);
178 if (lookupGravity_ < 0)
180 if (
db().time().foundObject<uniformDimensionedVectorField>(
"g"))
182 if (lookupGravity_ == -2)
185 <<
"Specifying the value of g in this boundary condition "
186 <<
"when g is available from the database is considered "
187 <<
"a fatal error to avoid the possibility of inconsistency"
203 const pointPatch& ptPatch = this->
patch();
206 bool firstIter =
false;
207 if (curTimeIndex_ != t.timeIndex())
210 curTimeIndex_ = t.timeIndex();
216 forcesDict.add(
"type", functionObjects::forces::typeName);
217 forcesDict.add(
"patches",
wordList(1, ptPatch.name()));
218 forcesDict.add(
"rhoInf", rhoInf_);
219 forcesDict.add(
"rho", rhoName_);
220 forcesDict.add(
"CofR", motion_.centreOfRotation());
222 functionObjects::forces
f(
"forces",
db(), forcesDict);
224 f.calcForcesMoments();
228 if (lookupGravity_ == 1)
242 ramp*(
f.forceEff() + motion_.mass()*g_),
243 ramp*(
f.momentEff() + motion_.mass()*(motion_.momentArm() ^ g_)),
250 motion_.transform(initialPoints_) - initialPoints_
261 os.writeEntry(
"rho", rhoName_);
263 if (rhoName_ ==
"rhoInf")
265 os.writeEntry(
"rhoInf", rhoInf_);
268 if (lookupGravity_ == 0 || lookupGravity_ == -2)
270 os.writeEntry(
"g", g_);
275 initialPoints_.writeEntry(
"initialPoints",
os);
286 sixDoFRigidBodyDisplacementPointPatchVectorField
Macros for easy insertion into run-time selection tables.
const uniformDimensionedVectorField & g
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const Mesh & mesh() const noexcept
Return const reference to mesh.
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
void autoMap(const FieldMapper &map, const bool applyFlip=true)
Map from self.
constexpr Field() noexcept
void rmap(const UList< Type > &mapF, const labelUList &mapAddressing)
1 to 1 reverse-map from the given field
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
label timeIndex() const noexcept
Return the current time index.
scalar deltaTValue() const noexcept
Return time step value.
scalar deltaT0Value() const noexcept
Return old time step value.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
entry * add(entry *entryPtr, bool mergeEntry=false)
Add a new entry.
const Type & value() const noexcept
Return const reference to value.
A FixedValue boundary condition for pointField.
fixedValuePointPatchField(const pointPatch &, const DimensionedField< vector, pointMesh > &)
Computes forces and moments over a given list of patches by integrating pressure and viscous forces a...
const Time & time() const noexcept
Return time registry.
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
const word & name() const noexcept
The patch name.
const Time & time() const
Return Time from polyMesh.
const pointPatch & patch() const noexcept
Return the patch.
const objectRegistry & db() const
The associated objectRegistry.
bool updated() const noexcept
True if the boundary condition has already been updated.
Foam::pointPatchFieldMapper.
Abstract base class for point-mesh patch fields.
virtual void write(Ostream &os) const
Write.
const DimensionedField< vector, pointMesh > & internalField() const noexcept
Basic pointPatch represents a set of points from the mesh.
Mesh consisting of general polyhedral cells.
Foam::sixDoFRigidBodyDisplacementPointPatchVectorField.
virtual void write(Ostream &) const
Write.
virtual void autoMap(const pointPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void rmap(const pointPatchField< vector > &, const labelList &)
Reverse map the given pointPatchField onto this pointPatchField.
sixDoFRigidBodyDisplacementPointPatchVectorField(const pointPatch &, const DimensionedField< vector, pointMesh > &)
Construct from patch and internal field.
void writeValueEntry(Ostream &os) const
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)
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
List< word > wordList
List of word.
List< label > labelList
A List of labels.
UniformDimensionedField< vector > uniformDimensionedVectorField
pointPatchField< vector > pointPatchVectorField
static constexpr const zero Zero
Global zero (0).
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
errorManipArg< error, int > exit(error &err, const int errNo=1)
#define makePointPatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete pointPatchField type and add to run-time tables Example, (pointPatchScalarField,...