In: src/postProcessing/functionObjects/forces, adding mechanisms for

restraints and constraints.

Moving sixDoFRigidBodyMotion back one directory.

Adding uncoupledSixDoFRigidBodyDisplacement to execute motion
specified by the six DoF motion state, but not applying surface
forces.  Useful for pre-displacing a mesh.

Adding constrainedSixDoFRigidBodyDisplacement to temporarily perform
Vorticity's lander simulation.  Will be removed when generalised
constraints are added.
This commit is contained in:
graham
2010-01-25 18:35:53 +00:00
parent 129f796211
commit c194087671
19 changed files with 1325 additions and 26 deletions

View File

@ -4,10 +4,21 @@ forces/forcesFunctionObject.C
forceCoeffs/forceCoeffs.C
forceCoeffs/forceCoeffsFunctionObject.C
sDoFRBM = pointPatchFields/derived/sixDoFRigidBodyMotion
$(sDoFRBM)/sixDoFRigidBodyMotion.C
$(sDoFRBM)/sixDoFRigidBodyMotionIO.C
$(sDoFRBM)/sixDoFRigidBodyMotionState.C
$(sDoFRBM)/sixDoFRigidBodyMotionStateIO.C
sDoFRBMR = $(sDoFRBM)/sixDoFRigidBodyMotionRestraint
$(sDoFRBMR)/sixDoFRigidBodyMotionRestraint/sixDoFRigidBodyMotionRestraint.C
$(sDoFRBMR)/sixDoFRigidBodyMotionRestraint/newSixDoFRigidBodyMotionRestraint.C
$(sDoFRBMR)/linearSpring/linearSpring.C
pointPatchFields/derived/sixDoFRigidBodyDisplacement/sixDoFRigidBodyDisplacementPointPatchVectorField.C
pointPatchFields/derived/sixDoFRigidBodyDisplacement/sixDoFRigidBodyMotion/sixDoFRigidBodyMotion.C
pointPatchFields/derived/sixDoFRigidBodyDisplacement/sixDoFRigidBodyMotion/sixDoFRigidBodyMotionIO.C
pointPatchFields/derived/sixDoFRigidBodyDisplacement/sixDoFRigidBodyMotion/sixDoFRigidBodyMotionState.C
pointPatchFields/derived/sixDoFRigidBodyDisplacement/sixDoFRigidBodyMotion/sixDoFRigidBodyMotionStateIO.C
pointPatchFields/derived/uncoupledSixDoFRigidBodyDisplacement/uncoupledSixDoFRigidBodyDisplacementPointPatchVectorField.C
pointPatchFields/derived/constrainedSixDoFRigidBodyDisplacement/constrainedSixDoFRigidBodyDisplacementPointPatchVectorField.C
LIB = $(FOAM_LIBBIN)/libforces

View File

@ -0,0 +1,206 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 1991-2009 OpenCFD Ltd.
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software; you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by the
Free Software Foundation; either version 2 of the License, or (at your
option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM; if not, write to the Free Software Foundation,
Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
\*---------------------------------------------------------------------------*/
#include "constrainedSixDoFRigidBodyDisplacementPointPatchVectorField.H"
#include "pointPatchFields.H"
#include "addToRunTimeSelectionTable.H"
#include "Time.H"
#include "fvMesh.H"
#include "volFields.H"
#include "uniformDimensionedFields.H"
#include "forces.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
constrainedSixDoFRigidBodyDisplacementPointPatchVectorField::
constrainedSixDoFRigidBodyDisplacementPointPatchVectorField
(
const pointPatch& p,
const DimensionedField<vector, pointMesh>& iF
)
:
fixedValuePointPatchField<vector>(p, iF),
motion_(),
p0_(p.localPoints()),
rhoInf_(1.0)
{}
constrainedSixDoFRigidBodyDisplacementPointPatchVectorField::
constrainedSixDoFRigidBodyDisplacementPointPatchVectorField
(
const pointPatch& p,
const DimensionedField<vector, pointMesh>& iF,
const dictionary& dict
)
:
fixedValuePointPatchField<vector>(p, iF, dict),
motion_(dict),
rhoInf_(readScalar(dict.lookup("rhoInf")))
{
if (!dict.found("value"))
{
updateCoeffs();
}
if (dict.found("p0"))
{
p0_ = vectorField("p0", dict , p.size());
}
else
{
p0_ = p.localPoints();
}
}
constrainedSixDoFRigidBodyDisplacementPointPatchVectorField::
constrainedSixDoFRigidBodyDisplacementPointPatchVectorField
(
const constrainedSixDoFRigidBodyDisplacementPointPatchVectorField& ptf,
const pointPatch& p,
const DimensionedField<vector, pointMesh>& iF,
const pointPatchFieldMapper& mapper
)
:
fixedValuePointPatchField<vector>(ptf, p, iF, mapper),
motion_(ptf.motion_),
p0_(ptf.p0_, mapper),
rhoInf_(ptf.rhoInf_)
{}
constrainedSixDoFRigidBodyDisplacementPointPatchVectorField::
constrainedSixDoFRigidBodyDisplacementPointPatchVectorField
(
const constrainedSixDoFRigidBodyDisplacementPointPatchVectorField& ptf,
const DimensionedField<vector, pointMesh>& iF
)
:
fixedValuePointPatchField<vector>(ptf, iF),
motion_(ptf.motion_),
p0_(ptf.p0_),
rhoInf_(ptf.rhoInf_)
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
void constrainedSixDoFRigidBodyDisplacementPointPatchVectorField::updateCoeffs()
{
if (this->updated())
{
return;
}
const polyMesh& mesh = this->dimensionedInternalField().mesh()();
const Time& t = mesh.time();
const pointPatch& ptPatch = this->patch();
// Patch force data is valid for the current positions, so
// calculate the forces on the motion object from this data, then
// update the positions
motion_.updatePosition(t.deltaTValue());
dictionary forcesDict;
forcesDict.add("patches", wordList(1, ptPatch.name()));
forcesDict.add("rhoInf", rhoInf_);
forcesDict.add("CofR", motion_.centreOfMass());
forces f("forces", db(), forcesDict);
forces::forcesMoments fm = f.calcForcesMoment();
// Get the forces on the patch faces at the current positions
vector gravity = vector::zero;
if (db().foundObject<uniformDimensionedVectorField>("g"))
{
uniformDimensionedVectorField g =
db().lookupObject<uniformDimensionedVectorField>("g");
gravity = g.value();
}
vector rotationAxis(0, 1, 0);
vector torque
(
(
(fm.second().first() + fm.second().second())
& rotationAxis
)
*rotationAxis
);
motion_.updateForce
(
vector::zero, // Force no centre of mass motion
torque, // Only rotation allowed around the unit rotationAxis
t.deltaTValue()
);
Field<vector>::operator=(motion_.generatePositions(p0_) - p0_);
fixedValuePointPatchField<vector>::updateCoeffs();
}
void constrainedSixDoFRigidBodyDisplacementPointPatchVectorField::write
(
Ostream& os
) const
{
pointPatchField<vector>::write(os);
motion_.write(os);
os.writeKeyword("rhoInf")
<< rhoInf_ << token::END_STATEMENT << nl;
p0_.writeEntry("p0", os);
writeEntry("value", os);
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
makePointPatchTypeField
(
pointPatchVectorField,
constrainedSixDoFRigidBodyDisplacementPointPatchVectorField
);
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace Foam
// ************************************************************************* //

View File

@ -0,0 +1,157 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 1991-2009 OpenCFD Ltd.
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software; you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by the
Free Software Foundation; either version 2 of the License, or (at your
option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM; if not, write to the Free Software Foundation,
Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
Class
Foam::constrainedSixDoFRigidBodyDisplacementPointPatchVectorField
Description
Foam::constrainedSixDoFRigidBodyDisplacementPointPatchVectorField
SourceFiles
constrainedSixDoFRigidBodyDisplacementPointPatchVectorField.C
\*---------------------------------------------------------------------------*/
#ifndef constrainedSixDoFRigidBodyDisplacementPointPatchVectorField_H
#define constrainedSixDoFRigidBodyDisplacementPointPatchVectorField_H
#include "fixedValuePointPatchField.H"
#include "sixDoFRigidBodyMotion.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
/*---------------------------------------------------------------------------*\
Class constrainedSixDoFRigidBodyDisplacementPointPatchVectorField Declaration
\*---------------------------------------------------------------------------*/
class constrainedSixDoFRigidBodyDisplacementPointPatchVectorField
:
public fixedValuePointPatchField<vector>
{
// Private data
//- Six dof motion object
sixDoFRigidBodyMotion motion_;
//- Reference positions of points on the patch
pointField p0_;
//- Reference density required by the forces object for
// incompressible calculations
scalar rhoInf_;
public:
//- Runtime type information
TypeName("constrainedSixDoFRigidBodyDisplacement");
// Constructors
//- Construct from patch and internal field
constrainedSixDoFRigidBodyDisplacementPointPatchVectorField
(
const pointPatch&,
const DimensionedField<vector, pointMesh>&
);
//- Construct from patch, internal field and dictionary
constrainedSixDoFRigidBodyDisplacementPointPatchVectorField
(
const pointPatch&,
const DimensionedField<vector, pointMesh>&,
const dictionary&
);
//- Construct by mapping given patchField<vector> onto a new patch
constrainedSixDoFRigidBodyDisplacementPointPatchVectorField
(
const constrainedSixDoFRigidBodyDisplacementPointPatchVectorField&,
const pointPatch&,
const DimensionedField<vector, pointMesh>&,
const pointPatchFieldMapper&
);
//- Construct and return a clone
virtual autoPtr<pointPatchField<vector> > clone() const
{
return autoPtr<pointPatchField<vector> >
(
new constrainedSixDoFRigidBodyDisplacementPointPatchVectorField
(
*this
)
);
}
//- Construct as copy setting internal field reference
constrainedSixDoFRigidBodyDisplacementPointPatchVectorField
(
const constrainedSixDoFRigidBodyDisplacementPointPatchVectorField&,
const DimensionedField<vector, pointMesh>&
);
//- Construct and return a clone setting internal field reference
virtual autoPtr<pointPatchField<vector> > clone
(
const DimensionedField<vector, pointMesh>& iF
) const
{
return autoPtr<pointPatchField<vector> >
(
new constrainedSixDoFRigidBodyDisplacementPointPatchVectorField
(
*this,
iF
)
);
}
// Member functions
// Evaluation functions
//- Update the coefficients associated with the patch field
virtual void updateCoeffs();
//- Write
virtual void write(Ostream&) const;
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //

View File

@ -160,6 +160,26 @@ void sixDoFRigidBodyDisplacementPointPatchVectorField::updateCoeffs()
t.deltaTValue()
);
// ----------------------------------------
// vector rotationAxis(0, 1, 0);
// vector torque
// (
// (
// (fm.second().first() + fm.second().second())
// & rotationAxis
// )
// *rotationAxis
// );
// motion_.updateForce
// (
// vector::zero, // Force no centre of mass motion
// torque, // Only rotation allowed around the unit rotationAxis
// t.deltaTValue()
// );
// ----------------------------------------
Field<vector>::operator=(motion_.generatePositions(p0_) - p0_);
fixedValuePointPatchField<vector>::updateCoeffs();

View File

@ -26,11 +26,29 @@ License
#include "sixDoFRigidBodyMotion.H"
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
void Foam::sixDoFRigidBodyMotion::applyRestraints()
{
forAll(restraints_, rI)
{
restraints_[rI].restraintForce();
}
}
void Foam::sixDoFRigidBodyMotion::applyConstraints()
{
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::sixDoFRigidBodyMotion::sixDoFRigidBodyMotion()
:
motionState_(),
restraints_(),
refCentreOfMass_(vector::zero),
momentOfInertia_(diagTensor::one*VSMALL),
mass_(VSMALL)
@ -59,6 +77,7 @@ Foam::sixDoFRigidBodyMotion::sixDoFRigidBodyMotion
pi,
tau
),
restraints_(),
refCentreOfMass_(refCentreOfMass),
momentOfInertia_(momentOfInertia),
mass_(mass)
@ -68,10 +87,15 @@ Foam::sixDoFRigidBodyMotion::sixDoFRigidBodyMotion
Foam::sixDoFRigidBodyMotion::sixDoFRigidBodyMotion(const dictionary& dict)
:
motionState_(dict),
restraints_(),
refCentreOfMass_(dict.lookupOrDefault("refCentreOfMass", centreOfMass())),
momentOfInertia_(dict.lookup("momentOfInertia")),
mass_(readScalar(dict.lookup("mass")))
{}
{
addRestraints(dict);
addConstraints(dict);
}
Foam::sixDoFRigidBodyMotion::sixDoFRigidBodyMotion
@ -80,6 +104,7 @@ Foam::sixDoFRigidBodyMotion::sixDoFRigidBodyMotion
)
:
motionState_(sDoFRBM.motionState()),
restraints_(sDoFRBM.restraints()),
refCentreOfMass_(sDoFRBM.refCentreOfMass()),
momentOfInertia_(sDoFRBM.momentOfInertia()),
mass_(sDoFRBM.mass())
@ -94,6 +119,49 @@ Foam::sixDoFRigidBodyMotion::~sixDoFRigidBodyMotion()
// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
void Foam::sixDoFRigidBodyMotion::addRestraints
(
const dictionary& dict
)
{
if (dict.found("restraints"))
{
const dictionary& restraintDict = dict.subDict("restraints");
label i = 0;
restraints_.setSize(restraintDict.size());
forAllConstIter(IDLList<entry>, restraintDict, iter)
{
if (iter().isDict())
{
Info<< "Adding restraint: " << iter().keyword() << endl;
restraints_.set
(
i,
sixDoFRigidBodyMotionRestraint::New(iter().dict())
);
}
i++;
}
restraints_.setSize(i);
}
}
void Foam::sixDoFRigidBodyMotion::addConstraints
(
const dictionary& dict
)
{
}
void Foam::sixDoFRigidBodyMotion::updatePosition
(
scalar deltaT
@ -149,7 +217,7 @@ void Foam::sixDoFRigidBodyMotion::updateForce
)
{
// Second leapfrog velocity adjust part, required after motion and
// force calculation part
// force calculation
if (Pstream::master())
{
@ -157,6 +225,10 @@ void Foam::sixDoFRigidBodyMotion::updateForce
tau() = (Q().T() & tauGlobal);
applyRestraints();
applyConstraints();
v() += 0.5*deltaT*a();
pi() += 0.5*deltaT*tau();
@ -173,30 +245,23 @@ void Foam::sixDoFRigidBodyMotion::updateForce
scalar deltaT
)
{
// Second leapfrog velocity adjust part, required after motion and
// force calculation part
vector a = vector::zero;
vector tau = vector::zero;
if (Pstream::master())
{
a() = vector::zero;
tau() = vector::zero;
forAll(positions, i)
{
const vector& f = forces[i];
a() += f/mass_;
a += f/mass_;
tau() += (positions[i] ^ (Q().T() & f));
tau += (positions[i] ^ (Q().T() & f));
}
v() += 0.5*deltaT*a();
pi() += 0.5*deltaT*tau();
}
Pstream::scatter(motionState_);
updateForce(a, tau, deltaT);
}

View File

@ -26,10 +26,10 @@ Class
Foam::sixDoFRigidBodyMotion
Description
Six degree of freedom motion for a rigid body. Angular momentum stored in
body fixed reference frame. Reference orientation of the body must align
with the cartesian axes such that the Inertia tensor is in principle
component form.
Six degree of freedom motion for a rigid body. Angular momentum
stored in body fixed reference frame. Reference orientation of
the body (where Q = I) must align with the cartesian axes such
that the Inertia tensor is in principle component form.
Symplectic motion as per:
@ -43,6 +43,9 @@ Description
url = {http://link.aip.org/link/?JCP/107/5840/1},
doi = {10.1063/1.474310}
Can add restraints (i.e. a spring) and constraints (i.e. motion
may only be on a plane).
SourceFiles
sixDoFRigidBodyMotionI.H
sixDoFRigidBodyMotion.C
@ -55,6 +58,8 @@ SourceFiles
#include "sixDoFRigidBodyMotionState.H"
#include "pointField.H"
#include "sixDoFRigidBodyMotionRestraint.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
@ -79,9 +84,15 @@ class sixDoFRigidBodyMotion
{
// Private data
// state data object
//- Motion state data object
sixDoFRigidBodyMotionState motionState_;
//- Restraints on the motion
PtrList<sixDoFRigidBodyMotionRestraint> restraints_;
// //- Constaints on the motion
// PtrList<sixDoFRigidBodyMotionConstraint> constraints_;
//- Centre of mass of reference state
point refCentreOfMass_;
@ -106,6 +117,12 @@ class sixDoFRigidBodyMotion
// frame z-axis by the given angle
inline tensor rotationTensorZ(scalar deltaT) const;
//- Apply the restraints to the object
void applyRestraints();
//- Apply the constraints to the object
void applyConstraints();
public:
@ -141,11 +158,23 @@ public:
// Member Functions
//- Add restraints to the motion, public to allow external
// addition of restraints after construction
void addRestraints(const dictionary& dict);
//- Add restraints to the motion, public to allow external
// addition of restraints after construction
void addConstraints(const dictionary& dict);
//- First leapfrog velocity adjust and motion part, required
// before force calculation
void updatePosition
(
scalar deltaT
);
//- Second leapfrog velocity adjust part, required after motion and
// force calculation
void updateForce
(
const vector& fGlobal,
@ -153,6 +182,8 @@ public:
scalar deltaT
);
//- Global forces supplied at locations, calculating net force
// and moment
void updateForce
(
const pointField& positions,
@ -160,6 +191,7 @@ public:
scalar deltaT
);
//- Transform the given pointField by the current motion state
tmp<pointField> generatePositions(const pointField& pts) const;
// Access
@ -167,6 +199,10 @@ public:
//- Return access to the motion state
inline const sixDoFRigidBodyMotionState& motionState() const;
//- Return access to the restraints
inline const PtrList<sixDoFRigidBodyMotionRestraint>&
restraints() const;
//- Return access to the centre of mass
inline const point& centreOfMass() const;

View File

@ -24,8 +24,6 @@ License
\*---------------------------------------------------------------------------*/
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
inline Foam::tensor
@ -71,6 +69,13 @@ Foam::sixDoFRigidBodyMotion::motionState() const
}
inline const Foam::PtrList<Foam::sixDoFRigidBodyMotionRestraint>&
Foam::sixDoFRigidBodyMotion::restraints() const
{
return restraints_;
}
inline const Foam::point& Foam::sixDoFRigidBodyMotion::centreOfMass() const
{
return motionState_.centreOfMass();

View File

@ -0,0 +1,87 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 1991-2009 OpenCFD Ltd.
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software; you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by the
Free Software Foundation; either version 2 of the License, or (at your
option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM; if not, write to the Free Software Foundation,
Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
\*---------------------------------------------------------------------------*/
#include "linearSpring.H"
#include "addToRunTimeSelectionTable.H"
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
namespace Foam
{
namespace sixDoFRigidBodyMotionRestraints
{
defineTypeNameAndDebug(linearSpring, 0);
addToRunTimeSelectionTable
(
sixDoFRigidBodyMotionRestraint,
linearSpring,
dictionary
);
};
};
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::sixDoFRigidBodyMotionRestraints::linearSpring::linearSpring
(
const dictionary& sDoFRBMRCoeffs
)
:
sixDoFRigidBodyMotionRestraint(sDoFRBMRCoeffs)
{
read(sDoFRBMRCoeffs);
}
// * * * * * * * * * * * * * * * * Destructors * * * * * * * * * * * * * * * //
Foam::sixDoFRigidBodyMotionRestraints::linearSpring::~linearSpring()
{}
// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
void
Foam::sixDoFRigidBodyMotionRestraints::linearSpring::restraintForce() const
{
}
bool Foam::sixDoFRigidBodyMotionRestraints::linearSpring::read
(
const dictionary& sDoFRBMRCoeffs
)
{
sixDoFRigidBodyMotionRestraint::read(sDoFRBMRCoeffs);
// sDoFRBMRCoeffs_.lookup("velocity") >> velocity_;
return true;
}
// ************************************************************************* //

View File

@ -0,0 +1,121 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 1991-2009 OpenCFD Ltd.
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software; you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by the
Free Software Foundation; either version 2 of the License, or (at your
option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM; if not, write to the Free Software Foundation,
Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
Class
Foam::sixDoFRigidBodyMotionRestraints::linearSpring
Description
sixDoFRigidBodyMotionRestraints model. Linear spring.
SourceFiles
linearSpring.C
\*---------------------------------------------------------------------------*/
#ifndef linearSpring_H
#define linearSpring_H
#include "sixDoFRigidBodyMotionRestraint.H"
#include "point.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
namespace sixDoFRigidBodyMotionRestraints
{
/*---------------------------------------------------------------------------*\
Class linearSpring Declaration
\*---------------------------------------------------------------------------*/
class linearSpring
:
public sixDoFRigidBodyMotionRestraint
{
// Private data
//- Anchor point, where the spring is attached to an immovable
// object
point anchor_;
//- Reference point of attachment to the solid body
point refAttachmentPt_;
//- Spring stiffness
scalar stiffness_;
//- Rest length - length of spring when no forces are applied to it
scalar restLength_;
public:
//- Runtime type information
TypeName("linearSpring");
// Constructors
//- Construct from components
linearSpring
(
const dictionary& sDoFRBMRCoeffs
);
//- Construct and return a clone
virtual autoPtr<sixDoFRigidBodyMotionRestraint> clone() const
{
return autoPtr<sixDoFRigidBodyMotionRestraint>
(
new linearSpring(*this)
);
}
// Destructor
virtual ~linearSpring();
// Member Functions
//- Return the restraint force and point of application
virtual void restraintForce() const;
//- Update properties from given dictionary
virtual bool read(const dictionary& sDoFRBMRCoeff);
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace solidBodyMotionFunctions
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //

View File

@ -0,0 +1,65 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 1991-2009 OpenCFD Ltd.
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software; you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by the
Free Software Foundation; either version 2 of the License, or (at your
option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM; if not, write to the Free Software Foundation,
Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
\*---------------------------------------------------------------------------*/
#include "sixDoFRigidBodyMotionRestraint.H"
// * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
Foam::autoPtr<Foam::sixDoFRigidBodyMotionRestraint>
Foam::sixDoFRigidBodyMotionRestraint::New(const dictionary& sDoFRBMRCoeffs)
{
word sixDoFRigidBodyMotionRestraintTypeName =
sDoFRBMRCoeffs.lookup("sixDoFRigidBodyMotionRestraint");
Info<< "Selecting sixDoFRigidBodyMotionRestraint function "
<< sixDoFRigidBodyMotionRestraintTypeName << endl;
dictionaryConstructorTable::iterator cstrIter =
dictionaryConstructorTablePtr_->find
(
sixDoFRigidBodyMotionRestraintTypeName
);
if (cstrIter == dictionaryConstructorTablePtr_->end())
{
FatalErrorIn
(
"sixDoFRigidBodyMotionRestraint::New"
"("
" const dictionary& sDoFRBMRCoeffs"
")"
) << "Unknown sixDoFRigidBodyMotionRestraint type "
<< sixDoFRigidBodyMotionRestraintTypeName << endl << endl
<< "Valid sixDoFRigidBodyMotionRestraints are : " << endl
<< dictionaryConstructorTablePtr_->sortedToc()
<< exit(FatalError);
}
return autoPtr<sixDoFRigidBodyMotionRestraint>(cstrIter()(sDoFRBMRCoeffs));
}
// ************************************************************************* //

View File

@ -0,0 +1,72 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 1991-2009 OpenCFD Ltd.
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software; you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by the
Free Software Foundation; either version 2 of the License, or (at your
option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM; if not, write to the Free Software Foundation,
Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
\*---------------------------------------------------------------------------*/
#include "sixDoFRigidBodyMotionRestraint.H"
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
defineTypeNameAndDebug(Foam::sixDoFRigidBodyMotionRestraint, 0);
defineRunTimeSelectionTable(Foam::sixDoFRigidBodyMotionRestraint, dictionary);
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::sixDoFRigidBodyMotionRestraint::sixDoFRigidBodyMotionRestraint
(
const dictionary& sDoFRBMRCoeffs
)
:
sDoFRBMRCoeffs_
(
sDoFRBMRCoeffs.subDict
(
word(sDoFRBMRCoeffs.lookup("sixDoFRigidBodyMotionRestraint"))
+ "Coeffs"
)
)
{}
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
Foam::sixDoFRigidBodyMotionRestraint::~sixDoFRigidBodyMotionRestraint()
{}
// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
bool Foam::sixDoFRigidBodyMotionRestraint::read
(
const dictionary& sDoFRBMRCoeffs
)
{
sDoFRBMRCoeffs_ = sDoFRBMRCoeffs.subDict(type() + "Coeffs");
return true;
}
// ************************************************************************* //

View File

@ -0,0 +1,133 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 1991-2009 OpenCFD Ltd.
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software; you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by the
Free Software Foundation; either version 2 of the License, or (at your
option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM; if not, write to the Free Software Foundation,
Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
Namespace
Foam::sixDoFRigidBodyMotionRestraints
Description
Namespace for six DoF motion restraints
Class
Foam::sixDoFRigidBodyMotionRestraint
Description
Base class for defining restraints for sixDoF motions
SourceFiles
sixDoFRigidBodyMotionRestraint.C
newDynamicFvMesh.C
\*---------------------------------------------------------------------------*/
#ifndef sixDoFRigidBodyMotionRestraint_H
#define sixDoFRigidBodyMotionRestraint_H
#include "Time.H"
#include "dictionary.H"
#include "autoPtr.H"
#include "runTimeSelectionTables.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
/*---------------------------------------------------------------------------*\
Class sixDoFRigidBodyMotionRestraint Declaration
\*---------------------------------------------------------------------------*/
class sixDoFRigidBodyMotionRestraint
{
protected:
// Protected data
dictionary sDoFRBMRCoeffs_;
public:
//- Runtime type information
TypeName("sixDoFRigidBodyMotionRestraint");
// Declare run-time constructor selection table
declareRunTimeSelectionTable
(
autoPtr,
sixDoFRigidBodyMotionRestraint,
dictionary,
(const dictionary& sDoFRBMRCoeffs),
(sDoFRBMRCoeffs)
);
// Constructors
//- Construct from the sDoFRBMRCoeffs dictionary and Time
sixDoFRigidBodyMotionRestraint
(
const dictionary& sDoFRBMRCoeffs
);
//- Construct and return a clone
virtual autoPtr<sixDoFRigidBodyMotionRestraint> clone() const = 0;
// Selectors
//- Select constructed from the sDoFRBMRCoeffs dictionary and Time
static autoPtr<sixDoFRigidBodyMotionRestraint> New
(
const dictionary& sDoFRBMRCoeffs
);
// Destructor
virtual ~sixDoFRigidBodyMotionRestraint();
// Member Functions
//- Return the restraint force and point of application
virtual void restraintForce() const = 0;
//- Update properties from given dictionary
virtual bool read(const dictionary& sDoFRBMRCoeffs) = 0;
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //

View File

@ -0,0 +1,163 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 1991-2009 OpenCFD Ltd.
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software; you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by the
Free Software Foundation; either version 2 of the License, or (at your
option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM; if not, write to the Free Software Foundation,
Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
\*---------------------------------------------------------------------------*/
#include "uncoupledSixDoFRigidBodyDisplacementPointPatchVectorField.H"
#include "pointPatchFields.H"
#include "addToRunTimeSelectionTable.H"
#include "Time.H"
#include "fvMesh.H"
#include "volFields.H"
#include "uniformDimensionedFields.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
uncoupledSixDoFRigidBodyDisplacementPointPatchVectorField::
uncoupledSixDoFRigidBodyDisplacementPointPatchVectorField
(
const pointPatch& p,
const DimensionedField<vector, pointMesh>& iF
)
:
fixedValuePointPatchField<vector>(p, iF),
motion_(),
p0_(p.localPoints()),
rhoInf_(1.0)
{}
uncoupledSixDoFRigidBodyDisplacementPointPatchVectorField::
uncoupledSixDoFRigidBodyDisplacementPointPatchVectorField
(
const pointPatch& p,
const DimensionedField<vector, pointMesh>& iF,
const dictionary& dict
)
:
fixedValuePointPatchField<vector>(p, iF, dict),
motion_(dict),
rhoInf_(readScalar(dict.lookup("rhoInf")))
{
if (!dict.found("value"))
{
updateCoeffs();
}
if (dict.found("p0"))
{
p0_ = vectorField("p0", dict , p.size());
}
else
{
p0_ = p.localPoints();
}
}
uncoupledSixDoFRigidBodyDisplacementPointPatchVectorField::
uncoupledSixDoFRigidBodyDisplacementPointPatchVectorField
(
const uncoupledSixDoFRigidBodyDisplacementPointPatchVectorField& ptf,
const pointPatch& p,
const DimensionedField<vector, pointMesh>& iF,
const pointPatchFieldMapper& mapper
)
:
fixedValuePointPatchField<vector>(ptf, p, iF, mapper),
motion_(ptf.motion_),
p0_(ptf.p0_, mapper),
rhoInf_(ptf.rhoInf_)
{}
uncoupledSixDoFRigidBodyDisplacementPointPatchVectorField::
uncoupledSixDoFRigidBodyDisplacementPointPatchVectorField
(
const uncoupledSixDoFRigidBodyDisplacementPointPatchVectorField& ptf,
const DimensionedField<vector, pointMesh>& iF
)
:
fixedValuePointPatchField<vector>(ptf, iF),
motion_(ptf.motion_),
p0_(ptf.p0_),
rhoInf_(ptf.rhoInf_)
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
void uncoupledSixDoFRigidBodyDisplacementPointPatchVectorField::updateCoeffs()
{
if (this->updated())
{
return;
}
const polyMesh& mesh = this->dimensionedInternalField().mesh()();
const Time& t = mesh.time();
motion_.updatePosition(t.deltaTValue());
// Do not modify the accelerations
motion_.updateForce(vector::zero, vector::zero, t.deltaTValue());
Field<vector>::operator=(motion_.generatePositions(p0_) - p0_);
fixedValuePointPatchField<vector>::updateCoeffs();
}
void uncoupledSixDoFRigidBodyDisplacementPointPatchVectorField::write
(
Ostream& os
) const
{
pointPatchField<vector>::write(os);
motion_.write(os);
os.writeKeyword("rhoInf")
<< rhoInf_ << token::END_STATEMENT << nl;
p0_.writeEntry("p0", os);
writeEntry("value", os);
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
makePointPatchTypeField
(
pointPatchVectorField,
uncoupledSixDoFRigidBodyDisplacementPointPatchVectorField
);
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace Foam
// ************************************************************************* //

View File

@ -0,0 +1,158 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 1991-2009 OpenCFD Ltd.
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software; you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by the
Free Software Foundation; either version 2 of the License, or (at your
option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM; if not, write to the Free Software Foundation,
Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
Class
Foam::uncoupledSixDoFRigidBodyDisplacementPointPatchVectorField
Description
Foam::uncoupledSixDoFRigidBodyDisplacementPointPatchVectorField
SourceFiles
uncoupledSixDoFRigidBodyDisplacementPointPatchVectorField.C
\*---------------------------------------------------------------------------*/
#ifndef uncoupledSixDoFRigidBodyDisplacementPointPatchVectorField_H
#define uncoupledSixDoFRigidBodyDisplacementPointPatchVectorField_H
#include "fixedValuePointPatchField.H"
#include "sixDoFRigidBodyMotion.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
/*---------------------------------------------------------------------------*\
Class uncoupledSixDoFRigidBodyDisplacementPointPatchVectorField Declaration
\*---------------------------------------------------------------------------*/
class uncoupledSixDoFRigidBodyDisplacementPointPatchVectorField
:
public fixedValuePointPatchField<vector>
{
// Private data
//- Six dof motion object
sixDoFRigidBodyMotion motion_;
//- Reference positions of points on the patch
pointField p0_;
//- Reference density required by the forces object for
// incompressible calculations. Retained here to give
// dictionary compatibility with other sixDoF patches.
scalar rhoInf_;
public:
//- Runtime type information
TypeName("uncoupledSixDoFRigidBodyDisplacement");
// Constructors
//- Construct from patch and internal field
uncoupledSixDoFRigidBodyDisplacementPointPatchVectorField
(
const pointPatch&,
const DimensionedField<vector, pointMesh>&
);
//- Construct from patch, internal field and dictionary
uncoupledSixDoFRigidBodyDisplacementPointPatchVectorField
(
const pointPatch&,
const DimensionedField<vector, pointMesh>&,
const dictionary&
);
//- Construct by mapping given patchField<vector> onto a new patch
uncoupledSixDoFRigidBodyDisplacementPointPatchVectorField
(
const uncoupledSixDoFRigidBodyDisplacementPointPatchVectorField&,
const pointPatch&,
const DimensionedField<vector, pointMesh>&,
const pointPatchFieldMapper&
);
//- Construct and return a clone
virtual autoPtr<pointPatchField<vector> > clone() const
{
return autoPtr<pointPatchField<vector> >
(
new uncoupledSixDoFRigidBodyDisplacementPointPatchVectorField
(
*this
)
);
}
//- Construct as copy setting internal field reference
uncoupledSixDoFRigidBodyDisplacementPointPatchVectorField
(
const uncoupledSixDoFRigidBodyDisplacementPointPatchVectorField&,
const DimensionedField<vector, pointMesh>&
);
//- Construct and return a clone setting internal field reference
virtual autoPtr<pointPatchField<vector> > clone
(
const DimensionedField<vector, pointMesh>& iF
) const
{
return autoPtr<pointPatchField<vector> >
(
new uncoupledSixDoFRigidBodyDisplacementPointPatchVectorField
(
*this,
iF
)
);
}
// Member functions
// Evaluation functions
//- Update the coefficients associated with the patch field
virtual void updateCoeffs();
//- Write
virtual void write(Ostream&) const;
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //