MAST3
Multidisciplinary-design Adaptation and Sensitivity Toolkit (MAST)
linear_strain_energy.hpp
Go to the documentation of this file.
1 /*
2 * MAST: Multidisciplinary-design Adaptation and Sensitivity Toolkit
3 * Copyright (C) 2013-2020 Manav Bhatia and MAST authors
4 *
5 * This library is free software; you can redistribute it and/or
6 * modify it under the terms of the GNU Lesser General Public
7 * License as published by the Free Software Foundation; either
8 * version 2.1 of the License, or (at your option) any later version.
9 *
10 * This library is distributed in the hope that it will be useful,
11 * but WITHOUT ANY WARRANTY; without even the implied warranty of
12 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
13 * Lesser General Public License for more details.
14 *
15 * You should have received a copy of the GNU Lesser General Public
16 * License along with this library; if not, write to the Free Software
17 * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
18 */
19 
20 #ifndef __mast_linear_continuum_strain_energy_h__
21 #define __mast_linear_continuum_strain_energy_h__
22 
23 // MAST includes
25 
26 namespace MAST {
27 namespace Physics {
28 namespace Elasticity {
29 namespace LinearContinuum {
30 
31 
32 template <typename FEVarType,
33  typename SectionPropertyType,
34  uint_t Dim,
35  typename ContextType>
36 class StrainEnergy {
37 
38 public:
39 
40  using scalar_t = typename FEVarType::scalar_t;
41  using basis_scalar_t = typename FEVarType::fe_shape_deriv_t::scalar_t;
42  using vector_t = typename Eigen::Matrix<scalar_t, Eigen::Dynamic, 1>;
43  using matrix_t = typename Eigen::Matrix<scalar_t, Eigen::Dynamic, Eigen::Dynamic>;
44  using fe_shape_deriv_t = typename FEVarType::fe_shape_deriv_t;
45  static const uint_t
47 
49  _property (nullptr),
50  _fe_var_data (nullptr)
51  { }
52 
53  virtual ~StrainEnergy() { }
54 
55  inline void
56  set_section_property(const SectionPropertyType& p) {
57 
58  Assert0(!_property, "Property already initialized.");
59 
60  _property = &p;
61  }
62 
63  inline void set_fe_var_data(const FEVarType& fe_data)
64  {
65  Assert0(!_fe_var_data, "FE data already initialized.");
66  _fe_var_data = &fe_data;
67  }
68 
69  inline uint_t n_dofs() const {
70 
71  Assert0(_fe_var_data, "FE data not initialized.");
72 
73  return Dim*_fe_var_data->get_fe_shape_data().n_basis();
74  }
75 
76  inline void compute(ContextType& c,
77  vector_t& res,
78  matrix_t* jac = nullptr) const {
79 
80  Assert0(_fe_var_data, "FE data not initialized.");
81  Assert0(_property, "Section property not initialized");
82 
83  const typename FEVarType::fe_shape_deriv_t
84  &fe = _fe_var_data->get_fe_shape_data();
85 
86  typename Eigen::Matrix<scalar_t, n_strain, 1>
87  epsilon,
88  stress;
89  vector_t
90  vec = vector_t::Zero(Dim*fe.n_basis());
91 
92  typename SectionPropertyType::value_t
93  mat;
94 
95  matrix_t
96  mat1 = matrix_t::Zero(n_strain, Dim*fe.n_basis()),
97  mat2 = matrix_t::Zero(Dim*fe.n_basis(), Dim*fe.n_basis());
98 
100  Bxmat;
101  Bxmat.reinit(n_strain, Dim, fe.n_basis());
102 
103 
104  for (uint_t i=0; i<fe.n_q_points(); i++) {
105 
106  c.qp = i;
107 
108  _property->value(c, mat);
110  <scalar_t, scalar_t, FEVarType, Dim>(*_fe_var_data, i, epsilon, Bxmat);
111  stress = mat * epsilon;
112  Bxmat.vector_mult_transpose(vec, stress);
113  res += fe.detJxW(i) * vec;
114 
115  if (jac) {
116 
117  Bxmat.left_multiply(mat1, mat);
118  Bxmat.right_multiply_transpose(mat2, mat1);
119  (*jac) += fe.detJxW(i) * mat2;
120  }
121  }
122  }
123 
124  template <typename ScalarFieldType>
125  inline void derivative(ContextType& c,
126  const ScalarFieldType& f,
127  vector_t& res,
128  matrix_t* jac = nullptr) const {
129 
130  Assert0(_fe_var_data, "FE data not initialized.");
131  Assert0(_property, "Section property not initialized");
132 
133  const typename FEVarType::fe_shape_deriv_t
134  &fe = _fe_var_data->get_fe_shape_data();
135 
136  typename Eigen::Matrix<scalar_t, n_strain, 1>
137  epsilon,
138  stress;
139  vector_t
140  vec = vector_t::Zero(Dim*fe.n_basis());
141 
142  typename SectionPropertyType::value_t
143  mat;
144  matrix_t
145  mat1 = matrix_t::Zero(n_strain, Dim*fe.n_basis()),
146  mat2 = matrix_t::Zero(Dim*fe.n_basis(), Dim*fe.n_basis());
147 
149  Bxmat;
150  Bxmat.reinit(n_strain, Dim, fe.n_basis());
151 
152 
153  for (uint_t i=0; i<fe.n_q_points(); i++) {
154 
155  c.qp = i;
156 
157  _property->derivative(c, f, mat);
159  <scalar_t, scalar_t, FEVarType, Dim>(*_fe_var_data, i, epsilon, Bxmat);
160  stress = mat * epsilon;
161  Bxmat.vector_mult_transpose(vec, stress);
162  res += fe.detJxW(i) * vec;
163 
164  if (jac) {
165 
166  Bxmat.left_multiply(mat1, mat);
167  Bxmat.right_multiply_transpose(mat2, mat1);
168  (*jac) += fe.detJxW(i) * mat2;
169  }
170  }
171  }
172 
173 
174 private:
175 
176 
177  const SectionPropertyType *_property;
178  const FEVarType *_fe_var_data;
179 };
180 
181 } // namespace LinearContinuum
182 } // namespace Elasticity
183 } // namespace Physics
184 } // namespace MAST
185 
186 #endif // __mast_linear_continuum_strain_energy_h__
void derivative(ContextType &c, const ScalarFieldType &f, vector_t &res, matrix_t *jac=nullptr) const
typename Eigen::Matrix< scalar_t, Eigen::Dynamic, Eigen::Dynamic > matrix_t
void compute(ContextType &c, vector_t &res, matrix_t *jac=nullptr) const
typename FEVarType::fe_shape_deriv_t::scalar_t basis_scalar_t
std::enable_if< Dim==2, void >::type strain(const FEVarType &fe_var, const uint_t qp, typename Eigen::Matrix< VarScalarType, 3, 1 > &epsilon, MAST::Numerics::FEMOperatorMatrix< NodalScalarType > &Bmat)
#define Assert0(cond, msg)
Definition: exceptions.hpp:134
unsigned int uint_t
typename Eigen::Matrix< scalar_t, Eigen::Dynamic, 1 > vector_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...