MAST3
Multidisciplinary-design Adaptation and Sensitivity Toolkit (MAST)
linear_acceleration.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_elasticity_continuum_linear_acceleration_kernel_h__
21 #define __mast_elasticity_continuum_linear_acceleration_kernel_h__
22 
23 // MAST includes
25 #include <mast/base/exceptions.hpp>
27 
28 
29 namespace MAST {
30 namespace Physics {
31 namespace Elasticity {
32 namespace LinearContinuum {
33 
34 
35 template <typename NodalScalarType,
36  typename VarScalarType,
37  typename FEVarType,
38  uint_t Dim>
39 inline void
40 displacement(const FEVarType& fe_var,
41  const uint_t qp,
42  typename Eigen::Matrix<VarScalarType, Dim, 1>& u,
44 
45  u.setZero();
46 
47  const typename FEVarType::fe_shape_deriv_t
48  &fe = fe_var.get_fe_shape_data();
49 
50  // make sure all matrices are the right size
51  Assert1(Bmat.m() == Dim,
52  Bmat.m(),
53  "Incompatible operator size");
54  Assert2(Bmat.n() == Dim*fe.n_basis(),
55  Bmat.n(), Dim*fe.n_basis(),
56  "Incompatible Operator size.");
57 
58 
59  // set the strain displacement relation
60  for (uint_t i=0; i<Dim; i++) {
61  Bmat.set_shape_function(i, i, fe.phi(qp));
62  u(i) = fe_var.u(qp, i);
63  }
64 }
65 
81 template <typename FEVarType,
82  typename DensityFieldType,
83  typename SectionAreaType,
84  uint_t Dim,
85  typename ContextType>
87 
88 public:
89 
90  using scalar_t = typename FEVarType::scalar_t;
91  using basis_scalar_t = typename FEVarType::fe_shape_deriv_t::scalar_t;
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>;
94 
96  _section (nullptr),
97  _density (nullptr),
98  _fe_var_data (nullptr)
99  { }
100 
101  virtual ~LinearAcceleration() { }
102 
103  inline void set_section_area(const SectionAreaType& s) { _section = &s;}
104 
105  inline void set_density(const DensityFieldType& rho) { _density = &rho;}
106 
107  inline void set_fe_var_data(const FEVarType& fe) { _fe_var_data = &fe;}
108 
109  inline uint_t n_dofs() const {
110 
111  Assert0(_fe_var_data, "FE data not initialized.");
112 
113  return Dim*_fe_var_data->get_fe_shape_data().n_basis();
114  }
115 
123  inline void compute(ContextType& c,
124  vector_t& res,
125  matrix_t* jac = nullptr) const {
126 
127  Assert0(_fe_var_data, "FE data not initialized.");
128  Assert0(_section, "Section property not initialized");
129  Assert0(_density, "Density not initialized");
130 
131  const typename FEVarType::fe_shape_deriv_t
132  &fe = _fe_var_data->get_fe_shape_data();
133 
134  typename Eigen::Matrix<scalar_t, Dim, 1>
135  u;
136  vector_t
137  vec = vector_t::Zero(Dim*fe.n_basis());
138 
139  scalar_t
140  rho = 0.,
141  sec = 0.;
142 
143  matrix_t
144  mat = matrix_t::Zero(Dim*fe.n_basis(), Dim*fe.n_basis());
145 
147  Bmat;
148  Bmat.reinit(Dim, Dim, fe.n_basis());
149 
150  for (uint_t i=0; i<fe.n_q_points(); i++) {
151 
152  c.qp = i;
153 
154  _density->value(c, rho);
155  _section->value(c, sec);
156 
158  <scalar_t, scalar_t, FEVarType, Dim>(*_fe_var_data, i, u, Bmat);
159  Bmat.vector_mult_transpose(vec, u);
160  res += fe.detJxW(i) * vec * rho * sec;
161 
162  if (jac) {
163 
164  Bmat.right_multiply_transpose(mat, Bmat);
165  (*jac) += fe.detJxW(i) * mat * rho * sec;
166  }
167  }
168  }
169 
170 
179  template <typename ScalarFieldType>
180  inline void derivative(ContextType& c,
181  const ScalarFieldType& f,
182  vector_t& res,
183  matrix_t* jac = nullptr) const {
184 
185  Assert0(_fe_var_data, "FE data not initialized.");
186  Assert0(_section, "Section property not initialized");
187  Assert0(_density, "Density not initialized");
188 
189  const typename FEVarType::fe_shape_deriv_t
190  &fe = _fe_var_data->get_fe_shape_data();
191 
192  typename Eigen::Matrix<scalar_t, Dim, 1>
193  u;
194  vector_t
195  vec = vector_t::Zero(Dim*fe.n_basis());
196 
197  scalar_t
198  rho = 0.,
199  drho = 0.,
200  sec = 0.,
201  dsec = 0.;
202 
203  matrix_t
204  mat = matrix_t::Zero(Dim*fe.n_basis(), Dim*fe.n_basis());
205 
207  Bmat;
208  Bmat.reinit(Dim, Dim, fe.n_basis());
209 
210  for (uint_t i=0; i<fe.n_q_points(); i++) {
211 
212  c.qp = i;
213 
214  _density->value(c, rho);
215  _density->derivative(c, f, drho);
216  _section->value(c, sec);
217  _section->derivative(c, f, dsec);
218 
220  <scalar_t, scalar_t, FEVarType, Dim>(*_fe_var_data, i, u, Bmat);
221  Bmat.vector_mult_transpose(vec, u);
222  res += fe.detJxW(i) * vec * (drho * sec + rho * dsec);
223 
224  if (jac) {
225 
226  Bmat.right_multiply_transpose(mat, Bmat);
227  (*jac) += fe.detJxW(i) * mat * (drho * sec + rho * dsec);
228  }
229  }
230  }
231 
232 private:
233 
234  const SectionAreaType *_section;
235  const DensityFieldType *_density;
236  const FEVarType *_fe_var_data;
237 };
238 
239 
240 } // namespace LinearContinuum
241 } // namespace Elasticity
242 } // namespace Physics
243 } // namespace MAST
244 
245 
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...
#define Assert1(cond, v1, msg)
Definition: exceptions.hpp:143
This class implements the discrete evaluation of the acceleration kernel defined as where...
void compute(ContextType &c, vector_t &res, matrix_t *jac=nullptr) const
Computes the residual of variational term and returns it in res.
void displacement(const FEVarType &fe_var, const uint_t qp, typename Eigen::Matrix< VarScalarType, Dim, 1 > &u, MAST::Numerics::FEMOperatorMatrix< NodalScalarType > &Bmat)
#define Assert0(cond, msg)
Definition: exceptions.hpp:134
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...
unsigned int uint_t
#define Assert2(cond, v1, v2, msg)
Definition: exceptions.hpp:152
typename Eigen::Matrix< scalar_t, Eigen::Dynamic, 1 > vector_t
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