MAST3
Multidisciplinary-design Adaptation and Sensitivity Toolkit (MAST)
linear_conduction_kernel.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_conduction_kernel_h__
21 #define __mast_conduction_kernel_h__
22 
23 // MAST includes
25 
26 namespace MAST {
27 namespace Physics {
28 namespace Conduction {
29 //namespace ConductionKernel {
30 
31 
32 template <typename FEVarType,
33  typename SectionPropertyType,
34  uint_t Dim,
35  typename ContextType,
36  bool IsotropicMaterial = SectionPropertyType::is_isotropic,
37  bool LinearMaterial = SectionPropertyType::is_linear>
38 class ConductionKernel { };
39 
40 
54 template <typename FEVarType,
55  typename SectionPropertyType,
56  uint_t Dim,
57  typename ContextType>
58 class ConductionKernel<FEVarType,
59  SectionPropertyType,
60  Dim,
61  ContextType,
62  true,
63  true> {
64 
65 public:
66 
67  using scalar_t = typename FEVarType::scalar_t;
68  using basis_scalar_t = typename FEVarType::fe_shape_deriv_t::scalar_t;
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>;
71  using fe_shape_deriv_t = typename FEVarType::fe_shape_deriv_t;
72 
74  _property (nullptr),
75  _fe_var_data (nullptr)
76  { }
77 
78  virtual ~ConductionKernel() { }
79 
80  inline void
81  set_section_property(const SectionPropertyType& p) {
82 
83  Assert0(!_property, "Property already initialized.");
84 
85  _property = &p;
86  }
87 
88  inline void set_fe_var_data(const FEVarType& fe_data)
89  {
90  Assert0(!_fe_var_data, "FE data already initialized.");
91  _fe_var_data = &fe_data;
92  }
93 
94  inline uint_t n_dofs() const {
95 
96  Assert0(_fe_var_data, "FE data not initialized.");
97 
98  return _fe_var_data->get_fe_shape_data().n_basis();
99  }
100 
108  inline void compute(ContextType& c,
109  vector_t& res,
110  matrix_t* jac = nullptr) const {
111 
112  Assert0(_fe_var_data, "FE data not initialized.");
113  Assert0(_property, "Section property not initialized");
114 
115  const typename FEVarType::fe_shape_deriv_t
116  &fe = _fe_var_data->get_fe_shape_data();
117 
118  typename Eigen::Matrix<scalar_t, Dim, 1>
119  grad;
120  vector_t
121  vec = vector_t::Zero(fe.n_basis());
122 
123  typename SectionPropertyType::value_t
124  mat;
125 
126  matrix_t
127  mat2 = matrix_t::Zero(fe.n_basis(), fe.n_basis());
128 
130  Bxmat;
131  Bxmat.reinit(Dim, 1, fe.n_basis());
132 
133 
134  for (uint_t i=0; i<fe.n_q_points(); i++) {
135 
136  c.qp = i;
137 
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;
143 
144  if (jac) {
145 
146  Bxmat.right_multiply_transpose(mat2, Bxmat);
147  (*jac) += fe.detJxW(i) * mat * mat2;
148  }
149  }
150  }
151 
159  template <typename ScalarFieldType>
160  inline void derivative(ContextType& c,
161  const ScalarFieldType& f,
162  vector_t& res,
163  matrix_t* jac = nullptr) const {
164 
165  Assert0(_fe_var_data, "FE data not initialized.");
166  Assert0(_property, "Section property not initialized");
167 
168  const typename FEVarType::fe_shape_deriv_t
169  &fe = _fe_var_data->get_fe_shape_data();
170 
171  typename Eigen::Matrix<scalar_t, Dim, 1>
172  grad;
173  vector_t
174  vec = vector_t::Zero(fe.n_basis());
175 
176  typename SectionPropertyType::value_t
177  mat;
178  matrix_t
179  mat2 = matrix_t::Zero(fe.n_basis(), fe.n_basis());
180 
182  Bxmat;
183  Bxmat.reinit(Dim, 1, fe.n_basis());
184 
185 
186  for (uint_t i=0; i<fe.n_q_points(); i++) {
187 
188  c.qp = i;
189 
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;
195 
196  if (jac) {
197 
198  Bxmat.right_multiply_transpose(mat2, Bxmat);
199  (*jac) += fe.detJxW(i) * mat * mat2;
200  }
201  }
202  }
203 
204 
205 private:
206 
207 
208  const SectionPropertyType *_property;
209  const FEVarType *_fe_var_data;
210 };
211 
212 //} // namespace ConductionKernel
213 } // namespace Conduction
214 } // namespace Physics
215 } // namespace MAST
216 
217 #endif // __mast_conduction_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 Assert0(cond, msg)
Definition: exceptions.hpp:134
unsigned int uint_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...
void gradient_operator(const FEVarType &fe_var, const uint_t qp, typename Eigen::Matrix< VarScalarType, Dim, 1 > &epsilon, MAST::Numerics::FEMOperatorMatrix< NodalScalarType > &Bmat)