20 #ifndef __mast_conduction_kernel_h__ 21 #define __mast_conduction_kernel_h__ 28 namespace Conduction {
32 template <
typename FEVarType,
33 typename SectionPropertyType,
36 bool IsotropicMaterial = SectionPropertyType::is_isotropic,
37 bool LinearMaterial = SectionPropertyType::is_linear>
54 template <
typename FEVarType,
55 typename SectionPropertyType,
69 using vector_t =
typename Eigen::Matrix<scalar_t, Eigen::Dynamic, 1>;
70 using matrix_t =
typename Eigen::Matrix<scalar_t, Eigen::Dynamic, Eigen::Dynamic>;
75 _fe_var_data (nullptr)
83 Assert0(!_property,
"Property already initialized.");
90 Assert0(!_fe_var_data,
"FE data already initialized.");
91 _fe_var_data = &fe_data;
96 Assert0(_fe_var_data,
"FE data not initialized.");
98 return _fe_var_data->get_fe_shape_data().n_basis();
112 Assert0(_fe_var_data,
"FE data not initialized.");
113 Assert0(_property,
"Section property not initialized");
115 const typename FEVarType::fe_shape_deriv_t
116 &fe = _fe_var_data->get_fe_shape_data();
118 typename Eigen::Matrix<scalar_t, Dim, 1>
121 vec = vector_t::Zero(fe.n_basis());
123 typename SectionPropertyType::value_t
127 mat2 = matrix_t::Zero(fe.n_basis(), fe.n_basis());
131 Bxmat.
reinit(Dim, 1, fe.n_basis());
134 for (
uint_t i=0; i<fe.n_q_points(); i++) {
138 _property->value(c, mat);
140 <
scalar_t, scalar_t, FEVarType, Dim>(*_fe_var_data, i, grad, Bxmat);
141 Bxmat.vector_mult_transpose(vec, grad);
142 res += fe.detJxW(i) * mat * vec;
146 Bxmat.right_multiply_transpose(mat2, Bxmat);
147 (*jac) += fe.detJxW(i) * mat * mat2;
159 template <
typename ScalarFieldType>
161 const ScalarFieldType& f,
165 Assert0(_fe_var_data,
"FE data not initialized.");
166 Assert0(_property,
"Section property not initialized");
168 const typename FEVarType::fe_shape_deriv_t
169 &fe = _fe_var_data->get_fe_shape_data();
171 typename Eigen::Matrix<scalar_t, Dim, 1>
174 vec = vector_t::Zero(fe.n_basis());
176 typename SectionPropertyType::value_t
179 mat2 = matrix_t::Zero(fe.n_basis(), fe.n_basis());
183 Bxmat.
reinit(Dim, 1, fe.n_basis());
186 for (
uint_t i=0; i<fe.n_q_points(); i++) {
190 _property->derivative(c, f, mat);
192 <
scalar_t, scalar_t, FEVarType, Dim>(*_fe_var_data, i, grad, Bxmat);
193 Bxmat.vector_mult_transpose(vec, grad);
194 res += fe.detJxW(i) * mat * vec;
198 Bxmat.right_multiply_transpose(mat2, Bxmat);
199 (*jac) += fe.detJxW(i) * mat * mat2;
217 #endif // __mast_conduction_kernel_h__ typename Eigen::Matrix< scalar_t, Eigen::Dynamic, 1 > vector_t
typename FEVarType::fe_shape_deriv_t fe_shape_deriv_t
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...
const SectionPropertyType * _property
typename FEVarType::fe_shape_deriv_t::scalar_t basis_scalar_t
typename Eigen::Matrix< scalar_t, Eigen::Dynamic, Eigen::Dynamic > matrix_t
#define Assert0(cond, msg)
void set_section_property(const SectionPropertyType &p)
typename FEVarType::scalar_t scalar_t
void compute(ContextType &c, vector_t &res, matrix_t *jac=nullptr) const
Computes the residual of variational term and returns it in res, and its Jacobian in jac if it is no...
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...
const FEVarType * _fe_var_data
virtual ~ConductionKernel()
void gradient_operator(const FEVarType &fe_var, const uint_t qp, typename Eigen::Matrix< VarScalarType, Dim, 1 > &epsilon, MAST::Numerics::FEMOperatorMatrix< NodalScalarType > &Bmat)
void set_fe_var_data(const FEVarType &fe_data)