20 #ifndef __mast_elasticity_continuum_linear_acceleration_kernel_h__ 21 #define __mast_elasticity_continuum_linear_acceleration_kernel_h__ 31 namespace Elasticity {
32 namespace LinearContinuum {
35 template <
typename NodalScalarType,
36 typename VarScalarType,
42 typename Eigen::Matrix<VarScalarType, Dim, 1>& u,
47 const typename FEVarType::fe_shape_deriv_t
48 &fe = fe_var.get_fe_shape_data();
53 "Incompatible operator size");
54 Assert2(Bmat.
n() == Dim*fe.n_basis(),
55 Bmat.
n(), Dim*fe.n_basis(),
56 "Incompatible Operator size.");
60 for (
uint_t i=0; i<Dim; i++) {
62 u(i) = fe_var.u(qp, i);
81 template <
typename FEVarType,
82 typename DensityFieldType,
83 typename SectionAreaType,
92 using vector_t =
typename Eigen::Matrix<scalar_t, Eigen::Dynamic, 1>;
93 using matrix_t =
typename Eigen::Matrix<scalar_t, Eigen::Dynamic, Eigen::Dynamic>;
131 const typename FEVarType::fe_shape_deriv_t
134 typename Eigen::Matrix<scalar_t, Dim, 1>
137 vec = vector_t::Zero(Dim*fe.n_basis());
144 mat = matrix_t::Zero(Dim*fe.n_basis(), Dim*fe.n_basis());
148 Bmat.
reinit(Dim, Dim, fe.n_basis());
150 for (
uint_t i=0; i<fe.n_q_points(); i++) {
159 Bmat.vector_mult_transpose(vec, u);
160 res += fe.detJxW(i) * vec * rho * sec;
164 Bmat.right_multiply_transpose(mat, Bmat);
165 (*jac) += fe.detJxW(i) * mat * rho * sec;
179 template <
typename ScalarFieldType>
181 const ScalarFieldType& f,
189 const typename FEVarType::fe_shape_deriv_t
192 typename Eigen::Matrix<scalar_t, Dim, 1>
195 vec = vector_t::Zero(Dim*fe.n_basis());
204 mat = matrix_t::Zero(Dim*fe.n_basis(), Dim*fe.n_basis());
208 Bmat.
reinit(Dim, Dim, fe.n_basis());
210 for (
uint_t i=0; i<fe.n_q_points(); i++) {
221 Bmat.vector_mult_transpose(vec, u);
222 res += fe.detJxW(i) * vec * (drho * sec + rho * dsec);
226 Bmat.right_multiply_transpose(mat, Bmat);
227 (*jac) += fe.detJxW(i) * mat * (drho * sec + rho * dsec);
246 #endif // __mast_elasticity_continuum_linear_acceleration_kernel_h__ void derivative(ContextType &c, const ScalarFieldType &f, vector_t &res, matrix_t *jac=nullptr) const
Computes the derivative of residual of variational term with respect to parameter and returns it in...
void set_fe_var_data(const FEVarType &fe)
const DensityFieldType * _density
void set_density(const DensityFieldType &rho)
#define Assert1(cond, v1, msg)
This class implements the discrete evaluation of the acceleration kernel defined as where...
virtual ~LinearAcceleration()
void compute(ContextType &c, vector_t &res, matrix_t *jac=nullptr) const
Computes the residual of variational term and returns it in res.
typename FEVarType::scalar_t scalar_t
const FEVarType * _fe_var_data
void displacement(const FEVarType &fe_var, const uint_t qp, typename Eigen::Matrix< VarScalarType, Dim, 1 > &u, MAST::Numerics::FEMOperatorMatrix< NodalScalarType > &Bmat)
const SectionAreaType * _section
#define Assert0(cond, msg)
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)
typename Eigen::Matrix< scalar_t, Eigen::Dynamic, 1 > vector_t
void set_section_area(const SectionAreaType &s)
typename FEVarType::fe_shape_deriv_t::scalar_t basis_scalar_t
void reinit(uint_t n_interpolated_vars, uint_t n_discrete_vars, uint_t n_discrete_dofs_per_var)
this initializes the operator for number of rows and variables, assuming that all variables has the s...
typename Eigen::Matrix< scalar_t, Eigen::Dynamic, Eigen::Dynamic > matrix_t