MAST3
Multidisciplinary-design Adaptation and Sensitivity Toolkit (MAST)
MAST::Physics::Elasticity Namespace Reference

Namespaces

 libMeshWrapper
 
 LinearContinuum
 
 MindlinPlate
 
 Plate
 

Classes

class  IsotropicMaterialStiffness
 
class  IsotropicMaterialStiffness< ScalarType, 2, ModulusType, PoissonType, ContextType >
 
class  IsotropicMaterialStiffness< ScalarType, 3, ModulusType, PoissonType, ContextType >
 
class  PlateBendingSectionProperty
 
class  PlateInertiaSectionProperty
 
class  ShellFacePressureLoad
 
class  SurfacePressureLoad
 
class  SurfaceTractionLoad
 This class implements the discrete evaluation of the surface traction kernel defined as

\[ - \int_{\Gamma_e} a \phi t \cdot \hat{n} ~d\Gamma, \]

where, $ \phi$ is the variation, $ t $ is the surface traction, $ a $ is the section thickness for 2D elements or section area for 1D elements, and $ \hat{n} $ is the surface normal. More...

 
class  Traction
 This class defines a data structure that can be used to define a parameterized surface traction vector, where each component of the Dim dimensional vector is scalar of type TractionScalarType, which is typically a ScalarConstant. More...
 

Functions

template<typename NodalScalarType , typename VarScalarType , typename FEVarType , uint_t Dim>
std::enable_if< Dim==2, void >::type green_lagrange_strain_operator (const FEVarType &fe_var, const uint_t qp, typename Eigen::Matrix< VarScalarType, 2, 2 > &E, typename Eigen::Matrix< VarScalarType, 2, 2 > &F, typename Eigen::Matrix< VarScalarType, 3, 1 > &epsilon, typename Eigen::Matrix< VarScalarType, 3, 2 > &mat_x, typename Eigen::Matrix< VarScalarType, 3, 2 > &mat_y, MAST::Numerics::FEMOperatorMatrix< NodalScalarType > &Bmat_lin, MAST::Numerics::FEMOperatorMatrix< NodalScalarType > &Bmat_nl_x, MAST::Numerics::FEMOperatorMatrix< NodalScalarType > &Bmat_nl_y, MAST::Numerics::FEMOperatorMatrix< NodalScalarType > &Bmat_nl_u, MAST::Numerics::FEMOperatorMatrix< NodalScalarType > &Bmat_nl_v)
 
template<typename ScalarType >
ScalarType shear_modulus (ScalarType E, ScalarType nu)
 
template<typename ScalarType >
ScalarType shear_modulus_derivative (ScalarType E, ScalarType nu, ScalarType dE, ScalarType dnu)
 
template<typename TractionType , uint_t Dim, typename ContextType , typename ScalarFieldType >
std::enable_if< Dim==1, void >::type traction_derivative (const ScalarFieldType &f, ContextType &c, const TractionType &t, typename TractionType::value_t &v)
 
template<typename TractionType , uint_t Dim, typename ContextType , typename ScalarFieldType >
std::enable_if< Dim==2, void >::type traction_derivative (const ScalarFieldType &f, ContextType &c, const TractionType &t, typename TractionType::value_t &v)
 
template<typename TractionType , uint_t Dim, typename ContextType , typename ScalarFieldType >
std::enable_if< Dim==3, void >::type traction_derivative (const ScalarFieldType &f, ContextType &c, const TractionType &t, typename TractionType::value_t &v)
 
template<typename TractionType , uint_t Dim, typename ContextType >
std::enable_if< Dim==1, void >::type traction_value (ContextType &c, const TractionType &t, typename TractionType::value_t &v)
 
template<typename TractionType , uint_t Dim, typename ContextType >
std::enable_if< Dim==2, void >::type traction_value (ContextType &c, const TractionType &t, typename TractionType::value_t &v)
 
template<typename TractionType , uint_t Dim, typename ContextType >
std::enable_if< Dim==3, void >::type traction_value (ContextType &c, const TractionType &t, typename TractionType::value_t &v)
 

Function Documentation

◆ green_lagrange_strain_operator()

template<typename NodalScalarType , typename VarScalarType , typename FEVarType , uint_t Dim>
std::enable_if<Dim == 2, void>::type MAST::Physics::Elasticity::green_lagrange_strain_operator ( const FEVarType &  fe_var,
const uint_t  qp,
typename Eigen::Matrix< VarScalarType, 2, 2 > &  E,
typename Eigen::Matrix< VarScalarType, 2, 2 > &  F,
typename Eigen::Matrix< VarScalarType, 3, 1 > &  epsilon,
typename Eigen::Matrix< VarScalarType, 3, 2 > &  mat_x,
typename Eigen::Matrix< VarScalarType, 3, 2 > &  mat_y,
MAST::Numerics::FEMOperatorMatrix< NodalScalarType > &  Bmat_lin,
MAST::Numerics::FEMOperatorMatrix< NodalScalarType > &  Bmat_nl_x,
MAST::Numerics::FEMOperatorMatrix< NodalScalarType > &  Bmat_nl_y,
MAST::Numerics::FEMOperatorMatrix< NodalScalarType > &  Bmat_nl_u,
MAST::Numerics::FEMOperatorMatrix< NodalScalarType > &  Bmat_nl_v 
)
inline

Definition at line 37 of file green_lagrange_strain.hpp.

48  {
49 
50 
51  epsilon.setZero();
52  mat_x.setZero();
53  mat_y.setZero();
54 
55  const typename FEVarType::fe_shape_data_type
56  &fe = fe_var.get_fe_shape_object();
57 
58  // make sure all matrices are the right size
59  Assert1(epsilon.size() == 3,
60  epsilon.size(),
61  "Strain vector for 2D continuum strain should be 3.");
62  Assert1(mat_x.rows() == 3,
63  mat_x.rows(),
64  "Incompatible matrix size.");
65  Assert1(mat_x.cols() == 2,
66  mat_x.cols(),
67  "Incompatible matrix size.");
68  Assert1(mat_y.rows() == 3,
69  mat_y.rows(),
70  "Incompatible matrix size.");
71  Assert1(mat_y.cols() == 2,
72  mat_y.cols(),
73  "Incompatible matrix size.");
74 
75  Assert1(Bmat_lin.m() == 3,
76  Bmat_lin.m(),
77  "Strain vector for 2D continuum strain should be 3");
78  Assert2(Bmat_lin.n() == 2*fe.n_basis(),
79  Bmat_lin.n(), 2*fe.n_basis(),
80  "Incompatible Operator size.");
81  Assert1(Bmat_nl_x.m() == 2,
82  Bmat_nl_x.m(),
83  "Incompatible matrix size.");
84  Assert2(Bmat_nl_x.n() == 2*fe.n_basis(),
85  Bmat_nl_x.n(), 2*fe.n_basis(),
86  "Incompatible matrix size.");
87  Assert1(Bmat_nl_y.m() == 2,
88  Bmat_nl_y.m(),
89  "Incompatible matrix size.");
90  Assert2(Bmat_nl_y.n() == 2*fe.n_basis(),
91  Bmat_nl_y.n(), 2*fe.n_basis(),
92  "Incompatible matrix size.");
93  Assert1(F.cols() == 2,
94  F.cols(),
95  "Incompatible matrix size.");
96  Assert1(F.rows() == 2,
97  F.rows(),
98  "Incompatible matrix size.");
99  Assert1(E.cols() == 2,
100  E.cols(),
101  "Incompatible matrix size.");
102  Assert1(E.rows() == 2,
103  E.rows(),
104  "Incompatible matrix size.");
105 
106 
107  // now set the shape function values
108  Bmat_lin.set_shape_function(0, 0, fe.dphi_dx(qp, 0)); // epsilon_xx = du/dx
109  Bmat_lin.set_shape_function(2, 1, fe.dphi_dx(qp, 0)); // gamma_xy = dv/dx + ...
110 
111  // nonlinear strain operator in x
112  Bmat_nl_x.set_shape_function(0, 0, fe.dphi_dx(qp, 0)); // du/dx
113  Bmat_nl_x.set_shape_function(1, 1, fe.dphi_dx(qp, 0)); // dv/dx
114 
115  // nonlinear strain operator in u
116  Bmat_nl_u.set_shape_function(0, 0, fe.dphi_dx(qp, 0)); // du/dx
117  Bmat_nl_v.set_shape_function(0, 1, fe.dphi_dx(qp, 0)); // dv/dx
118 
119  // dN/dy
120  Bmat_lin.set_shape_function(1, 1, fe.dphi_dx(qp, 1)); // epsilon_yy = dv/dy
121  Bmat_lin.set_shape_function(2, 0, fe.dphi_dx(qp, 1)); // gamma_xy = du/dy + ...
122 
123  // nonlinear strain operator in y
124  Bmat_nl_y.set_shape_function(0, 0, fe.dphi_dx(qp, 1)); // du/dy
125  Bmat_nl_y.set_shape_function(1, 1, fe.dphi_dx(qp, 1)); // dv/dy
126 
127  // nonlinear strain operator in v
128  Bmat_nl_u.set_shape_function(1, 0, fe.dphi_dx(qp, 1)); // du/dy
129  Bmat_nl_v.set_shape_function(1, 1, fe.dphi_dx(qp, 1)); // dv/dy
130 
131  // prepare the deformation gradient matrix
132  F.row(0) = fe_var.du_dx(qp, 0);
133  F.row(1) = fe_var.du_dx(qp, 1);
134 
135  // this calculates the Green-Lagrange strain in the reference config
136  E = 0.5*(F + F.transpose() + F.transpose() * F);
137 
138  // now, add this to the strain vector
139  epsilon(0) = E(0,0);
140  epsilon(1) = E(1,1);
141  epsilon(2) = E(0,1) + E(1,0);
142 
143  // now initialize the matrices with strain components
144  // that multiply the Bmat_nl terms
145  mat_x(0, 0) = fe_var.du_dx(qp, 0, 0);
146  mat_x(0, 1) = fe_var.du_dx(qp, 1, 0);
147  mat_x(2, 0) = fe_var.du_dx(qp, 0, 1);
148  mat_x(2, 1) = fe_var.du_dx(qp, 1, 1);
149 
150  mat_y(1, 0) = fe_var.du_dx(qp, 0, 1);
151  mat_y(1, 1) = fe_var.du_dx(qp, 1, 1);
152  mat_y(2, 0) = fe_var.du_dx(qp, 0, 0);
153  mat_y(2, 1) = fe_var.du_dx(qp, 1, 0);
154 }
#define Assert1(cond, v1, msg)
Definition: exceptions.hpp:143
void set_shape_function(uint_t interpolated_var, uint_t discrete_var, const VecType &shape_func)
sets the shape function values for the block corresponding to interpolated_var and discrete_var...
#define Assert2(cond, v1, v2, msg)
Definition: exceptions.hpp:152

◆ shear_modulus()

template<typename ScalarType >
ScalarType MAST::Physics::Elasticity::shear_modulus ( ScalarType  E,
ScalarType  nu 
)
inline

Definition at line 33 of file isotropic_stiffness.hpp.

33 { return E/2./(1.+nu);}

◆ shear_modulus_derivative()

template<typename ScalarType >
ScalarType MAST::Physics::Elasticity::shear_modulus_derivative ( ScalarType  E,
ScalarType  nu,
ScalarType  dE,
ScalarType  dnu 
)
inline

Definition at line 37 of file isotropic_stiffness.hpp.

39 { return dE/2./(1.+nu) - E/2./pow(1.+nu,2) * dnu;}

◆ traction_derivative() [1/3]

template<typename TractionType , uint_t Dim, typename ContextType , typename ScalarFieldType >
std::enable_if<Dim == 1, void>::type MAST::Physics::Elasticity::traction_derivative ( const ScalarFieldType &  f,
ContextType &  c,
const TractionType &  t,
typename TractionType::value_t &  v 
)
inline

Definition at line 79 of file traction_load.hpp.

82  {
83 
84  Assert2(v.size() == Dim, v.size(), Dim, "Incorrect size");
85 
86  v(0) = t.get_scalar_for_dim(0).derivative(c, f);
87 }
#define Assert2(cond, v1, v2, msg)
Definition: exceptions.hpp:152

◆ traction_derivative() [2/3]

template<typename TractionType , uint_t Dim, typename ContextType , typename ScalarFieldType >
std::enable_if<Dim == 2, void>::type MAST::Physics::Elasticity::traction_derivative ( const ScalarFieldType &  f,
ContextType &  c,
const TractionType &  t,
typename TractionType::value_t &  v 
)
inline

Definition at line 96 of file traction_load.hpp.

99  {
100 
101  Assert2(v.size() == Dim, v.size(), Dim, "Incorrect size");
102 
103  v(0) = t.get_scalar_for_dim(0).derivative(c, f);
104  v(1) = t.get_scalar_for_dim(1).derivative(c, f);
105 }
#define Assert2(cond, v1, v2, msg)
Definition: exceptions.hpp:152

◆ traction_derivative() [3/3]

template<typename TractionType , uint_t Dim, typename ContextType , typename ScalarFieldType >
std::enable_if<Dim == 3, void>::type MAST::Physics::Elasticity::traction_derivative ( const ScalarFieldType &  f,
ContextType &  c,
const TractionType &  t,
typename TractionType::value_t &  v 
)
inline

Definition at line 114 of file traction_load.hpp.

117  {
118 
119  Assert2(v.size() == Dim, v.size(), Dim, "Incorrect size");
120 
121  v(0) = t.get_scalar_for_dim(0).derivative(c, f);
122  v(1) = t.get_scalar_for_dim(1).derivative(c, f);
123  v(2) = t.get_scalar_for_dim(2).derivative(c, f);
124 }
#define Assert2(cond, v1, v2, msg)
Definition: exceptions.hpp:152

◆ traction_value() [1/3]

template<typename TractionType , uint_t Dim, typename ContextType >
std::enable_if<Dim == 1, void>::type MAST::Physics::Elasticity::traction_value ( ContextType &  c,
const TractionType &  t,
typename TractionType::value_t &  v 
)
inline

Definition at line 34 of file traction_load.hpp.

36  {
37 
38  Assert2(v.size() == Dim, v.size(), Dim, "Incorrect size");
39 
40  v(0) = t.get_scalar_for_dim(0).value(c);
41 }
#define Assert2(cond, v1, v2, msg)
Definition: exceptions.hpp:152

◆ traction_value() [2/3]

template<typename TractionType , uint_t Dim, typename ContextType >
std::enable_if<Dim == 2, void>::type MAST::Physics::Elasticity::traction_value ( ContextType &  c,
const TractionType &  t,
typename TractionType::value_t &  v 
)
inline

Definition at line 47 of file traction_load.hpp.

49  {
50 
51  Assert2(v.size() == Dim, v.size(), Dim, "Incorrect size");
52 
53  v(0) = t.get_scalar_for_dim(0).value(c);
54  v(1) = t.get_scalar_for_dim(1).value(c);
55 }
#define Assert2(cond, v1, v2, msg)
Definition: exceptions.hpp:152

◆ traction_value() [3/3]

template<typename TractionType , uint_t Dim, typename ContextType >
std::enable_if<Dim == 3, void>::type MAST::Physics::Elasticity::traction_value ( ContextType &  c,
const TractionType &  t,
typename TractionType::value_t &  v 
)
inline

Definition at line 61 of file traction_load.hpp.

63  {
64 
65  Assert2(v.size() == Dim, v.size(), Dim, "Incorrect size");
66 
67  v(0) = t.get_scalar_for_dim(0).value(c);
68  v(1) = t.get_scalar_for_dim(1).value(c);
69  v(2) = t.get_scalar_for_dim(2).value(c);
70 }
#define Assert2(cond, v1, v2, msg)
Definition: exceptions.hpp:152