libMesh
|
This class implements a constitutive formulation for an Neo-Hookean elastic solid in terms of the current configuration. More...
#include <nonlinear_neohooke_cc.h>
Public Member Functions | |
NonlinearNeoHookeCurrentConfig (const std::vector< std::vector< RealGradient >> &dphi_in, GetPot &args, bool calculate_linearized_stiffness_in) | |
void | init_for_qp (VectorValue< Gradient > &grad_u, unsigned int qp) |
Initialize the class for the given displacement gradient at the specified quadrature point. More... | |
void | get_residual (DenseVector< Real > &residuum, unsigned int &i) |
Return the residual vector for the current state. More... | |
void | get_linearized_stiffness (DenseMatrix< Real > &stiffness, unsigned int &i, unsigned int &j) |
Return the stiffness matrix for the current state. More... | |
Public Attributes | |
bool | calculate_linearized_stiffness |
Flag to indicate if it is necessary to calculate values for stiffness matrix during initialization. More... | |
Private Member Functions | |
void | build_b_0_mat (int i, DenseMatrix< Real > &b_l_mat) |
void | calculate_stress () |
void | calculate_tangent () |
Static Private Member Functions | |
static void | tensor_to_voigt (const RealTensor &tensor, DenseVector< Real > &vec) |
Private Attributes | |
unsigned int | current_qp |
const std::vector< std::vector< RealGradient > > & | dphi |
DenseMatrix< Real > | C_mat |
Real | E |
Real | nu |
RealTensor | F |
RealTensor | S |
RealTensor | tau |
RealTensor | sigma |
DenseMatrix< Real > | B_L |
DenseMatrix< Real > | B_K |
This class implements a constitutive formulation for an Neo-Hookean elastic solid in terms of the current configuration.
This implementation is suitable for computing large deformations. See e.g. "Nonlinear finite element methods" (P. Wriggers, 2008, Springer) for details.
Definition at line 40 of file nonlinear_neohooke_cc.h.
|
inline |
Definition at line 43 of file nonlinear_neohooke_cc.h.
|
private |
Definition at line 141 of file nonlinear_neohooke_cc.C.
References current_qp, and dphi.
Referenced by get_linearized_stiffness(), and get_residual().
|
private |
Definition at line 73 of file nonlinear_neohooke_cc.C.
References E, F, libMesh::TypeTensor< T >::inverse(), nu, libMesh::Real, S, sigma, and tau.
Referenced by init_for_qp().
|
private |
Definition at line 52 of file nonlinear_neohooke_cc.C.
References C_mat, E, F, nu, libMesh::Real, and libMesh::DenseMatrix< T >::resize().
Referenced by init_for_qp().
void NonlinearNeoHookeCurrentConfig::get_linearized_stiffness | ( | DenseMatrix< Real > & | stiffness, |
unsigned int & | i, | ||
unsigned int & | j | ||
) |
Return the stiffness matrix for the current state.
Definition at line 118 of file nonlinear_neohooke_cc.C.
References B_K, B_L, build_b_0_mat(), C_mat, current_qp, dphi, F, libMesh::Real, libMesh::DenseMatrix< T >::resize(), libMesh::DenseMatrix< T >::right_multiply(), libMesh::DenseMatrix< T >::right_multiply_transpose(), and sigma.
Referenced by SolidSystem::element_time_derivative().
void NonlinearNeoHookeCurrentConfig::get_residual | ( | DenseVector< Real > & | residuum, |
unsigned int & | i | ||
) |
Return the residual vector for the current state.
Definition at line 93 of file nonlinear_neohooke_cc.C.
References B_L, build_b_0_mat(), libMesh::DenseVector< Real >, libMesh::DenseMatrix< T >::resize(), sigma, tensor_to_voigt(), and libMesh::DenseMatrix< T >::vector_mult().
Referenced by SolidSystem::element_time_derivative().
void NonlinearNeoHookeCurrentConfig::init_for_qp | ( | VectorValue< Gradient > & | grad_u, |
unsigned int | qp | ||
) |
Initialize the class for the given displacement gradient at the specified quadrature point.
Definition at line 25 of file nonlinear_neohooke_cc.C.
References calculate_linearized_stiffness, calculate_stress(), calculate_tangent(), current_qp, F, libMesh::libmesh_real(), S, libMesh::TOLERANCE, and libMesh::TypeTensor< T >::zero().
Referenced by SolidSystem::element_time_derivative().
|
staticprivate |
Definition at line 106 of file nonlinear_neohooke_cc.C.
Referenced by get_residual().
|
private |
Definition at line 93 of file nonlinear_neohooke_cc.h.
Referenced by get_linearized_stiffness().
|
private |
Definition at line 92 of file nonlinear_neohooke_cc.h.
Referenced by get_linearized_stiffness(), and get_residual().
|
private |
Definition at line 88 of file nonlinear_neohooke_cc.h.
Referenced by calculate_tangent(), and get_linearized_stiffness().
bool NonlinearNeoHookeCurrentConfig::calculate_linearized_stiffness |
Flag to indicate if it is necessary to calculate values for stiffness matrix during initialization.
Definition at line 77 of file nonlinear_neohooke_cc.h.
Referenced by init_for_qp().
|
private |
Definition at line 85 of file nonlinear_neohooke_cc.h.
Referenced by build_b_0_mat(), get_linearized_stiffness(), and init_for_qp().
|
private |
Definition at line 86 of file nonlinear_neohooke_cc.h.
Referenced by build_b_0_mat(), and get_linearized_stiffness().
|
private |
Definition at line 89 of file nonlinear_neohooke_cc.h.
Referenced by calculate_stress(), and calculate_tangent().
|
private |
Definition at line 91 of file nonlinear_neohooke_cc.h.
Referenced by calculate_stress(), calculate_tangent(), get_linearized_stiffness(), and init_for_qp().
|
private |
Definition at line 90 of file nonlinear_neohooke_cc.h.
Referenced by calculate_stress(), and calculate_tangent().
|
private |
Definition at line 91 of file nonlinear_neohooke_cc.h.
Referenced by calculate_stress(), and init_for_qp().
|
private |
Definition at line 91 of file nonlinear_neohooke_cc.h.
Referenced by calculate_stress(), get_linearized_stiffness(), and get_residual().
|
private |
Definition at line 91 of file nonlinear_neohooke_cc.h.
Referenced by calculate_stress().