CG  Version 25
Public Types | Public Member Functions | Public Attributes | Static Public Attributes | Protected Types | Protected Member Functions | Protected Attributes | List of all members
RigidBodyMotion Class Reference

#include <RigidBodyMotion.h>

Public Types

enum  PositionConstraintEnum { positionHasNoConstraint =0, positionConstrainedToAPlane, positionConstrainedToALine, positionIsFixed }
 
enum  RotationConstraintEnum { rotationHasNoConstraint =0, rotationConstrainedToAPlane, rotationConstrainedToALine, rotationIsFixed }
 
enum  TimeSteppingMethodEnum { leapFrogTrapezoidal =0, improvedEuler, implicitRungeKutta }
 
enum  TwilightZoneTypeEnum { polynomialTwilightZone =0, trigonometricTwilightZone }
 
enum  BodyForceTypeEnum { timePolynomialBodyForce =0, timeFunctionBodyForce }
 
enum  { defaultOrderOfAccuracy =-1234 }
 

Public Member Functions

 RigidBodyMotion (int numberOfDimensions=3)
 
 ~RigidBodyMotion ()
 
bool centerOfMassHasBeenInitialized () const
 
int correct (real t, const RealArray &force, const RealArray &torque, RealArray &xCM, RealArray &rotation)
 Correct the solution at time t using new values for the forces and torques.
 
int correct (real t, const RealArray &force, const RealArray &torque, const RealArray &A11, const RealArray &A12, const RealArray &A21, const RealArray &A22, RealArray &xCM, RealArray &rotation)
 Correct the solution at time t using new values for the forces and torques. This version takes added mass matrices.
 
int get (const GenericDataBase &dir, const aString &name)
 
int getAcceleration (real t, RealArray &vCM) const
 
int getAddedMassMatrices (const real t, RealArray &A11, RealArray &A12, RealArray &A21, RealArray &A22) const
 : Evaluate the added mass matrices at time t.
 
int getAngularAcceleration (real t, RealArray &omegaDot) const
 
int getAngularVelocities (real t, RealArray &omega) const
 
int getAxesOfInertia (real t, RealArray &axesOfInertia) const
 
int getBodyForces (const real t, RealArray &bodyForce, RealArray &bodyTorque) const
 : Evaluate the body force and torque.
 
bool getCorrectionHasConverged ()
 Return true if the correction steps for moving grids have converged. This is usually only an issue for "light" bodies.
 
int getCoordinates (real t, RealArray &xCM=Overture::nullRealArray(), RealArray &vCM=Overture::nullRealArray(), RealArray &rotation=Overture::nullRealArray(), RealArray &omega=Overture::nullRealArray(), RealArray &omegaDot=Overture::nullRealArray(), RealArray &aCM=Overture::nullRealArray(), RealArray &axesOfInertia=Overture::nullRealArray()) const
 
real getDensity () const
 
int getExactSolution (const int deriv, const real t, RealArray &xe, RealArray &ve, RealArray &we) const
 Evaluate the exact solution (for twilightZone)
 
int getForce (const real t, RealArray &fv, RealArray &gv) const
 : Return the force and torque at time t (interpolate in time if necessary).
 
int getInitialConditions (RealArray &x0=Overture::nullRealArray(), RealArray &v0=Overture::nullRealArray(), RealArray &w0=Overture::nullRealArray(), RealArray &axesOfInertia=Overture::nullRealArray()) const
 
real getMass () const
 
real getMaximumRelativeCorrection ()
 Return the maximum relative change in the moving grid correction scheme. This is usually only an issue for "light" bodies.
 
RealArray getMomentsOfInertia () const
 
int getPosition (real t, RealArray &xCM) const
 
int getPointAcceleration (real t, const RealArray &p, RealArray &ap) const
 
int getPointVelocity (real t, const RealArray &p, RealArray &vp) const
 
int getPointTransformationMatrix (real t, RealArray &rDotRt=Overture::nullRealArray(), RealArray &rDotDotRt=Overture::nullRealArray()) const
 
real getTimeStepEstimate () const
 
aString getTimeSteppingMethodName () const
 
int getVelocity (real t, RealArray &vCM) const
 
int getRotationMatrix (real t, RealArray &r, RealArray &rDot=Overture::nullRealArray(), RealArray &rDotDot=Overture::nullRealArray()) const
 
int includeAddedMass (bool trueOrFalse=true)
 Indicate whether added mass matrices will be used (and provided by the user).
 
int integrate (real tf, const RealArray &force, const RealArray &torque, real t, RealArray &xCM, RealArray &rotation)
 Integrate the equations of motion from time t0 to time t using force and moment information at time t0.
 
int integrate (real tf, const RealArray &force, const RealArray &torque, real t, const RealArray &A11, const RealArray &A12, const RealArray &A21, const RealArray &A22, RealArray &xCM, RealArray &rotation)
 Integrate the equations of motion from time t0 to time t using force and moment information at time t0. This version takes added mass matrices.
 
bool massHasBeenInitialized () const
 
int momentumTransfer (real t, RealArray &v)
 
int put (GenericDataBase &dir, const aString &name) const
 
int reset ()
 Reset the rigid body to it's initial state. (e.g. remove saved solutions etc.)
 
int setInitialCentreOfMass (const RealArray &x0)
 
int setInitialConditions (real t0=0., const RealArray &x0=Overture::nullRealArray(), const RealArray &v0=Overture::nullRealArray(), const RealArray &w0=Overture::nullRealArray(), const RealArray &axesOfInertia=Overture::nullRealArray())
 
void setMass (const real totalMass)
 
int setNewtonTolerance (real tol)
 Set the convergence tolerance for the Newton iteration of implicit solvers.
 
int setPastTimeForcing (real t, const RealArray &force, const RealArray &torque)
 Supply forcing at "negative" times for startup.
 
int setPastTimeForcing (real t, const RealArray &force, const RealArray &torque, const RealArray &A11, const RealArray &A12, const RealArray &A21, const RealArray &A22)
 Supply forcing at "negative" times for startup.
 
int setPositionConstraint (PositionConstraintEnum constraint, const RealArray &values)
 
int setProperties (real totalMass, const RealArray &momentsOfInertia, const int numberOfDimensions=3)
 
int setRotationConstraint (RotationConstraintEnum constraint, const RealArray &values)
 
int setTwilightZone (bool trueOrFalse, TwilightZoneTypeEnum type=trigonometricTwilightZone)
 Turn on (or off) twilight-zone forcing.
 
int setTimeSteppingMethod (const TimeSteppingMethodEnum method, int orderOfAccuracy=defaultOrderOfAccuracy)
 Choose the time stepping method. Optionally set the order of accuracy.
 
int update (GenericGraphicsInterface &gi)
 
bool useAddedMass () const
 Return true if the added mass matrices are being used.
 

Public Attributes

enum
RigidBodyMotion::TimeSteppingMethodEnum 
timeSteppingMethod
 

Static Public Attributes

static int debug = 0
 

Protected Types

enum  ConstraintApplicationEnum { applyConstraintsToPosition =0, applyConstraintsToVelocity, applyConstraintsToForces }
 

Protected Member Functions

int applyConstraints (const int step, const ConstraintApplicationEnum option)
 
int buildBodyForceOptionsDialog (DialogData &dialog)
 Build the plot options dialog.
 
int dirkImplicitSolve (const real dt, const real aii, const real tc, const RealArray &yv, const RealArray &yv0, RealArray &kv)
 : Solve the implicit DIRK equation for kv. This is a protected routine.
 
int getBodyForceOption (const aString &answer, DialogData &dialog)
 : Look for a plot option in the string "answer"
 
int getForceInternal (const real t, RealArray &fv, RealArray &gv, RealArray *pA=NULL) const
 : Return the force and torque at time t (interpolate in time if necessary). This is a protected routine. The force and torque do NOT include the added mass terms.
 
int takeStepImplicitRungeKutta (const real t0, const real dt)
 Take a time step with the implicit Runge-Kutta method. This is a protected routine.
 
int takeStepLeapFrog (const real t0, const real dt)
 Take a time step with the leap-frog (predictor). This is a protected routine.
 
int takeStepTrapezoid (const real t0, const real dt)
 Take a time step with the trapezoidal rule (corrector). This is a protected routine.
 
int takeStepExtrapolation (const real t0, const real dt)
 Take a (predictor) time step by extrapolation in time. This can be used with the DIRK schemes, for example.
 

Protected Attributes

real mass
 
int numberOfDimensions
 
int current
 
int numberOfSteps
 
int numberSaved
 
int maximumNumberToSave
 
int initialConditionsGiven
 
real density
 
real damping
 
RealArray mI
 
RealArray e
 
RealArray e0
 
RealArray x
 
RealArray v
 
RealArray x0
 
RealArray v0
 
RealArray w0
 
RealArray f
 
RealArray g
 
RealArray w
 
RealArray time
 
RealArray bodyForceCoeff
 
RealArray bodyTorqueCoeff
 
const Range R
 
PositionConstraintEnum positionConstraint
 
RotationConstraintEnum rotationConstraint
 
RealArray constraintValues
 
bool relaxCorrectionSteps
 
bool correctionHasConverged
 
real maximumRelativeCorrection
 
real correctionAbsoluteToleranceForce
 
real correctionRelativeToleranceForce
 
real correctionRelaxationParameterForce
 
real correctionAbsoluteToleranceTorque
 
real correctionRelativeToleranceTorque
 
real correctionRelaxationParameterTorque
 
bool twilightZone
 
TwilightZoneTypeEnum twilightZoneType
 
FILE * logFile
 
BodyForceTypeEnum bodyForceType
 
DataBase dbase
 

Member Enumeration Documentation

anonymous enum
Enumerator
defaultOrderOfAccuracy 
Enumerator
timePolynomialBodyForce 
timeFunctionBodyForce 
Enumerator
applyConstraintsToPosition 
applyConstraintsToVelocity 
applyConstraintsToForces 
Enumerator
positionHasNoConstraint 
positionConstrainedToAPlane 
positionConstrainedToALine 
positionIsFixed 
Enumerator
rotationHasNoConstraint 
rotationConstrainedToAPlane 
rotationConstrainedToALine 
rotationIsFixed 
Enumerator
leapFrogTrapezoidal 
improvedEuler 
implicitRungeKutta 
Enumerator
polynomialTwilightZone 
trigonometricTwilightZone 

Constructor & Destructor Documentation

RigidBodyMotion::RigidBodyMotion ( int  numberOfDimensions = 3)
RigidBodyMotion::~RigidBodyMotion ( )

Member Function Documentation

int RigidBodyMotion::applyConstraints ( const int  step,
const ConstraintApplicationEnum  option 
)
protected
int RigidBodyMotion::buildBodyForceOptionsDialog ( DialogData &  dialog)
protected

Build the plot options dialog.

Parameters
dialog(input) : graphics dialog to use.

References assert(), bodyForceCoeff, bodyForceType, and bodyTorqueCoeff.

Referenced by update().

bool RigidBodyMotion::centerOfMassHasBeenInitialized ( ) const

References initialConditionsGiven.

Referenced by update().

int RigidBodyMotion::correct ( real  t,
const RealArray &  force,
const RealArray &  torque,
RealArray &  xCM,
RealArray &  rotation 
)

Correct the solution at time t using new values for the forces and torques.

This routine might be called by a predictor-corrector type algorithm.

Parameters
t(input) : time (corresponding to the end of the time step).
force(0:2)(input) : components of the total force at time t
torque(0:2)(input) : components of the torque (about the center of mass) at time t
xCM(0:2)(output) : new positions of the centre of mass at time t
rotation(0:2,0:2)(output) : rotation matrix at time t. This matrix will be orthonormal (explicitly enforced by this routine).

Referenced by main(), MovingGrids::rigidBodyMotion(), and TestRigidBody::solve().

int RigidBodyMotion::correct ( real  t,
const RealArray &  force,
const RealArray &  torque,
const RealArray &  A11,
const RealArray &  A12,
const RealArray &  A21,
const RealArray &  A22,
RealArray &  xCM,
RealArray &  rotation 
)

Correct the solution at time t using new values for the forces and torques. This version takes added mass matrices.

This routine might be called by a predictor-corrector type algorithm.

Parameters
t(input) : time (corresponding to the end of the time step).
force(0:2)(input) : components of the total force at time t
torque(0:2)(input) : components of the torque (about the center of mass) at time t
A11,A12,A21,A22(input) : added matrices of size (0:2,0:2). These are used if the use of added mass matrices has been turned on with a call to includeAddedMass.
xCM(0:2)(output) : new positions of the centre of mass at time t
rotation(0:2,0:2)(output) : rotation matrix at time t. This matrix will be orthonormal (explicitly enforced by this routine).

References applyConstraints(), applyConstraintsToForces, assert(), correctionAbsoluteToleranceForce, correctionAbsoluteToleranceTorque, correctionHasConverged, correctionRelativeToleranceForce, correctionRelativeToleranceTorque, correctionRelaxationParameterForce, correctionRelaxationParameterTorque, current, dbase, debug, e, e0, f, g, getBodyForces(), implicitRungeKutta, improvedEuler, includeAddedMass(), k, leapFrogTrapezoidal, logFile, maximumNumberToSave, maximumRelativeCorrection, n, numberOfSteps, OV_ABORT(), printF(), R, relaxCorrectionSteps, takeStepImplicitRungeKutta(), takeStepTrapezoid(), time, timeSteppingMethod, v, w, and x.

int RigidBodyMotion::dirkImplicitSolve ( const real  dt,
const real  aii,
const real  tc,
const RealArray &  yv,
const RealArray &  yv0,
RealArray &  kv 
)
protected

: Solve the implicit DIRK equation for kv. This is a protected routine.

This code is based on the matlab code in ??

References arrayToState(), dbase, debug, dot(), getAddedMassMatrices(), getCrossProductMatrix(), getExactSolution(), getForceInternal(), j, k, logFile, m, mass, mI, mr, mult(), orderOfAccuracy, printF(), r, R, solve(), trans(), twilightZone, and vvn.

Referenced by takeStepImplicitRungeKutta().

int RigidBodyMotion::get ( const GenericDataBase &  dir,
const aString &  name 
)
int RigidBodyMotion::getAcceleration ( real  t,
RealArray &  vCM 
) const
int RigidBodyMotion::getAddedMassMatrices ( const real  t,
RealArray &  A11,
RealArray &  A12,
RealArray &  A21,
RealArray &  A22 
) const

: Evaluate the added mass matrices at time t.

Interpolate from the sequence of saved values in time.

Parameters
t(input) : evaluate matrices at this time.
A11,A12,A21,A22(output) : added mass matrices.

References assert(), c1, c2, current, dbase, includeAddedMass(), k, maximumNumberToSave, numberSaved, orderOfAccuracy, OV_ABORT(), printF(), R, and time.

Referenced by dirkImplicitSolve(), getForce(), and getForceInternal().

int RigidBodyMotion::getAngularAcceleration ( real  t,
RealArray &  omegaDot 
) const
int RigidBodyMotion::getAngularVelocities ( real  t,
RealArray &  omega 
) const
int RigidBodyMotion::getAxesOfInertia ( real  t,
RealArray &  axesOfInertia 
) const
int RigidBodyMotion::getBodyForceOption ( const aString &  answer,
DialogData &  dialog 
)
protected

: Look for a plot option in the string "answer"

Parameters
answer(input) : check this command
Returns
return 1 if the command was found, 0 otherwise.

References assert(), axis, bodyForceCoeff, bodyForceType, bodyTorqueCoeff, dbase, OV_ABORT(), P, printF(), timeFunctionBodyForce, timePolynomialBodyForce, and TimeFunction::update().

Referenced by update().

int RigidBodyMotion::getBodyForces ( const real  t,
RealArray &  bodyForce,
RealArray &  bodyTorque 
) const

: Evaluate the body force and torque.

Parameters
t(input) : evaluate at this time.
bodyForce,bodyTorque(output) : body force and boy torque at time t

References axis, bodyForceCoeff, bodyForceType, bodyTorqueCoeff, dbase, TimeFunction::eval(), f, g, OV_ABORT(), printF(), R, timeFunctionBodyForce, and timePolynomialBodyForce.

Referenced by correct(), and integrate().

int RigidBodyMotion::getCoordinates ( real  t,
RealArray &  xCM = Overture::nullRealArray(),
RealArray &  vCM = Overture::nullRealArray(),
RealArray &  rotation = Overture::nullRealArray(),
RealArray &  omega = Overture::nullRealArray(),
RealArray &  omegaDot = Overture::nullRealArray(),
RealArray &  aCM = Overture::nullRealArray(),
RealArray &  axesOfInertia = Overture::nullRealArray() 
) const
bool RigidBodyMotion::getCorrectionHasConverged ( )

Return true if the correction steps for moving grids have converged. This is usually only an issue for "light" bodies.

References correctionHasConverged.

Referenced by MovingGrids::rigidBodyMotion().

real RigidBodyMotion::getDensity ( ) const
int RigidBodyMotion::getExactSolution ( const int  deriv,
const real  t,
RealArray &  xe,
RealArray &  ve,
RealArray &  we 
) const

Evaluate the exact solution (for twilightZone)

Parameters
deriv(input) : evaluate this derivative of the exact solution
t(input) : time
x,v,w(output)

References dbase, j, OV_ABORT(), polynomialTwilightZone, printF(), trigonometricTwilightZone, twilightZone, twilightZoneType, and xe.

Referenced by dirkImplicitSolve(), TestRigidBody::getExactSolution(), getForceInternal(), and TestRigidBody::initialConditions().

int RigidBodyMotion::getForce ( const real  t,
RealArray &  fv,
RealArray &  gv 
) const

: Return the force and torque at time t (interpolate in time if necessary).

Parameters
t(input) : evaluate the force and torque at this time.
fv,gv(output) : force and torque

References assert(), c1, c2, current, dbase, f, g, getAddedMassMatrices(), getAngularVelocities(), getVelocity(), includeAddedMass(), maximumNumberToSave, mult(), numberSaved, orderOfAccuracy, OV_ABORT(), printF(), R, time, v0, and w0.

Referenced by MovingGrids::rigidBodyMotion().

int RigidBodyMotion::getForceInternal ( const real  t,
RealArray &  fv,
RealArray &  gv,
RealArray *  pA = NULL 
) const
protected

: Return the force and torque at time t (interpolate in time if necessary). This is a protected routine. The force and torque do NOT include the added mass terms.

Parameters
t(input) : evaluate the force and torque at this time.
fv,gv(output) : force and torque
pA(input) : for twilightzone forcing supply A (the angular momentum matrix, h=A*w)

References assert(), c1, c2, current, dbase, f, g, getAddedMassMatrices(), getCrossProductMatrix(), getExactSolution(), includeAddedMass(), mass, maximumNumberToSave, mult(), numberSaved, orderOfAccuracy, OV_ABORT(), printF(), R, time, and twilightZone.

Referenced by dirkImplicitSolve(), takeStepLeapFrog(), and takeStepTrapezoid().

int RigidBodyMotion::getInitialConditions ( RealArray &  x0 = Overture::nullRealArray(),
RealArray &  v0 = Overture::nullRealArray(),
RealArray &  w0 = Overture::nullRealArray(),
RealArray &  axesOfInertia = Overture::nullRealArray() 
) const

References R.

Referenced by MovingGrids::moveGrids().

real RigidBodyMotion::getMass ( ) const
real RigidBodyMotion::getMaximumRelativeCorrection ( )

Return the maximum relative change in the moving grid correction scheme. This is usually only an issue for "light" bodies.

References maximumRelativeCorrection.

RealArray RigidBodyMotion::getMomentsOfInertia ( ) const

References mI.

int RigidBodyMotion::getPointAcceleration ( real  t,
const RealArray &  p,
RealArray &  ap 
) const

References R.

int RigidBodyMotion::getPointTransformationMatrix ( real  t,
RealArray &  rDotRt = Overture::nullRealArray(),
RealArray &  rDotDotRt = Overture::nullRealArray() 
) const
int RigidBodyMotion::getPointVelocity ( real  t,
const RealArray &  p,
RealArray &  vp 
) const

References R.

int RigidBodyMotion::getPosition ( real  t,
RealArray &  xCM 
) const
int RigidBodyMotion::getRotationMatrix ( real  t,
RealArray &  r,
RealArray &  rDot = Overture::nullRealArray(),
RealArray &  rDotDot = Overture::nullRealArray() 
) const

References ee, n, omega, R, and r.

Referenced by MovingGrids::moveGrids().

real RigidBodyMotion::getTimeStepEstimate ( ) const
aString RigidBodyMotion::getTimeSteppingMethodName ( ) const
int RigidBodyMotion::getVelocity ( real  t,
RealArray &  vCM 
) const
int RigidBodyMotion::includeAddedMass ( bool  trueOrFalse = true)

Indicate whether added mass matrices will be used (and provided by the user).

Parameters
trueOrFalse(input) : set to true if added mass matrices will be provided.

References dbase.

Referenced by correct(), getAddedMassMatrices(), getForce(), getForceInternal(), integrate(), setPastTimeForcing(), TestRigidBody::solve(), and takeStepLeapFrog().

int RigidBodyMotion::integrate ( real  t0,
const RealArray &  force,
const RealArray &  torque,
real  t,
RealArray &  xCM,
RealArray &  rotation 
)

Integrate the equations of motion from time t0 to time t using force and moment information at time t0.

Parameters
t0(input) : starting time for the time step.
force(0:2)(input) : components of the total force at time t0
torque(0:2)(input) : components of the torque (about the center of mass) at time t0.
t(input) : end time for the time step.
xCM(0:2)(output) : position of the centre of mass at time t
rotation(0:2,0:2)(output) : rotation matrix at time t. This matrix will be orthonormal (explicitly enforced by this routine).

Referenced by main(), MovingGrids::rigidBodyMotion(), and TestRigidBody::solve().

int RigidBodyMotion::integrate ( real  t0,
const RealArray &  force,
const RealArray &  torque,
real  t,
const RealArray &  A11,
const RealArray &  A12,
const RealArray &  A21,
const RealArray &  A22,
RealArray &  xCM,
RealArray &  rotation 
)

Integrate the equations of motion from time t0 to time t using force and moment information at time t0. This version takes added mass matrices.

Parameters
t0(input) : starting time for the time step.
force(0:2)(input) : components of the total force at time t0
torque(0:2)(input) : components of the torque (about the center of mass) at time t0.
t(input) : end time for the time step.
A11,A12,A21,A22(input) : added matrices of size (0:2,0:2). These are used if the use of added mass matrices has been turned on with a call to includeAddedMass.
xCM(0:2)(output) : position of the centre of mass at time t
rotation(0:2,0:2)(output) : rotation matrix at time t. This matrix will be orthonormal (explicitly enforced by this routine).

References applyConstraints(), applyConstraintsToForces, c1, c2, current, dbase, debug, e, e0, f, g, getBodyForces(), implicitRungeKutta, improvedEuler, includeAddedMass(), k, leapFrogTrapezoidal, logFile, maximumNumberToSave, mult(), n, numberOfSteps, numberSaved, orderOfAccuracy, OV_ABORT(), printF(), R, takeStepImplicitRungeKutta(), takeStepLeapFrog(), time, timeSteppingMethod, v, w, and x.

bool RigidBodyMotion::massHasBeenInitialized ( ) const
int RigidBodyMotion::momentumTransfer ( real  t,
RealArray &  v 
)

References current, debug, maximumNumberToSave, printF(), R, time, v, and x.

Referenced by detectCollisions().

int RigidBodyMotion::put ( GenericDataBase &  dir,
const aString &  name 
) const

References dir, e, f, v, and x.

int RigidBodyMotion::reset ( )

Reset the rigid body to it's initial state. (e.g. remove saved solutions etc.)

References bodyForceCoeff, bodyTorqueCoeff, current, dbase, e, e0, f, g, initialConditionsGiven, maximumNumberToSave, mI, numberOfSteps, numberSaved, time, v, v0, w, w0, x, and x0.

Referenced by TestRigidBody::solve().

int RigidBodyMotion::setInitialCentreOfMass ( const RealArray &  x0)

References initialConditionsGiven, R, x, and x0.

Referenced by MovingGrids::update(), and update().

int RigidBodyMotion::setInitialConditions ( real  t0 = 0.,
const RealArray &  x0 = Overture::nullRealArray(),
const RealArray &  v0 = Overture::nullRealArray(),
const RealArray &  w0 = Overture::nullRealArray(),
const RealArray &  axesOfInertia = Overture::nullRealArray() 
)
void RigidBodyMotion::setMass ( const real  totalMass)
int RigidBodyMotion::setNewtonTolerance ( real  tol)

Set the convergence tolerance for the Newton iteration of implicit solvers.

If the Newton is converging quadratically then the actual error in the Newton iteration should be tol*tol. Thus setting tol=1.e-5 will likely give an error of 1.e-10.

References dbase.

Referenced by main().

int RigidBodyMotion::setPastTimeForcing ( real  t,
const RealArray &  force,
const RealArray &  torque 
)

Supply forcing at "negative" times for startup.

The fourth-order scheme requires values at t=-dt.

Parameters
t(input) : a past time value (e.g. t=-dt) , torque (input) : force and torque at time t

Referenced by TestRigidBody::solve().

int RigidBodyMotion::setPastTimeForcing ( real  t,
const RealArray &  force,
const RealArray &  torque,
const RealArray &  A11,
const RealArray &  A12,
const RealArray &  A21,
const RealArray &  A22 
)

Supply forcing at "negative" times for startup.

The fourth-order scheme requires values at t=-dt.

Parameters
t(input) : a past time value (e.g. t=-dt)
force,torque(input) : force and torque at time t
A11,A12,A21,A22(input) : added mass matrices. (only used if added mass matrices are being used)

References current, dbase, f, g, includeAddedMass(), maximumNumberToSave, OV_ABORT(), printF(), R, and time.

int RigidBodyMotion::setPositionConstraint ( PositionConstraintEnum  constraint,
const RealArray &  values 
)
int RigidBodyMotion::setProperties ( real  totalMass,
const RealArray &  momentsOfInertia,
const int  numberOfDimensions = 3 
)
int RigidBodyMotion::setRotationConstraint ( RotationConstraintEnum  constraint,
const RealArray &  values 
)
int RigidBodyMotion::setTimeSteppingMethod ( const TimeSteppingMethodEnum  method,
int  orderOfAccuracy = defaultOrderOfAccuracy 
)

Choose the time stepping method. Optionally set the order of accuracy.

Parameters
method(input) : use this time stepping method.
orderOfAccuracy(input) : optionally choose the order of accuracy. Currently only applies to implicit Runge Kutta.

References dbase, defaultOrderOfAccuracy, and timeSteppingMethod.

Referenced by TestRigidBody::solve().

int RigidBodyMotion::setTwilightZone ( bool  trueOrFalse,
TwilightZoneTypeEnum  type = trigonometricTwilightZone 
)

Turn on (or off) twilight-zone forcing.

Parameters
trueOrFalse(input) : turn on or off
type(input) : twilight-zone type

References dbase, R, twilightZone, and twilightZoneType.

Referenced by TestRigidBody::initialConditions().

int RigidBodyMotion::takeStepExtrapolation ( const real  t0,
const real  dt 
)
protected

Take a (predictor) time step by extrapolation in time. This can be used with the DIRK schemes, for example.

Parameters
t0: current time.
dt: time step.

References c1, c2, current, dbase, e, logFile, maximumNumberToSave, numberSaved, orderOfAccuracy, R, time, v, w, and x.

int RigidBodyMotion::takeStepImplicitRungeKutta ( const real  t0,
const real  dt 
)
protected

Take a time step with the implicit Runge-Kutta method. This is a protected routine.

Parameters
t0(input) : current time (at the start of the step)
dt: time step.

References a11, applyConstraints(), applyConstraintsToPosition, assert(), c1, c2, ci, current, dbase, dirkImplicitSolve(), dot(), e, j, m, maximumNumberToSave, n, numberOfDimensions, orderOfAccuracy, OV_ABORT(), R, v, w, and x.

Referenced by correct(), and integrate().

int RigidBodyMotion::takeStepLeapFrog ( const real  t0,
const real  dt 
)
protected

Take a time step with the leap-frog (predictor). This is a protected routine.

Parameters
t0: current time.
dt: time step.

References applyConstraints(), applyConstraintsToPosition, axis, current, damping, dbase, debug, dot(), e, f, g, getCrossProductMatrix(), getForceInternal(), improvedEuler, includeAddedMass(), j, leapFrogTrapezoidal, m, mass, maximumNumberToSave, mI, mult(), n, numberOfDimensions, numberOfSteps, OV_ABORT(), printF(), R, timeSteppingMethod, trans(), twilightZone, v, v0, w, and x.

Referenced by integrate().

int RigidBodyMotion::takeStepTrapezoid ( const real  t0,
const real  dt 
)
protected

Take a time step with the trapezoidal rule (corrector). This is a protected routine.

Parameters
t0(input) : current time (at the start of the step)
dt: time step.

References applyConstraints(), applyConstraintsToPosition, assert(), axis, current, dot(), e, f, g, getCrossProductMatrix(), getForceInternal(), j, m, mass, maximumNumberToSave, mI, mult(), n, numberOfDimensions, R, trans(), twilightZone, v, w, and x.

Referenced by correct().

int RigidBodyMotion::update ( GenericGraphicsInterface &  gi)
bool RigidBodyMotion::useAddedMass ( ) const

Return true if the added mass matrices are being used.

References dbase.

Referenced by MovingGrids::rigidBodyMotion().

Member Data Documentation

RealArray RigidBodyMotion::bodyForceCoeff
protected
BodyForceTypeEnum RigidBodyMotion::bodyForceType
protected
RealArray RigidBodyMotion::bodyTorqueCoeff
protected
RealArray RigidBodyMotion::constraintValues
protected
real RigidBodyMotion::correctionAbsoluteToleranceForce
protected

Referenced by correct(), get(), RigidBodyMotion(), and update().

real RigidBodyMotion::correctionAbsoluteToleranceTorque
protected

Referenced by correct(), get(), RigidBodyMotion(), and update().

bool RigidBodyMotion::correctionHasConverged
protected
real RigidBodyMotion::correctionRelativeToleranceForce
protected

Referenced by correct(), get(), RigidBodyMotion(), and update().

real RigidBodyMotion::correctionRelativeToleranceTorque
protected

Referenced by correct(), get(), RigidBodyMotion(), and update().

real RigidBodyMotion::correctionRelaxationParameterForce
protected

Referenced by correct(), get(), RigidBodyMotion(), and update().

real RigidBodyMotion::correctionRelaxationParameterTorque
protected

Referenced by correct(), get(), RigidBodyMotion(), and update().

int RigidBodyMotion::current
protected
real RigidBodyMotion::damping
protected
DataBase RigidBodyMotion::dbase
mutableprotected
int RigidBodyMotion::debug = 0
static
real RigidBodyMotion::density
protected
RealArray RigidBodyMotion::e
protected
RealArray RigidBodyMotion::e0
protected
RealArray RigidBodyMotion::f
protected
RealArray RigidBodyMotion::g
protected
int RigidBodyMotion::initialConditionsGiven
protected
FILE* RigidBodyMotion::logFile
protected
real RigidBodyMotion::mass
protected
int RigidBodyMotion::maximumNumberToSave
protected
real RigidBodyMotion::maximumRelativeCorrection
protected
RealArray RigidBodyMotion::mI
protected
int RigidBodyMotion::numberOfDimensions
protected
int RigidBodyMotion::numberOfSteps
protected
int RigidBodyMotion::numberSaved
protected
PositionConstraintEnum RigidBodyMotion::positionConstraint
protected
const Range RigidBodyMotion::R
protected
bool RigidBodyMotion::relaxCorrectionSteps
protected

Referenced by correct(), get(), RigidBodyMotion(), and update().

RotationConstraintEnum RigidBodyMotion::rotationConstraint
protected
RealArray RigidBodyMotion::time
protected
enum RigidBodyMotion::TimeSteppingMethodEnum RigidBodyMotion::timeSteppingMethod
bool RigidBodyMotion::twilightZone
protected
TwilightZoneTypeEnum RigidBodyMotion::twilightZoneType
protected
RealArray RigidBodyMotion::v
protected
RealArray RigidBodyMotion::v0
protected
RealArray RigidBodyMotion::w
protected
RealArray RigidBodyMotion::w0
protected
RealArray RigidBodyMotion::x
protected
RealArray RigidBodyMotion::x0
protected

The documentation for this class was generated from the following files: