IBAMR
An adaptive and distributed-memory parallel implementation of the immersed boundary (IB) method
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
IBAMR::FEMechanicsExplicitIntegrator Class Reference

Class FEMechanicsExplicitIntegrator is an implementation of the abstract base class FEMechanicsBase that provides a simple explicit elastodynamics time integrator with an interface that is similar to IBFEMethod to facilitate model re-use. More...

#include </home/runner/work/IBAMR/IBAMR/include/ibamr/FEMechanicsExplicitIntegrator.h>

Inheritance diagram for IBAMR::FEMechanicsExplicitIntegrator:
Inheritance graph
[legend]

Public Member Functions

 FEMechanicsExplicitIntegrator (const std::string &object_name, const SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > &input_db, libMesh::MeshBase *mesh, bool register_for_restart=true, const std::string &restart_read_dirname="", unsigned int restart_restore_number=0)
 Constructor for a single-part model.
 
 FEMechanicsExplicitIntegrator (const std::string &object_name, const SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > &input_db, const std::vector< libMesh::MeshBase * > &meshes, bool register_for_restart=true, const std::string &restart_read_dirname="", unsigned int restart_restore_number=0)
 Constructor for a multi-part model.
 
 FEMechanicsExplicitIntegrator ()=delete
 Deleted default constructor.
 
 FEMechanicsExplicitIntegrator (const FEMechanicsExplicitIntegrator &from)=delete
 Deleted copy constructor.
 
FEMechanicsExplicitIntegratoroperator= (const FEMechanicsExplicitIntegrator &that)=delete
 Deleted assignment operator.
 
 ~FEMechanicsExplicitIntegrator () override=default
 Defaulted destructor.
 
void preprocessIntegrateData (double current_time, double new_time, int num_cycles) override
 
void postprocessIntegrateData (double current_time, double new_time, int num_cycles) override
 
void forwardEulerStep (double current_time, double new_time)
 
void modifiedEulerStep (double current_time, double new_time)
 
void backwardEulerStep (double current_time, double new_time)
 
void midpointStep (double current_time, double new_time)
 
void trapezoidalStep (double current_time, double new_time)
 
void modifiedTrapezoidalStep (double current_time, double new_time)
 
void SSPRK3Step (double current_time, double new_time, unsigned int n_stages=4)
 
void computeLagrangianForce (double data_time)
 
void computeLagrangianForce (libMesh::PetscVector< double > &F_vec, libMesh::PetscVector< double > &X_vec, libMesh::PetscVector< double > &U_vec, libMesh::PetscVector< double > *P_vec, double data_time, unsigned int part)
 
void putToDatabase (SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > db) override
 

Protected Member Functions

void doInitializeFEEquationSystems () override
 
void doInitializeFESystemVectors () override
 
void doInitializeFEData (bool use_present_data) override
 
void doForwardEulerStep (libMesh::PetscVector< double > &X_new_vec, libMesh::PetscVector< double > &U_new_vec, libMesh::PetscVector< double > *P_new_vec, libMesh::PetscVector< double > &X_current_vec, libMesh::PetscVector< double > &U_current_vec, libMesh::PetscVector< double > *P_current_vec, libMesh::PetscVector< double > &F_current_vec, libMesh::PetscVector< double > *dP_dt_current_vec, double current_time, double new_time, unsigned int part)
 

Protected Attributes

std::vector< double > d_rhos
 Structure mass densities.
 
std::vector< double > d_etas
 Structure damping coefficients.
 

Detailed Description

Class FEMechanicsExplicitIntegrator is an implementation of the abstract base class FEMechanicsBase that provides a simple explicit elastodynamics time integrator with an interface that is similar to IBFEMethod to facilitate model re-use.

See also
IBFEMethod

Member Function Documentation

◆ backwardEulerStep()

void IBAMR::FEMechanicsExplicitIntegrator::backwardEulerStep ( double  current_time,
double  new_time 
)

Advance the structural velocities and positions using backward Euler.

NOTE: This is not a true backward Euler method; we simply use the most recent approximations to the updated velocity and position to compute the force at the end of the time step, and then use that force to increment the velocity and position via backward Euler.

◆ computeLagrangianForce() [1/2]

void IBAMR::FEMechanicsExplicitIntegrator::computeLagrangianForce ( double  data_time)

Compute the Lagrangian force at the specified time within the current time interval.

◆ computeLagrangianForce() [2/2]

void IBAMR::FEMechanicsExplicitIntegrator::computeLagrangianForce ( libMesh::PetscVector< double > &  F_vec,
libMesh::PetscVector< double > &  X_vec,
libMesh::PetscVector< double > &  U_vec,
libMesh::PetscVector< double > *  P_vec,
double  data_time,
unsigned int  part 
)

Compute the Lagrangian force at the specified time within the current time interval using the provided data vectors.

◆ doForwardEulerStep()

void IBAMR::FEMechanicsExplicitIntegrator::doForwardEulerStep ( libMesh::PetscVector< double > &  X_new_vec,
libMesh::PetscVector< double > &  U_new_vec,
libMesh::PetscVector< double > *  P_new_vec,
libMesh::PetscVector< double > &  X_current_vec,
libMesh::PetscVector< double > &  U_current_vec,
libMesh::PetscVector< double > *  P_current_vec,
libMesh::PetscVector< double > &  F_current_vec,
libMesh::PetscVector< double > *  dP_dt_current_vec,
double  current_time,
double  new_time,
unsigned int  part 
)
protected

Perform a forward Euler step.

◆ doInitializeFEData()

void IBAMR::FEMechanicsExplicitIntegrator::doInitializeFEData ( bool  use_present_data)
overrideprotected

Do the actual work in reinitializeFEData and initializeFEData. if use_present_data is true then the current content of the solution vectors is used: more exactly, the coordinates and velocities (computed by initializeCoordinates and initializeVelocity) are considered as being up to date, as is the direct forcing kinematic data.

◆ doInitializeFEEquationSystems()

void IBAMR::FEMechanicsExplicitIntegrator::doInitializeFEEquationSystems ( )
overrideprotected

Do the actual work in initializeFEEquationSystems.

◆ doInitializeFESystemVectors()

void IBAMR::FEMechanicsExplicitIntegrator::doInitializeFESystemVectors ( )
overrideprotected

Do the actual work of setting up libMesh system vectors.

◆ forwardEulerStep()

void IBAMR::FEMechanicsExplicitIntegrator::forwardEulerStep ( double  current_time,
double  new_time 
)

Advance the structural velocities and positions using forward Euler.

◆ midpointStep()

void IBAMR::FEMechanicsExplicitIntegrator::midpointStep ( double  current_time,
double  new_time 
)

Advance the structural velocities and positions using the explicit midpoint rule.

◆ modifiedEulerStep()

void IBAMR::FEMechanicsExplicitIntegrator::modifiedEulerStep ( double  current_time,
double  new_time 
)

Advance the structural velocities and positions using "modified" Euler.

NOTE: This scheme first advances the velocity using forward Euler, and then uses the updated velocity to advance the position of the structure.

◆ modifiedTrapezoidalStep()

void IBAMR::FEMechanicsExplicitIntegrator::modifiedTrapezoidalStep ( double  current_time,
double  new_time 
)

Advance the structural velocities and positions using a modified explicit trapezoidal rule.

NOTE: In this scheme, the first stage uses a modified Euler scheme (

See also
modifiedEulerStep), and the second stage uses a modified explicit trapezoidal rule, in which the final approximation to the updated velocity is used instead of the intermediate updated approximation used in the standard SSP RK2 scheme.

◆ postprocessIntegrateData()

void IBAMR::FEMechanicsExplicitIntegrator::postprocessIntegrateData ( double  current_time,
double  new_time,
int  num_cycles 
)
override

Method to clean up data following call(s) to integrateHierarchy().

◆ preprocessIntegrateData()

void IBAMR::FEMechanicsExplicitIntegrator::preprocessIntegrateData ( double  current_time,
double  new_time,
int  num_cycles 
)
override

Method to prepare to advance data from current_time to new_time.

◆ putToDatabase()

void IBAMR::FEMechanicsExplicitIntegrator::putToDatabase ( SAMRAI::tbox::Pointer< SAMRAI::tbox::Database db)
override

Write out object state to the given database.

◆ SSPRK3Step()

void IBAMR::FEMechanicsExplicitIntegrator::SSPRK3Step ( double  current_time,
double  new_time,
unsigned int  n_stages = 4 
)

Advance the structural velocities and positions using the explicit SSP RK3 scheme.

The implementation supports either 3- or 4-stage variants. The default is the 4-stage version, which has a larger stability region.

◆ trapezoidalStep()

void IBAMR::FEMechanicsExplicitIntegrator::trapezoidalStep ( double  current_time,
double  new_time 
)

Advance the structural velocities and positions using the explicit SSP RK2 scheme, which is the explicit trapezoidal rule.


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