IBAMR  IBAMR version 0.19.
Classes | Public Types | Public Member Functions | Static Public Attributes | Protected Member Functions | Static Protected Member Functions | Protected Attributes | Private Member Functions | List of all members
IBAMR::FEMechanicsBase Class Reference

Class FEMechanicsBase provides core finite element mechanics functionality and data management. More...

#include <ibamr/FEMechanicsBase.h>

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

Classes

struct  CoordinateMappingFcnData
 
struct  InitialVelocityFcnData
 
struct  LagBodyForceFcnData
 
struct  LagSurfaceForceFcnData
 
struct  LagSurfacePressureFcnData
 
struct  PK1StressFcnData
 

Public Types

using CoordinateMappingFcnPtr = void(*)(libMesh::Point &x, const libMesh::Point &X, void *ctx)
 
using InitialVelocityFcnPtr = void(*)(libMesh::VectorValue< double > &U0, const libMesh::Point &X0, void *ctx)
 
using PK1StressFcnPtr = IBTK::TensorMeshFcnPtr
 
using LagBodyForceFcnPtr = IBTK::VectorMeshFcnPtr
 
using LagSurfacePressureFcnPtr = IBTK::ScalarSurfaceFcnPtr
 
using LagSurfaceForceFcnPtr = IBTK::VectorSurfaceFcnPtr
 
using VolumetricEnergyDerivativeFcn = double(*)(double)
 

Public Member Functions

 FEMechanicsBase (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)
 
 FEMechanicsBase (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)
 
 FEMechanicsBase ()=delete
 
 FEMechanicsBase (const FEMechanicsBase &from)=delete
 
FEMechanicsBaseoperator= (const FEMechanicsBase &that)=delete
 
virtual ~FEMechanicsBase () override
 
libMesh::EquationSystems * getEquationSystems (unsigned int part=0) const
 
const std::string & getCurrentCoordinatesSystemName () const
 
const std::string & getDisplacementSystemName () const
 
const std::string & getForceSystemName () const
 
const std::string & getPressureSystemName () const
 
const std::string & getVelocitySystemName () const
 
std::shared_ptr< IBTK::FEDatagetFEData (unsigned int part=0) const
 
virtual void registerInitialCoordinateMappingFunction (const CoordinateMappingFcnData &data, unsigned int part=0)
 
CoordinateMappingFcnData getInitialCoordinateMappingFunction (unsigned int part=0) const
 
virtual void registerInitialVelocityFunction (const InitialVelocityFcnData &data, unsigned int part=0)
 
InitialVelocityFcnData getInitialVelocityFunction (unsigned int part=0) const
 
virtual void registerPK1StressFunction (const PK1StressFcnData &data, unsigned int part=0)
 
std::vector< PK1StressFcnDatagetPK1StressFunction (unsigned int part=0) const
 
virtual void registerLagBodyForceFunction (const LagBodyForceFcnData &data, unsigned int part=0)
 
LagBodyForceFcnData getLagBodyForceFunction (unsigned int part=0) const
 
virtual void registerLagSurfacePressureFunction (const LagSurfacePressureFcnData &data, unsigned int part=0)
 
LagSurfacePressureFcnData getLagSurfacePressureFunction (unsigned int part=0) const
 
virtual void registerLagSurfaceForceFunction (const LagSurfaceForceFcnData &data, unsigned int part=0)
 
LagSurfaceForceFcnData getLagSurfaceForceFunction (unsigned int part=0) const
 
virtual void registerStaticPressurePart (PressureProjectionType projection_type=CONSISTENT_PROJECTION, VolumetricEnergyDerivativeFcn dU_dJ_fcn=nullptr, unsigned int part=0)
 
virtual void registerDynamicPressurePart (PressureProjectionType projection_type=CONSISTENT_PROJECTION, VolumetricEnergyDerivativeFcn d2U_dJ2_fcn=nullptr, unsigned int part=0)
 
bool partHasPressure (unsigned int part=0)
 
virtual void preprocessIntegrateData (double current_time, double new_time, int num_cycles)
 
virtual void postprocessIntegrateData (double current_time, double new_time, int num_cycles)
 
virtual void initializeFEEquationSystems ()
 
virtual void initializeFEData ()
 
virtual void reinitializeFEData ()
 
virtual void putToDatabase (SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > db) override
 
virtual void writeFEDataToRestartFile (const std::string &restart_dump_dirname, unsigned int time_step_number)
 

Static Public Attributes

static const std::string COORDS_SYSTEM_NAME
 
static const std::string COORD_MAPPING_SYSTEM_NAME
 
static const std::string FORCE_SYSTEM_NAME
 
static const std::string PRESSURE_SYSTEM_NAME
 
static const std::string VELOCITY_SYSTEM_NAME
 

Protected Member Functions

virtual void doInitializeFEEquationSystems ()
 
virtual void doInitializeFESystemVectors ()
 
virtual void doInitializeFEData (bool use_present_data)
 
void computeStaticPressure (libMesh::PetscVector< double > &P_vec, libMesh::PetscVector< double > &X_vec, double data_time, unsigned int part)
 Compute the static pressure field. More...
 
void computeDynamicPressureRateOfChange (libMesh::PetscVector< double > &dP_dt_vec, libMesh::PetscVector< double > &X_vec, libMesh::PetscVector< double > &U_vec, double data_time, unsigned int part)
 Compute the dynamic pressure time rate of change. More...
 
virtual void assembleInteriorForceDensityRHS (libMesh::PetscVector< double > &F_rhs_vec, libMesh::PetscVector< double > &X_vec, libMesh::PetscVector< double > *P_vec, double data_time, unsigned int part)
 
virtual void initializeCoordinates (unsigned int part)
 
virtual void updateCoordinateMapping (unsigned int part)
 
virtual void initializeVelocity (unsigned int part)
 
virtual std::string getLibMeshRestartFileName (const std::string &restart_dump_dirname, unsigned int time_step_number, unsigned int part, const std::string &extension) const
 

Static Protected Member Functions

static void setup_system_vectors (libMesh::EquationSystems *equation_systems, const std::vector< std::string > &system_names, const std::vector< std::string > &vector_names)
 
static void setup_system_vector (libMesh::System &system, const std::string &vector_name)
 

Protected Attributes

SAMRAI::tbox::Pointer< SAMRAI::tbox::Databased_fe_projector_db
 
bool d_do_log = false
 
double d_current_time = std::numeric_limits<double>::quiet_NaN()
 
double d_new_time = std::numeric_limits<double>::quiet_NaN()
 
double d_half_time = std::numeric_limits<double>::quiet_NaN()
 
std::vector< libMesh::MeshBase * > d_meshes
 
libMesh::CouplingMatrix d_diagonal_system_coupling
 
std::vector< std::unique_ptr< libMesh::EquationSystems > > d_equation_systems
 
const std::string d_current_coordinates_system_name
 System names for key variables. More...
 
const std::string d_displacement_system_name
 
const std::string d_force_system_name
 
const std::string d_pressure_system_name
 
const std::string d_velocity_system_name
 
std::vector< std::shared_ptr< IBTK::FEData > > d_fe_data
 FEData objects provide key FE data management. More...
 
std::vector< std::shared_ptr< IBTK::FEProjector > > d_fe_projectors
 FEProjector objects provide L2 projection functionality. More...
 
std::vector< libMesh::ExplicitSystem * > d_X_systems
 
std::vector< libMesh::ExplicitSystem * > d_U_systems
 
std::vector< libMesh::ExplicitSystem * > d_F_systems
 
std::vector< libMesh::ExplicitSystem * > d_P_systems
 
std::unique_ptr< IBTK::LibMeshSystemVectorsd_X_vecs
 
std::unique_ptr< IBTK::LibMeshSystemVectorsd_U_vecs
 
std::unique_ptr< IBTK::LibMeshSystemVectorsd_F_vecs
 
std::unique_ptr< IBTK::LibMeshSystemVectorsd_P_vecs
 
bool d_fe_equation_systems_initialized = false
 
bool d_fe_data_initialized = false
 
LibmeshPartitionerType d_libmesh_partitioner_type = LIBMESH_DEFAULT
 
bool d_libmesh_use_amr = false
 
std::vector< libMesh::Order > d_fe_order_position
 
std::vector< libMesh::Order > d_fe_order_force
 
std::vector< libMesh::Order > d_fe_order_pressure
 
std::vector< libMesh::FEFamily > d_fe_family_position
 
std::vector< libMesh::FEFamily > d_fe_family_force
 
std::vector< libMesh::FEFamily > d_fe_family_pressure
 
std::vector< libMesh::QuadratureType > d_default_quad_type_stress
 
std::vector< libMesh::QuadratureType > d_default_quad_type_force
 
std::vector< libMesh::QuadratureType > d_default_quad_type_pressure
 
std::vector< libMesh::Order > d_default_quad_order_stress
 
std::vector< libMesh::Order > d_default_quad_order_force
 
std::vector< libMesh::Order > d_default_quad_order_pressure
 
bool d_use_consistent_mass_matrix = true
 
bool d_allow_rules_with_negative_weights = true
 
bool d_include_normal_stress_in_weak_form = false
 
bool d_include_tangential_stress_in_weak_form = false
 
bool d_include_normal_surface_forces_in_weak_form = true
 
bool d_include_tangential_surface_forces_in_weak_form = true
 
std::vector< CoordinateMappingFcnDatad_coordinate_mapping_fcn_data
 
std::vector< InitialVelocityFcnDatad_initial_velocity_fcn_data
 
std::vector< std::vector< PK1StressFcnData > > d_PK1_stress_fcn_data
 
std::vector< LagBodyForceFcnDatad_lag_body_force_fcn_data
 
std::vector< LagSurfacePressureFcnDatad_lag_surface_pressure_fcn_data
 
std::vector< LagSurfaceForceFcnDatad_lag_surface_force_fcn_data
 
double d_static_pressure_kappa = 0.0
 
double d_static_pressure_stab_param = 0.0
 
bool d_has_static_pressure_parts = false
 
std::vector< boold_static_pressure_part
 
std::vector< PressureProjectionTyped_static_pressure_proj_type
 
std::vector< VolumetricEnergyDerivativeFcnd_static_pressure_dU_dJ_fcn
 
double d_dynamic_pressure_kappa = 0.0
 
double d_dynamic_pressure_stab_param = 0.0
 
bool d_has_dynamic_pressure_parts = false
 
std::vector< boold_dynamic_pressure_part
 
std::vector< PressureProjectionTyped_dynamic_pressure_proj_type
 
std::vector< VolumetricEnergyDerivativeFcnd_dynamic_pressure_d2U_dJ2_fcn
 
std::string d_object_name
 
bool d_registered_for_restart
 
std::string d_libmesh_restart_read_dir
 
unsigned int d_libmesh_restart_restore_number
 
std::string d_libmesh_restart_file_extension
 

Private Member Functions

void commonConstructor (const std::string &object_name, const SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > &input_db, const std::vector< libMesh::MeshBase * > &meshes, bool register_for_restart, const std::string &restart_read_dirname, unsigned int restart_restore_number)
 
void getFromInput (const SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > &db, bool is_from_restart)
 
void getFromRestart ()
 

Detailed Description

Parameters read from the input database

  1. FEProjector: Input database passed along to the object responsible for computing projections onto the finite element space. See the documentation of IBTK::FEProjector for more information.

Member Typedef Documentation

◆ CoordinateMappingFcnPtr

using IBAMR::FEMechanicsBase::CoordinateMappingFcnPtr = void (*)(libMesh::Point& x, const libMesh::Point& X, void* ctx)

Typedef specifying interface for coordinate mapping function.

◆ InitialVelocityFcnPtr

using IBAMR::FEMechanicsBase::InitialVelocityFcnPtr = void (*)(libMesh::VectorValue<double>& U0, const libMesh::Point& X0, void* ctx)

Typedef specifying interface for initial velocity specification function.

◆ PK1StressFcnPtr

Typedef specifying interface for PK1 stress tensor function.

◆ LagBodyForceFcnPtr

Typedef specifying interface for Lagrangian body force distribution function.

◆ LagSurfacePressureFcnPtr

Typedef specifying interface for Lagrangian pressure force distribution function.

◆ LagSurfaceForceFcnPtr

Typedef specifying interface for Lagrangian surface force distribution function.

◆ VolumetricEnergyDerivativeFcn

Function signature for specifying the energy functional that determines the pressure.

Constructor & Destructor Documentation

◆ FEMechanicsBase() [1/4]

IBAMR::FEMechanicsBase::FEMechanicsBase ( 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.

◆ FEMechanicsBase() [2/4]

IBAMR::FEMechanicsBase::FEMechanicsBase ( 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.

◆ FEMechanicsBase() [3/4]

IBAMR::FEMechanicsBase::FEMechanicsBase ( )
delete

Deleted default constructor.

◆ FEMechanicsBase() [4/4]

IBAMR::FEMechanicsBase::FEMechanicsBase ( const FEMechanicsBase from)
delete

Deleted copy constructor.

◆ ~FEMechanicsBase()

virtual IBAMR::FEMechanicsBase::~FEMechanicsBase ( )
overridevirtual

Destructor.

Member Function Documentation

◆ operator=()

FEMechanicsBase& IBAMR::FEMechanicsBase::operator= ( const FEMechanicsBase that)
delete

Deleted assignment operator.

◆ getEquationSystems()

libMesh::EquationSystems* IBAMR::FEMechanicsBase::getEquationSystems ( unsigned int  part = 0) const

Return a pointer to the equations system object for the specified part.

◆ getCurrentCoordinatesSystemName()

const std::string& IBAMR::FEMechanicsBase::getCurrentCoordinatesSystemName ( ) const
inline

Return the current coordinates system name.

◆ getDisplacementSystemName()

const std::string& IBAMR::FEMechanicsBase::getDisplacementSystemName ( ) const
inline

Return the displacement system name.

◆ getForceSystemName()

const std::string& IBAMR::FEMechanicsBase::getForceSystemName ( ) const
inline

Return the force system name.

◆ getPressureSystemName()

const std::string& IBAMR::FEMechanicsBase::getPressureSystemName ( ) const
inline

Return the pressure system name.

◆ getVelocitySystemName()

const std::string& IBAMR::FEMechanicsBase::getVelocitySystemName ( ) const
inline

Return the velocity system name.

◆ getFEData()

std::shared_ptr<IBTK::FEData> IBAMR::FEMechanicsBase::getFEData ( unsigned int  part = 0) const

Return a pointer to the FEData object for the specified part.

◆ registerInitialCoordinateMappingFunction()

virtual void IBAMR::FEMechanicsBase::registerInitialCoordinateMappingFunction ( const CoordinateMappingFcnData data,
unsigned int  part = 0 
)
virtual

Register the (optional) function used to initialize the physical coordinates from the Lagrangian coordinates.

Note
If no function is provided, the initial physical coordinates are taken to be the same as the Lagrangian coordinate system, i.e., the initial coordinate mapping is assumed to be the identity mapping.

◆ getInitialCoordinateMappingFunction()

CoordinateMappingFcnData IBAMR::FEMechanicsBase::getInitialCoordinateMappingFunction ( unsigned int  part = 0) const

Get the initial coordinate mapping function data.

◆ registerInitialVelocityFunction()

virtual void IBAMR::FEMechanicsBase::registerInitialVelocityFunction ( const InitialVelocityFcnData data,
unsigned int  part = 0 
)
virtual

Register the (optional) function used to initialize the velocity of the solid mesh.

Note
If no function is provided, the initial velocity is taken to be zero.

◆ getInitialVelocityFunction()

InitialVelocityFcnData IBAMR::FEMechanicsBase::getInitialVelocityFunction ( unsigned int  part = 0) const

Get the initial velocity function data.

◆ registerPK1StressFunction()

virtual void IBAMR::FEMechanicsBase::registerPK1StressFunction ( const PK1StressFcnData data,
unsigned int  part = 0 
)
virtual

Register the (optional) function to compute the first Piola-Kirchhoff stress tensor, used to compute the forces on the Lagrangian finite element mesh.

Note
It is possible to register multiple PK1 stress functions with this class. This is intended to be used to implement selective reduced integration.

◆ getPK1StressFunction()

std::vector<PK1StressFcnData> IBAMR::FEMechanicsBase::getPK1StressFunction ( unsigned int  part = 0) const

Get the PK1 stress function data.

◆ registerLagBodyForceFunction()

virtual void IBAMR::FEMechanicsBase::registerLagBodyForceFunction ( const LagBodyForceFcnData data,
unsigned int  part = 0 
)
virtual

Register the (optional) function to compute body force distributions on the Lagrangian finite element mesh.

Note
It is NOT possible to register multiple body force functions with this class.

◆ getLagBodyForceFunction()

LagBodyForceFcnData IBAMR::FEMechanicsBase::getLagBodyForceFunction ( unsigned int  part = 0) const

Get the Lagrangian body force function data.

◆ registerLagSurfacePressureFunction()

virtual void IBAMR::FEMechanicsBase::registerLagSurfacePressureFunction ( const LagSurfacePressureFcnData data,
unsigned int  part = 0 
)
virtual

Register the (optional) function to compute surface pressure distributions on the Lagrangian finite element mesh.

Note
It is NOT possible to register multiple pressure functions with this class.

◆ getLagSurfacePressureFunction()

LagSurfacePressureFcnData IBAMR::FEMechanicsBase::getLagSurfacePressureFunction ( unsigned int  part = 0) const

Get the Lagrangian surface pressure function data.

◆ registerLagSurfaceForceFunction()

virtual void IBAMR::FEMechanicsBase::registerLagSurfaceForceFunction ( const LagSurfaceForceFcnData data,
unsigned int  part = 0 
)
virtual

Register the (optional) function to compute surface force distributions on the Lagrangian finite element mesh.

Note
It is NOT possible to register multiple surface force functions with this class.

◆ getLagSurfaceForceFunction()

LagSurfaceForceFcnData IBAMR::FEMechanicsBase::getLagSurfaceForceFunction ( unsigned int  part = 0) const

Get the Lagrangian surface force function data.

◆ registerStaticPressurePart()

virtual void IBAMR::FEMechanicsBase::registerStaticPressurePart ( PressureProjectionType  projection_type = CONSISTENT_PROJECTION,
VolumetricEnergyDerivativeFcn  dU_dJ_fcn = nullptr,
unsigned int  part = 0 
)
virtual

Indicate that a part should include a static pressure.

The pressure is determined via (P, Q) = (U'(J), Q), using either a consistent or lumped mass matrix, or via a locally stabilized projection of the form (P, Q) + epsilon (P - Pi P, Q - Pi Q) = (U'(J), Q), in which P is the pressure and Q is an arbitrary test function.

Users can provide a function to evaluate U'(J). If no function is provided, we default to using U(J) = -kappa (J ln(J) − J + 1), so that U'(J) = -kappa ln J. (Ref: C.H. Liu, G. Hofstetter, H.A. Mang, 3D finite element analysis of rubber-like materials at finite strains, Eng. Comput. 11 (2) (1994) 111–128.)

The sign convention used in the implementation generates a PK1 stress of the form PP = -J P FF^{-T}.

Note
The same part cannot have both static and dynamic pressures.
See also
registerDynamicPressurePart

Reimplemented in IBAMR::IBFEMethod.

◆ registerDynamicPressurePart()

virtual void IBAMR::FEMechanicsBase::registerDynamicPressurePart ( PressureProjectionType  projection_type = CONSISTENT_PROJECTION,
VolumetricEnergyDerivativeFcn  d2U_dJ2_fcn = nullptr,
unsigned int  part = 0 
)
virtual

Indicate that a part should include a dynamic pressure.

The pressure is determined via (P-dot, Q) = (J U''(J) (FF : Grad U), Q), using either a consistent or lumped mass matrix, or via a locally stabilized projection of the form (P, Q) + epsilon (P - Pi P, Q - Pi Q) = (U'(J), Q), in which P is the pressure and Q is an arbitrary test function.

Users can provide a function to evaluate U''(J). If no function is provided, we default to using U(J) = -kappa (J ln(J) − J + 1), so that U'(J) = -kappa ln J and U''(J) = -kappa J^{-1}. (Ref: C.H. Liu, G. Hofstetter, H.A. Mang, 3D finite element analysis of rubber-like materials at finite strains, Eng. Comput. 11 (2) (1994) 111–128.)

The sign convention used in the implementation generates a PK1 stress of the form PP = -J P FF^{-T}.

Note
The same part cannot have both static and dynamic pressures.
See also
registerStaticPressurePart

◆ partHasPressure()

bool IBAMR::FEMechanicsBase::partHasPressure ( unsigned int  part = 0)
inline

Indicate whether there is a static or dynamic pressure associated with the specified FE mesh part.

◆ preprocessIntegrateData()

virtual void IBAMR::FEMechanicsBase::preprocessIntegrateData ( double  current_time,
double  new_time,
int  num_cycles 
)
virtual

Method to prepare to advance data from current_time to new_time.

Reimplemented in IBAMR::IBFEMethod, and IBAMR::FEMechanicsExplicitIntegrator.

◆ postprocessIntegrateData()

virtual void IBAMR::FEMechanicsBase::postprocessIntegrateData ( double  current_time,
double  new_time,
int  num_cycles 
)
virtual

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

Reimplemented in IBAMR::IBFEMethod, and IBAMR::FEMechanicsExplicitIntegrator.

◆ initializeFEEquationSystems()

virtual void IBAMR::FEMechanicsBase::initializeFEEquationSystems ( )
virtual

Initialize the FE equation systems objects. This method must be called prior to calling initializeFEData().

◆ initializeFEData()

virtual void IBAMR::FEMechanicsBase::initializeFEData ( )
virtual

Initialize FE data. This method must be called prior to calling IBHierarchyIntegrator::initializePatchHierarchy().

◆ reinitializeFEData()

virtual void IBAMR::FEMechanicsBase::reinitializeFEData ( )
virtual

Reinitialize FE data by calling reinit on each part's EquationSystem, reassembling the system matrices, and setting boundary conditions.

◆ putToDatabase()

virtual void IBAMR::FEMechanicsBase::putToDatabase ( SAMRAI::tbox::Pointer< SAMRAI::tbox::Database db)
overridevirtual

Write out object state to the given database.

Implements SAMRAI::tbox::Serializable.

Reimplemented in IBAMR::IBFEMethod, and IBAMR::FEMechanicsExplicitIntegrator.

◆ writeFEDataToRestartFile()

virtual void IBAMR::FEMechanicsBase::writeFEDataToRestartFile ( const std::string &  restart_dump_dirname,
unsigned int  time_step_number 
)
virtual

For technical reasons this class does not use SAMRAI's RestartManager, so restart files must be separately written for the FE objects. This function saves the solutions to the defined EquationSystems in an xdr file in restart_dump_dirname for each FE part. An example snippet is included below to show the distinct FE restart data saving step. The data will then be automatically read back into the system along with the RestartManager data during restart.

if (dump_restart_data && (iteration_num % restart_dump_interval == 0 || last_step))
{
RestartManager::getManager()->writeRestartFile(restart_dump_dirname, iteration_num);
fe_mechanics_base->writeFEDataToRestartFile(restart_dump_dirname, iteration_num);
}

◆ doInitializeFEEquationSystems()

virtual void IBAMR::FEMechanicsBase::doInitializeFEEquationSystems ( )
protectedvirtual

Do the actual work in initializeFEEquationSystems.

Reimplemented in IBAMR::IBFEMethod, and IBAMR::FEMechanicsExplicitIntegrator.

◆ doInitializeFESystemVectors()

virtual void IBAMR::FEMechanicsBase::doInitializeFESystemVectors ( )
protectedvirtual

Do the actual work of setting up libMesh system vectors.

Reimplemented in IBAMR::IBFEMethod, and IBAMR::FEMechanicsExplicitIntegrator.

◆ doInitializeFEData()

virtual void IBAMR::FEMechanicsBase::doInitializeFEData ( bool  use_present_data)
protectedvirtual

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.

Reimplemented in IBAMR::IBFEMethod, and IBAMR::FEMechanicsExplicitIntegrator.

◆ computeStaticPressure()

void IBAMR::FEMechanicsBase::computeStaticPressure ( libMesh::PetscVector< double > &  P_vec,
libMesh::PetscVector< double > &  X_vec,
double  data_time,
unsigned int  part 
)
protected

◆ computeDynamicPressureRateOfChange()

void IBAMR::FEMechanicsBase::computeDynamicPressureRateOfChange ( libMesh::PetscVector< double > &  dP_dt_vec,
libMesh::PetscVector< double > &  X_vec,
libMesh::PetscVector< double > &  U_vec,
double  data_time,
unsigned int  part 
)
protected

◆ assembleInteriorForceDensityRHS()

virtual void IBAMR::FEMechanicsBase::assembleInteriorForceDensityRHS ( libMesh::PetscVector< double > &  F_rhs_vec,
libMesh::PetscVector< double > &  X_vec,
libMesh::PetscVector< double > *  P_vec,
double  data_time,
unsigned int  part 
)
protectedvirtual

Assemble the RHS for the interior elastic density, possibly splitting off the normal component of the transmission force along the physical boundary of the Lagrangian structure.

◆ initializeCoordinates()

virtual void IBAMR::FEMechanicsBase::initializeCoordinates ( unsigned int  part)
protectedvirtual

Initialize the physical coordinates using the supplied coordinate mapping function. If no function is provided, the initial coordinates are taken to be the Lagrangian coordinates.

◆ updateCoordinateMapping()

virtual void IBAMR::FEMechanicsBase::updateCoordinateMapping ( unsigned int  part)
protectedvirtual

Compute dX = x - X, useful mainly for visualization purposes.

◆ initializeVelocity()

virtual void IBAMR::FEMechanicsBase::initializeVelocity ( unsigned int  part)
protectedvirtual

Initialize the velocity field using the supplied initial velocity specification function. If no function is provided, the initial velocity is taken to be zero.

◆ getLibMeshRestartFileName()

virtual std::string IBAMR::FEMechanicsBase::getLibMeshRestartFileName ( const std::string &  restart_dump_dirname,
unsigned int  time_step_number,
unsigned int  part,
const std::string &  extension 
) const
protectedvirtual

Get the libMesh restart file name.

◆ setup_system_vectors()

static void IBAMR::FEMechanicsBase::setup_system_vectors ( libMesh::EquationSystems *  equation_systems,
const std::vector< std::string > &  system_names,
const std::vector< std::string > &  vector_names 
)
staticprotected

Convenience function to setup system vectors and, if necessary, convert PARALLEL vectors into GHOSTED vectors for a collection of Systems.

Deprecated:
use IBTK::setup_system_vectors instead.

◆ setup_system_vector()

static void IBAMR::FEMechanicsBase::setup_system_vector ( libMesh::System &  system,
const std::string &  vector_name 
)
staticprotected

Convenience function to setup a system vector and, if necessary, convert a PARALLEL vector into a GHOSTED vector.

Deprecated:
use IBTK::setup_system_vector instead.

◆ commonConstructor()

void IBAMR::FEMechanicsBase::commonConstructor ( const std::string &  object_name,
const SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > &  input_db,
const std::vector< libMesh::MeshBase * > &  meshes,
bool  register_for_restart,
const std::string &  restart_read_dirname,
unsigned int  restart_restore_number 
)
private

Implementation of class constructor.

◆ getFromInput()

void IBAMR::FEMechanicsBase::getFromInput ( const SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > &  db,
bool  is_from_restart 
)
private

Read input values from a given database.

◆ getFromRestart()

void IBAMR::FEMechanicsBase::getFromRestart ( )
private

Read object state from the restart file and initialize class data members.

Member Data Documentation

◆ COORDS_SYSTEM_NAME

const std::string IBAMR::FEMechanicsBase::COORDS_SYSTEM_NAME
static

◆ COORD_MAPPING_SYSTEM_NAME

const std::string IBAMR::FEMechanicsBase::COORD_MAPPING_SYSTEM_NAME
static

◆ FORCE_SYSTEM_NAME

const std::string IBAMR::FEMechanicsBase::FORCE_SYSTEM_NAME
static

◆ PRESSURE_SYSTEM_NAME

const std::string IBAMR::FEMechanicsBase::PRESSURE_SYSTEM_NAME
static

◆ VELOCITY_SYSTEM_NAME

const std::string IBAMR::FEMechanicsBase::VELOCITY_SYSTEM_NAME
static

◆ d_fe_projector_db

SAMRAI::tbox::Pointer<SAMRAI::tbox::Database> IBAMR::FEMechanicsBase::d_fe_projector_db
protected

Cached input databases.

◆ d_do_log

bool IBAMR::FEMechanicsBase::d_do_log = false
protected

Indicates whether the integrator should output logging messages.

◆ d_current_time

double IBAMR::FEMechanicsBase::d_current_time = std::numeric_limits<double>::quiet_NaN()
protected

◆ d_new_time

double IBAMR::FEMechanicsBase::d_new_time = std::numeric_limits<double>::quiet_NaN()
protected

◆ d_half_time

double IBAMR::FEMechanicsBase::d_half_time = std::numeric_limits<double>::quiet_NaN()
protected

◆ d_meshes

std::vector<libMesh::MeshBase*> IBAMR::FEMechanicsBase::d_meshes
protected

Meshes provided to this object. These are set up and managed outside this class. These meshes are modified by FEMechanicsBase since this class creates several libMesh Systems (and hence stores DoF information in these meshes).

◆ d_diagonal_system_coupling

libMesh::CouplingMatrix IBAMR::FEMechanicsBase::d_diagonal_system_coupling
protected

The libMesh Systems set up by this system (for example, for velocity projection) consist of one variable per spatial component. By default, libMesh assumes that all variables in a given System couple to each other which, since we only ever solve projection problems in this class, is not the case. Hence we can save some memory by explicitly informing libMesh that the variables in a system only couple to themselves by providing a diagonal coupling matrix to each System.

◆ d_equation_systems

std::vector<std::unique_ptr<libMesh::EquationSystems> > IBAMR::FEMechanicsBase::d_equation_systems
protected

EquationSystems objects, one per part. These contain the actual matrices and solution vectors for each relevant libMesh system.

◆ d_current_coordinates_system_name

const std::string IBAMR::FEMechanicsBase::d_current_coordinates_system_name
protected

◆ d_displacement_system_name

const std::string IBAMR::FEMechanicsBase::d_displacement_system_name
protected

◆ d_force_system_name

const std::string IBAMR::FEMechanicsBase::d_force_system_name
protected

◆ d_pressure_system_name

const std::string IBAMR::FEMechanicsBase::d_pressure_system_name
protected

◆ d_velocity_system_name

const std::string IBAMR::FEMechanicsBase::d_velocity_system_name
protected

◆ d_fe_data

std::vector<std::shared_ptr<IBTK::FEData> > IBAMR::FEMechanicsBase::d_fe_data
protected

◆ d_fe_projectors

std::vector<std::shared_ptr<IBTK::FEProjector> > IBAMR::FEMechanicsBase::d_fe_projectors
protected

◆ d_X_systems

std::vector<libMesh::ExplicitSystem*> IBAMR::FEMechanicsBase::d_X_systems
protected

Vectors of pointers to the systems for each part (for position, velocity, force, and pressure).

◆ d_U_systems

std::vector<libMesh::ExplicitSystem*> IBAMR::FEMechanicsBase::d_U_systems
protected

◆ d_F_systems

std::vector<libMesh::ExplicitSystem*> IBAMR::FEMechanicsBase::d_F_systems
protected

◆ d_P_systems

std::vector<libMesh::ExplicitSystem*> IBAMR::FEMechanicsBase::d_P_systems
protected

◆ d_X_vecs

std::unique_ptr<IBTK::LibMeshSystemVectors> IBAMR::FEMechanicsBase::d_X_vecs
protected

Object managing access to libMesh system vectors for the position.

◆ d_U_vecs

std::unique_ptr<IBTK::LibMeshSystemVectors> IBAMR::FEMechanicsBase::d_U_vecs
protected

Object managing access to libMesh system vectors for the velocity.

◆ d_F_vecs

std::unique_ptr<IBTK::LibMeshSystemVectors> IBAMR::FEMechanicsBase::d_F_vecs
protected

Object managing access to libMesh system vectors for the force.

◆ d_P_vecs

std::unique_ptr<IBTK::LibMeshSystemVectors> IBAMR::FEMechanicsBase::d_P_vecs
protected

Object managing access to libMesh system vectors for the pressure.

◆ d_fe_equation_systems_initialized

bool IBAMR::FEMechanicsBase::d_fe_equation_systems_initialized = false
protected

Whether or not the libMesh equation systems objects have been initialized (i.e., whether or not initializeFEEquationSystems has been called).

◆ d_fe_data_initialized

bool IBAMR::FEMechanicsBase::d_fe_data_initialized = false
protected

Whether or not all finite element data (including that initialized by initializeFEEquationSystems), such system matrices, is available.

◆ d_libmesh_partitioner_type

LibmeshPartitionerType IBAMR::FEMechanicsBase::d_libmesh_partitioner_type = LIBMESH_DEFAULT
protected

Type of partitioner to use. See the main documentation of this class for more information.

◆ d_libmesh_use_amr

bool IBAMR::FEMechanicsBase::d_libmesh_use_amr = false
protected

Whether or not to use AMR in the finite element discretization. This feature is not yet implemented and currently defaults to false.

◆ d_fe_order_position

std::vector<libMesh::Order> IBAMR::FEMechanicsBase::d_fe_order_position
protected

Method parameters.

◆ d_fe_order_force

std::vector<libMesh::Order> IBAMR::FEMechanicsBase::d_fe_order_force
protected

◆ d_fe_order_pressure

std::vector<libMesh::Order> IBAMR::FEMechanicsBase::d_fe_order_pressure
protected

◆ d_fe_family_position

std::vector<libMesh::FEFamily> IBAMR::FEMechanicsBase::d_fe_family_position
protected

◆ d_fe_family_force

std::vector<libMesh::FEFamily> IBAMR::FEMechanicsBase::d_fe_family_force
protected

◆ d_fe_family_pressure

std::vector<libMesh::FEFamily> IBAMR::FEMechanicsBase::d_fe_family_pressure
protected

◆ d_default_quad_type_stress

std::vector<libMesh::QuadratureType> IBAMR::FEMechanicsBase::d_default_quad_type_stress
protected

◆ d_default_quad_type_force

std::vector<libMesh::QuadratureType> IBAMR::FEMechanicsBase::d_default_quad_type_force
protected

◆ d_default_quad_type_pressure

std::vector<libMesh::QuadratureType> IBAMR::FEMechanicsBase::d_default_quad_type_pressure
protected

◆ d_default_quad_order_stress

std::vector<libMesh::Order> IBAMR::FEMechanicsBase::d_default_quad_order_stress
protected

◆ d_default_quad_order_force

std::vector<libMesh::Order> IBAMR::FEMechanicsBase::d_default_quad_order_force
protected

◆ d_default_quad_order_pressure

std::vector<libMesh::Order> IBAMR::FEMechanicsBase::d_default_quad_order_pressure
protected

◆ d_use_consistent_mass_matrix

bool IBAMR::FEMechanicsBase::d_use_consistent_mass_matrix = true
protected

◆ d_allow_rules_with_negative_weights

bool IBAMR::FEMechanicsBase::d_allow_rules_with_negative_weights = true
protected

◆ d_include_normal_stress_in_weak_form

bool IBAMR::FEMechanicsBase::d_include_normal_stress_in_weak_form = false
protected

◆ d_include_tangential_stress_in_weak_form

bool IBAMR::FEMechanicsBase::d_include_tangential_stress_in_weak_form = false
protected

◆ d_include_normal_surface_forces_in_weak_form

bool IBAMR::FEMechanicsBase::d_include_normal_surface_forces_in_weak_form = true
protected

◆ d_include_tangential_surface_forces_in_weak_form

bool IBAMR::FEMechanicsBase::d_include_tangential_surface_forces_in_weak_form = true
protected

◆ d_coordinate_mapping_fcn_data

std::vector<CoordinateMappingFcnData> IBAMR::FEMechanicsBase::d_coordinate_mapping_fcn_data
protected

Functions used to compute the initial coordinates of the Lagrangian mesh.

◆ d_initial_velocity_fcn_data

std::vector<InitialVelocityFcnData> IBAMR::FEMechanicsBase::d_initial_velocity_fcn_data
protected

Functions used to compute the initial coordinates of the Lagrangian mesh.

◆ d_PK1_stress_fcn_data

std::vector<std::vector<PK1StressFcnData> > IBAMR::FEMechanicsBase::d_PK1_stress_fcn_data
protected

Functions used to compute the first Piola-Kirchhoff stress tensor.

◆ d_lag_body_force_fcn_data

std::vector<LagBodyForceFcnData> IBAMR::FEMechanicsBase::d_lag_body_force_fcn_data
protected

Functions used to compute additional body and surface forces on the Lagrangian mesh.

◆ d_lag_surface_pressure_fcn_data

std::vector<LagSurfacePressureFcnData> IBAMR::FEMechanicsBase::d_lag_surface_pressure_fcn_data
protected

◆ d_lag_surface_force_fcn_data

std::vector<LagSurfaceForceFcnData> IBAMR::FEMechanicsBase::d_lag_surface_force_fcn_data
protected

◆ d_static_pressure_kappa

double IBAMR::FEMechanicsBase::d_static_pressure_kappa = 0.0
protected

Data related to handling static pressures.

◆ d_static_pressure_stab_param

double IBAMR::FEMechanicsBase::d_static_pressure_stab_param = 0.0
protected

◆ d_has_static_pressure_parts

bool IBAMR::FEMechanicsBase::d_has_static_pressure_parts = false
protected

◆ d_static_pressure_part

std::vector<bool> IBAMR::FEMechanicsBase::d_static_pressure_part
protected

◆ d_static_pressure_proj_type

std::vector<PressureProjectionType> IBAMR::FEMechanicsBase::d_static_pressure_proj_type
protected

◆ d_static_pressure_dU_dJ_fcn

std::vector<VolumetricEnergyDerivativeFcn> IBAMR::FEMechanicsBase::d_static_pressure_dU_dJ_fcn
protected

◆ d_dynamic_pressure_kappa

double IBAMR::FEMechanicsBase::d_dynamic_pressure_kappa = 0.0
protected

Data related to handling dynamic pressures.

◆ d_dynamic_pressure_stab_param

double IBAMR::FEMechanicsBase::d_dynamic_pressure_stab_param = 0.0
protected

◆ d_has_dynamic_pressure_parts

bool IBAMR::FEMechanicsBase::d_has_dynamic_pressure_parts = false
protected

◆ d_dynamic_pressure_part

std::vector<bool> IBAMR::FEMechanicsBase::d_dynamic_pressure_part
protected

◆ d_dynamic_pressure_proj_type

std::vector<PressureProjectionType> IBAMR::FEMechanicsBase::d_dynamic_pressure_proj_type
protected

◆ d_dynamic_pressure_d2U_dJ2_fcn

std::vector<VolumetricEnergyDerivativeFcn> IBAMR::FEMechanicsBase::d_dynamic_pressure_d2U_dJ2_fcn
protected

◆ d_object_name

std::string IBAMR::FEMechanicsBase::d_object_name
protected

The object name is used as a handle to databases stored in restart files and for error reporting purposes.

◆ d_registered_for_restart

bool IBAMR::FEMechanicsBase::d_registered_for_restart
protected

A boolean value indicating whether the class is registered with the restart database.

◆ d_libmesh_restart_read_dir

std::string IBAMR::FEMechanicsBase::d_libmesh_restart_read_dir
protected

Directory and time step number to use when restarting.

◆ d_libmesh_restart_restore_number

unsigned int IBAMR::FEMechanicsBase::d_libmesh_restart_restore_number
protected

◆ d_libmesh_restart_file_extension

std::string IBAMR::FEMechanicsBase::d_libmesh_restart_file_extension
protected

Restart file type for libMesh equation systems (e.g. xda or xdr).


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