walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers Class Reference

Detailed Description

Particular implementation of the collision resolution for the hard contacts.

The following code example illustrates the setup of the solver:

cr::HCSITS hcsits(globalBodyStorage, forest->getBlockStoragePointer(), storageID, ccdID, fcdID);
hcsits.setMaxIterations ( uint_c(10) );
hcsits.setRelaxationParameter ( real_t(0.7) );
hcsits.setErrorReductionParameter ( real_t(0.8) );
hcsits.setGlobalLinearAcceleration( Vec3(0,0,-1) );

#include <HCSITS.h>

+ Inheritance diagram for walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers:

Classes

struct  BodyCache
 
struct  ContactCache
 

Public Types

enum  RelaxationModel {
  InelasticFrictionlessContact, ApproximateInelasticCoulombContactByDecoupling, ApproximateInelasticCoulombContactByOrthogonalProjections, InelasticCoulombContactByDecoupling,
  InelasticCoulombContactByOrthogonalProjections, InelasticGeneralizedMaximumDissipationContact
}
 

Public Member Functions

Constructor
 HardContactSemiImplicitTimesteppingSolvers (const shared_ptr< BodyStorage > &globalBodyStorage, const shared_ptr< BlockStorage > &blockStorage, domain_decomposition::BlockDataID storageID, domain_decomposition::BlockDataID ccdID, domain_decomposition::BlockDataID fcdID, WcTimingTree *tt=nullptr)
 Constructor of the CollisionSystem class. More...
 
Destructor
 ~HardContactSemiImplicitTimesteppingSolvers () override
 Destructor of the CollisionSystem class. More...
 
Get functions
real_t getMaximumPenetration () const override
 Returns the maximum penetration depth found in the last collision detection. More...
 
size_t getNumberOfContacts () const override
 Returns the number of contacts found by the last collision detection. More...
 
size_t getNumberOfContactsTreated () const override
 
const std::map< IBlockID::IDType, ContactCachegetContactCache () const
 
real_t getSpeedLimitFactor () const
 
size_t getMaxIterations () const
 
real_t getOverRelaxationParameter () const
 
real_t getRelaxationParameter () const
 
real_t getErrorReductionParameter () const
 
RelaxationModel getRelaxationModel () const
 
Set functions
void setOverRelaxationParameter (real_t omega)
 Sets the relaxation parameter for boundary bodies. More...
 
void setRelaxationParameter (real_t f)
 Sets the relaxation parameter for boundary bodies. More...
 
void setMaxIterations (size_t n)
 Sets the maximum number of iterations performed by the iterative solver. More...
 
void setRelaxationModel (RelaxationModel relaxationModel)
 Sets the relaxation model used by the iterative solver. More...
 
void setErrorReductionParameter (real_t erp)
 Sets the error reduction parameter. More...
 
void setAbortThreshold (real_t threshold)
 Sets the threshold of movement of a particle during collision resolution. More...
 
void setSpeedLimiter (bool active, const real_t speedLimitFactor=real_t(0.0))
 Activates/Deactivates the speed limiter and sets the limit. More...
 
Query functions
bool isSyncRequired () const
 Returns if a synchronization is required before the next time step. More...
 
bool isSyncRequiredLocally () const
 Returns if a synchronization is required by the local process before the next time step. More...
 
bool isSpeedLimiterActive () const
 Returns if speed limiter is currently active and working. More...
 
void operator() (const real_t dt)
 forwards to timestep Convenience operator to make class a functor. More...
 
void timestep (const real_t dt) override
 Advances the simulation dt seconds. More...
 
- Public Member Functions inherited from walberla::pe::cr::ICR
 ICR ()
 
virtual ~ICR ()=default
 
void setGlobalLinearAcceleration (const Vec3 &acc)
 Sets the global linear acceleration. More...
 
const Vec3getGlobalLinearAcceleration () const
 

Private Member Functions

Simulation functions
void resolveContacts (const Contacts &contacts, real_t dt)
 
real_t relaxInelasticFrictionlessContacts (real_t dtinv, HardContactSemiImplicitTimesteppingSolvers::ContactCache &contactCache, HardContactSemiImplicitTimesteppingSolvers::BodyCache &bodyCache)
 Relaxes all contacts once. More...
 
real_t relaxApproximateInelasticCoulombContactsByDecoupling (real_t dtinv, HardContactSemiImplicitTimesteppingSolvers::ContactCache &contactCache, HardContactSemiImplicitTimesteppingSolvers::BodyCache &bodyCache)
 Relaxes all contacts once. More...
 
real_t relaxInelasticCoulombContactsByDecoupling (real_t dtinv, HardContactSemiImplicitTimesteppingSolvers::ContactCache &contactCache, HardContactSemiImplicitTimesteppingSolvers::BodyCache &bodyCache)
 Relaxes all contacts once. More...
 
real_t relaxInelasticCoulombContactsByOrthogonalProjections (real_t dtinv, bool approximate, HardContactSemiImplicitTimesteppingSolvers::ContactCache &contactCache, HardContactSemiImplicitTimesteppingSolvers::BodyCache &bodyCache)
 Relaxes all contacts once. More...
 
real_t relaxInelasticGeneralizedMaximumDissipationContacts (real_t dtinv, HardContactSemiImplicitTimesteppingSolvers::ContactCache &contactCache, HardContactSemiImplicitTimesteppingSolvers::BodyCache &bodyCache)
 Relaxes all contacts once. More...
 
Communication functions
void parseVelocityCorrection (mpi::RecvBuffer &rb, BodyStorage &bodyStorage, BodyCache &bodyCache)
 
void parseVelocityCorrectionShadow (mpi::RecvBuffer &rb, BodyStorage &bodyStorage, BodyCache &bodyCache)
 
void synchronizeVelocities ()
 TODO. More...
 
Time-integration functions
void initializeVelocityCorrections (BodyID body, Vec3 &dv, Vec3 &dw, real_t dt) const
 Calculates the initial velocity corrections of a given body. More...
 
void integratePositions (BodyID body, Vec3 v, Vec3 w, real_t dt) const
 Time integration of the position and orientation of a given body. More...
 
Utility functions
bool checkUpdateFlags ()
 
- Private Member Functions inherited from walberla::NonCopyable
 NonCopyable ()=default
 
 ~NonCopyable ()=default
 

Private Attributes

std::map< IBlockID::IDType, BodyCacheblockToBodyCache_
 
std::map< IBlockID::IDType, ContactCacheblockToContactCache_
 
Member variables
shared_ptr< BodyStorageglobalBodyStorage_
 
shared_ptr< BlockStorageblockStorage_
 
domain_decomposition::BlockDataID storageID_
 
domain_decomposition::BlockDataID ccdID_
 
domain_decomposition::BlockDataID fcdID_
 
WcTimingTreett_
 
real_t erp_
 The error reduction parameter (0 <= erp_ <= 1). More...
 
size_t maxIterations_
 Maximum number of iterations. More...
 
size_t iteration_
 
size_t maxSubIterations_
 Maximum number of iterations of iterative solvers in the one-contact problem. More...
 
real_t abortThreshold_
 If L-infinity iterate difference drops below this threshold the iteration is aborted. More...
 
RelaxationModel relaxationModel_
 The method used to relax unilateral contacts. More...
 
real_t overRelaxationParam_
 Parameter specifying the convergence speed for orthogonal projection models. More...
 
real_t relaxationParam_
 Parameter specifying underrelaxation of velocity corrections for boundary bodies. More...
 
real_t maximumPenetration_
 
size_t numContacts_
 
size_t numContactsTreated_
 
bool speedLimiterActive_
 is the speed limiter active? More...
 
real_t speedLimitFactor_
 what multiple of boundingbox edge length is the body allowed to travel in one timestep More...
 
bool requireSync_
 Flag indicating whether this process requires a synchronization prior to the next time step. More...
 

Member Enumeration Documentation

◆ RelaxationModel

Enumerator
InelasticFrictionlessContact 
ApproximateInelasticCoulombContactByDecoupling 
ApproximateInelasticCoulombContactByOrthogonalProjections 
InelasticCoulombContactByDecoupling 
InelasticCoulombContactByOrthogonalProjections 
InelasticGeneralizedMaximumDissipationContact 

Constructor & Destructor Documentation

◆ HardContactSemiImplicitTimesteppingSolvers()

walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::HardContactSemiImplicitTimesteppingSolvers ( const shared_ptr< BodyStorage > &  globalBodyStorage,
const shared_ptr< BlockStorage > &  blockStorage,
domain_decomposition::BlockDataID  storageID,
domain_decomposition::BlockDataID  ccdID,
domain_decomposition::BlockDataID  fcdID,
WcTimingTree tt = nullptr 
)
inlineexplicit

Constructor of the CollisionSystem class.

◆ ~HardContactSemiImplicitTimesteppingSolvers()

walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::~HardContactSemiImplicitTimesteppingSolvers ( )
inlineoverride

Destructor of the CollisionSystem class.

Member Function Documentation

◆ checkUpdateFlags()

bool walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::checkUpdateFlags ( )
private

◆ getContactCache()

const std::map<IBlockID::IDType, ContactCache> walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::getContactCache ( ) const
inline

◆ getErrorReductionParameter()

real_t walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::getErrorReductionParameter ( ) const
inline

◆ getMaximumPenetration()

real_t walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::getMaximumPenetration ( ) const
inlineoverridevirtual

Returns the maximum penetration depth found in the last collision detection.

Returns
The maximum penetration depth found in the last collision detection.

Only contacts treated on the local process are considered.

Reimplemented from walberla::pe::cr::ICR.

◆ getMaxIterations()

size_t walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::getMaxIterations ( ) const
inline

◆ getNumberOfContacts()

size_t walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::getNumberOfContacts ( ) const
inlineoverridevirtual

Returns the number of contacts found by the last collision detection.

Returns
The number of contacts found by the last collision detection.

Only contacts treated on the local process are counted.

Reimplemented from walberla::pe::cr::ICR.

◆ getNumberOfContactsTreated()

size_t walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::getNumberOfContactsTreated ( ) const
inlineoverridevirtual

Reimplemented from walberla::pe::cr::ICR.

◆ getOverRelaxationParameter()

real_t walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::getOverRelaxationParameter ( ) const
inline

◆ getRelaxationModel()

RelaxationModel walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::getRelaxationModel ( ) const
inline

◆ getRelaxationParameter()

real_t walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::getRelaxationParameter ( ) const
inline

◆ getSpeedLimitFactor()

real_t walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::getSpeedLimitFactor ( ) const
inline

◆ initializeVelocityCorrections()

void walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::initializeVelocityCorrections ( BodyID  body,
Vec3 dv,
Vec3 dw,
real_t  dt 
) const
inlineprivate

Calculates the initial velocity corrections of a given body.

Parameters
bodyThe body whose velocities to time integrate
dvOn return the initial linear velocity correction.
dwOn return the initial angular velocity correction.
dtThe time step size.
Returns
void

Calculates the velocity corrections effected by external forces and torques in an explicit Euler time integration of the velocities of the given body. For fixed objects the velocity corrections are set to zero. External forces and torques are reset if indicated by the settings.

◆ integratePositions()

void walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::integratePositions ( BodyID  body,
Vec3  v,
Vec3  w,
real_t  dt 
) const
inlineprivate

Time integration of the position and orientation of a given body.

Parameters
bodyThe body whose position and orientation to time integrate
vThe linear velocity to use for time integration of the position.
wThe angular velocity to use for time integration of the orientation.
dtThe time step size.
Returns
void

Performs an Euler time integration of the positions of the given body. Velocities are damped if indicated by the settings and stored back in the body properties. The bounding box is recalculated and it is redetermined whether the body is awake or not. Also the data structure tracking the contacts attached to the body are cleared and

◆ isSpeedLimiterActive()

bool walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::isSpeedLimiterActive ( ) const
inline

Returns if speed limiter is currently active and working.

Returns
status of the speed limiter

◆ isSyncRequired()

bool walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::isSyncRequired ( ) const
inline

Returns if a synchronization is required before the next time step.

Returns
True if a synchronization is strictly required, false if a synchronization is possibly required for correct dynamics but not enforced.

Insertion or removal of bodies from the simulation can require a subsequent synchronization before performing the next time step. If e.g. no other process obtained a shadow copy of a body to be removed then a synchronization is not enforced. However, if a neighbor has a shadow copy a synchronization is required. Changing e.g. velocities or positions can lead to inconsistent descriptions of bodies across process boundaries but synchronization is not enforced. In this case it is the users obligation to synchronize whenever necessary.

WARNING: This query function uses an expensive allreduce MPI operation to determine the result!

◆ isSyncRequiredLocally()

bool walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::isSyncRequiredLocally ( ) const
inline

Returns if a synchronization is required by the local process before the next time step.

Returns
True if a synchronization is strictly required by the local process, false otherwise.

Insertion or removal of bodies from the simulation can require a subsequent synchronization before performing the next time step. If e.g. no other process obtained a shadow copy of a body to be removed then a synchronization is not enforced. However, if a neighbor has a shadow copy a synchronization is required. Changing e.g. velocities or positions can lead to inconsistent descriptions of bodies across process boundaries but synchronization is not enforced. In this case it is the users obligation to synchronize whenever necessary.

◆ operator()()

void walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::operator() ( const real_t  dt)
inline

forwards to timestep Convenience operator to make class a functor.

◆ parseVelocityCorrection()

void walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::parseVelocityCorrection ( mpi::RecvBuffer rb,
BodyStorage bodyStorage,
BodyCache bodyCache 
)
inlineprivate

◆ parseVelocityCorrectionShadow()

void walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::parseVelocityCorrectionShadow ( mpi::RecvBuffer rb,
BodyStorage bodyStorage,
BodyCache bodyCache 
)
inlineprivate

◆ relaxApproximateInelasticCoulombContactsByDecoupling()

real_t walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::relaxApproximateInelasticCoulombContactsByDecoupling ( real_t  dtinv,
HardContactSemiImplicitTimesteppingSolvers::ContactCache contactCache,
HardContactSemiImplicitTimesteppingSolvers::BodyCache bodyCache 
)
inlineprivate

Relaxes all contacts once.

The contact model is for inelastic unilateral contacts with approximate Coulomb friction.

Returns
The largest variation of contact impulses in the L-infinity norm.

This function is to be called from resolveContacts(). Separating contacts are preferred over other solutions if valid. Static solutions are preferred over dynamic solutions. Dynamic solutions are computed by decoupling the normal from the frictional components. That is for a dynamic contact the normal component is relaxed first followed by the frictional components. The determination of the frictional components does not perform any subiterations and guarantees that the friction partially opposes slip.

◆ relaxInelasticCoulombContactsByDecoupling()

real_t walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::relaxInelasticCoulombContactsByDecoupling ( real_t  dtinv,
HardContactSemiImplicitTimesteppingSolvers::ContactCache contactCache,
HardContactSemiImplicitTimesteppingSolvers::BodyCache bodyCache 
)
inlineprivate

Relaxes all contacts once.

The contact model is for inelastic unilateral contacts with Coulomb friction.

Returns
The largest variation of contact impulses in the L-infinity norm.

This function is to be called from resolveContacts(). Separating contacts are preferred over other solutions if valid. Static solutions are preferred over dynamic solutions. Dynamic solutions are computed by decoupling the normal from the frictional components. That is for a dynamic contact the normal component is relaxed first followed by the frictional components. How much the frictional components directly oppose slip as required by the Coulomb friction model depends on the number of subiterations performed. If no subiterations are performed the friction is guaranteed to be at least partially dissipative.

◆ relaxInelasticCoulombContactsByOrthogonalProjections()

real_t walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::relaxInelasticCoulombContactsByOrthogonalProjections ( real_t  dtinv,
bool  approximate,
HardContactSemiImplicitTimesteppingSolvers::ContactCache contactCache,
HardContactSemiImplicitTimesteppingSolvers::BodyCache bodyCache 
)
inlineprivate

Relaxes all contacts once.

The contact model is for inelastic unilateral contacts with Coulomb friction.

Parameters
dtinvThe inverse of the current time step.
approximateUse the approximate model showing bouncing.
Returns
The largest variation of contact impulses in the L-infinity norm.

This function is to be called from resolveContacts(). The iterative method to solve the contact problem is e.g. described in the article "A matrix-free cone complementarity approach for solving large-scale, nonsmooth, rigid body dynamics" by A. Tasora and M. Anitescu in Computer Methods in Applied Mechanics and Engineering (Volume 200, Issues 5–8, 15 January 2011, Pages 439-453). The contact model is purely quadratic and convergence should be good but depends on a parameter. The one-contact problem has a unique solution. The frictional reactions for a dynamic contact converge to those that directly oppose slip. However, the contact is not perfectly inelastic for dynamic contacts but bounces. These vertical motions tend to go to zero for smaller time steps and can be interpreted as exaggerated vertical motions coming from micro asperities (see "Optimization-based simulation of nonsmooth rigid multibody dynamics" by M. Anitescu in Mathematical Programming (Volume 105, Issue 1, January 2006, Pages 113-143). These motions can be prevented by a small change in the iteration proposed in "The bipotential method: a constructive approach to design the complete contact law with friction and improved numerical algorithms" by G. De Saxce and Z-Q. Feng in Mathematical and Computer Modelling (Volume 28, Issue 4, 1998, Pages 225-245). Which iteration is used is controlled with the approximate parameter.

◆ relaxInelasticFrictionlessContacts()

real_t walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::relaxInelasticFrictionlessContacts ( real_t  dtinv,
HardContactSemiImplicitTimesteppingSolvers::ContactCache contactCache,
HardContactSemiImplicitTimesteppingSolvers::BodyCache bodyCache 
)
inlineprivate

Relaxes all contacts once.

The contact model is for inelastic unilateral contacts without friction.

Returns
The largest variation of contact impulses in the L-infinity norm.

This function is to be called from resolveContacts(). Separating contacts are preferred over persisting solutions if valid.

◆ relaxInelasticGeneralizedMaximumDissipationContacts()

real_t walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::relaxInelasticGeneralizedMaximumDissipationContacts ( real_t  dtinv,
HardContactSemiImplicitTimesteppingSolvers::ContactCache contactCache,
HardContactSemiImplicitTimesteppingSolvers::BodyCache bodyCache 
)
inlineprivate

Relaxes all contacts once.

The contact model is for inelastic unilateral contacts with the generalized maximum dissipation principle for friction.

Returns
The largest variation of contact impulses in the L-infinity norm.

This function is to be called from resolveContacts(). Dynamic solutions are computed by minimizing the kinetic energy along the intersection of the plane of maximum compression and the friction cone.

◆ resolveContacts()

void walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::resolveContacts ( const Contacts contacts,
real_t  dt 
)
private

◆ setAbortThreshold()

void walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::setAbortThreshold ( real_t  threshold)
inline

Sets the threshold of movement of a particle during collision resolution.

Parameters
thresholdIf movement is smaller than threshold, col. resolution is stopped.
Returns
void

◆ setErrorReductionParameter()

void walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::setErrorReductionParameter ( real_t  erp)
inline

Sets the error reduction parameter.

Parameters
erpThe error reduction parameter.
Returns
void

If body shapes overlap by x at a contact then the contact resolution aims to remove erp*x of the overlap. Thus the error reduction parameter must be between 0 and 1. 0 corresponds to no error reduction and is the default. 1 corresponds to full error reduction. Note that error reduction (constraint stabilization) introduces additional energy to the system.

◆ setMaxIterations()

void walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::setMaxIterations ( size_t  n)
inline

Sets the maximum number of iterations performed by the iterative solver.

Parameters
nThe maximum number of iterations.
Returns
void

◆ setOverRelaxationParameter()

void walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::setOverRelaxationParameter ( real_t  omega)
inline

Sets the relaxation parameter for boundary bodies.

Parameters
omegaThe overrelaxation parameter.
Returns
void

The overrelaxation parameter omega is only used when the relaxation model is one of

  • ApproximateInelasticCoulombContactByOrthogonalProjections
  • InelasticCoulombContactByOrthogonalProjections

It is used to control the convergence of the model. Large values show faster convergence, but they can also lead to divergence ("exploding" particles). The default value is 1.0.

◆ setRelaxationModel()

void walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::setRelaxationModel ( RelaxationModel  relaxationModel)
inline

Sets the relaxation model used by the iterative solver.

Parameters
relaxationModelThe relaxation model to be used by the iterative solver.
Returns
void

◆ setRelaxationParameter()

void walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::setRelaxationParameter ( real_t  f)
inline

Sets the relaxation parameter for boundary bodies.

Parameters
fThe relaxation parameter.
Returns
void

The iterative solvers are a mixture of non-linear Gauss-Seidel and Jacobi solvers. This might require underrelaxation. The parameter must be positive. Note that for dilute systems the solver might need stronger underrelaxation (smaller f) than for dense systems.

◆ setSpeedLimiter()

void walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::setSpeedLimiter ( bool  active,
const real_t  speedLimitFactor = real_t(0.0) 
)
inline

Activates/Deactivates the speed limiter and sets the limit.

Parameters
activeactivate/deactivate speed limiter
speedLimitFactorsize of bounding box will be multiplied by this factor to get the maximal distance a body is allowed to travel within one timestep
Returns
void

◆ synchronizeVelocities()

void walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::synchronizeVelocities ( )
inlineprivate

TODO.

Returns
void

◆ timestep()

void walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::timestep ( const real_t  dt)
inlineoverridevirtual

Advances the simulation dt seconds.

Runs the simulation for one timestep.

Parameters
dtSize of the time step.
Returns
void

Implements walberla::pe::cr::ICR.

Member Data Documentation

◆ abortThreshold_

real_t walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::abortThreshold_
private

If L-infinity iterate difference drops below this threshold the iteration is aborted.

◆ blockStorage_

shared_ptr<BlockStorage> walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::blockStorage_
private

◆ blockToBodyCache_

std::map<IBlockID::IDType, BodyCache> walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::blockToBodyCache_
private

◆ blockToContactCache_

std::map<IBlockID::IDType, ContactCache> walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::blockToContactCache_
private

◆ ccdID_

domain_decomposition::BlockDataID walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::ccdID_
private

◆ erp_

real_t walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::erp_
private

The error reduction parameter (0 <= erp_ <= 1).

◆ fcdID_

domain_decomposition::BlockDataID walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::fcdID_
private

◆ globalBodyStorage_

shared_ptr<BodyStorage> walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::globalBodyStorage_
private

◆ iteration_

size_t walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::iteration_
private

◆ maximumPenetration_

real_t walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::maximumPenetration_
private

◆ maxIterations_

size_t walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::maxIterations_
private

Maximum number of iterations.

◆ maxSubIterations_

size_t walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::maxSubIterations_
private

Maximum number of iterations of iterative solvers in the one-contact problem.

◆ numContacts_

size_t walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::numContacts_
private

◆ numContactsTreated_

size_t walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::numContactsTreated_
private

◆ overRelaxationParam_

real_t walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::overRelaxationParam_
private

Parameter specifying the convergence speed for orthogonal projection models.

◆ relaxationModel_

RelaxationModel walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::relaxationModel_
private

The method used to relax unilateral contacts.

◆ relaxationParam_

real_t walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::relaxationParam_
private

Parameter specifying underrelaxation of velocity corrections for boundary bodies.

◆ requireSync_

bool walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::requireSync_
private

Flag indicating whether this process requires a synchronization prior to the next time step.

◆ speedLimiterActive_

bool walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::speedLimiterActive_
private

is the speed limiter active?

◆ speedLimitFactor_

real_t walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::speedLimitFactor_
private

what multiple of boundingbox edge length is the body allowed to travel in one timestep

◆ storageID_

domain_decomposition::BlockDataID walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::storageID_
private

◆ tt_

WcTimingTree* walberla::pe::cr::HardContactSemiImplicitTimesteppingSolvers::tt_
private

The documentation for this class was generated from the following files:
@ InelasticGeneralizedMaximumDissipationContact
Definition: HCSITS.h:108
math::Vector3< real_t > Vec3
Definition: Types.h:46
HardContactSemiImplicitTimesteppingSolvers HCSITS
Definition: HCSITS.h:273
float real_t
Definition: DataTypes.h:167