MAST3
Multidisciplinary-design Adaptation and Sensitivity Toolkit (MAST)
source_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_source_load_h__
21 #define __mast_conduction_source_load_h__
22 
23 // MAST includes
25 #include <mast/base/exceptions.hpp>
26 
27 namespace MAST {
28 namespace Physics {
29 namespace Conduction {
30 
31 template <typename ScalarType,
32  typename SectionAreaType,
33  typename SourceLoadFieldType,
34  typename ContextType,
35  uint_t Dim>
36 typename std::enable_if<Dim<3, ScalarType>::type
37 source_load_multiplier(const SourceLoadFieldType *f,
38  const SectionAreaType *s,
39  ContextType &c) {
40  Assert0(f, "Invalid pointer");
41  Assert0(s, "Invalid pointer");
42  return f->value(c) * s->value(c);
43 }
44 
45 template <typename ScalarType,
46  typename SectionAreaType,
47  typename SourceLoadFieldType,
48  typename ContextType,
49  uint_t Dim>
50 typename std::enable_if<Dim==3, ScalarType>::type
51 source_load_multiplier(const SourceLoadFieldType *f,
52  const SectionAreaType *s,
53  ContextType &c) {
54  Assert0(f, "Invalid pointer");
55  Assert0(!s, "Pointer must be nullptr");
56  return f->value(c);
57 }
58 
59 template <typename ScalarType,
60  typename SectionAreaType,
61  typename SourceLoadFieldType,
62  typename ContextType,
63  typename ScalarFieldType,
64  uint_t Dim>
65 typename std::enable_if<Dim<3, ScalarType>::type
66 source_load_derivative_multiplier(const SourceLoadFieldType *f,
67  const SectionAreaType *s,
68  ContextType &c,
69  const ScalarFieldType &p) {
70 
71  Assert0(f, "Invalid pointer");
72  Assert0(s, "Invalid pointer");
73  return (f->value(c) * s->derivative(c, p) +
74  s->value(c) * f->derivative(c, p));
75 }
76 
77 template <typename ScalarType,
78  typename SectionAreaType,
79  typename SourceLoadFieldType,
80  typename ContextType,
81  typename ScalarFieldType,
82  uint_t Dim>
83 typename std::enable_if<Dim==3, ScalarType>::type
84 source_load_derivative_multiplier(const SourceLoadFieldType *f,
85  const SectionAreaType *s,
86  ContextType &c,
87  const ScalarFieldType &p) {
88 
89  Assert0(f, "Invalid pointer");
90  Assert0(!s, "Pointer must be nullptr");
91  return f->derivative(c);
92 }
93 
94 
110 template <typename FEVarType,
111  typename SourceLoadFieldType,
112  typename SectionAreaType,
113  uint_t Dim,
114  typename ContextType>
116 
117 public:
118 
119  using scalar_t = typename FEVarType::scalar_t;
120  using basis_scalar_t = typename FEVarType::fe_shape_deriv_t::scalar_t;
121  using vector_t = typename Eigen::Matrix<scalar_t, Eigen::Dynamic, 1>;
122  using matrix_t = typename Eigen::Matrix<scalar_t, Eigen::Dynamic, Eigen::Dynamic>;
123 
125  _section (nullptr),
126  _load (nullptr),
127  _fe_var_data (nullptr)
128  { }
129 
130  virtual ~SourceHeatLoad() { }
131 
135  inline void set_section_area(const SectionAreaType& s) {
136 
137  Assert1(Dim < 3, Dim, "SectionAreaType only used for 1D and 2D elements");
138 
139  _section = &s;
140  }
141 
142  inline void set_source(const SourceLoadFieldType& s) { _load = &s;}
143 
144  inline void set_fe_var_data(const FEVarType& fe) { _fe_var_data = &fe;}
145 
146  inline uint_t n_dofs() const {
147 
148  Assert0(_fe_var_data, "FE data not initialized.");
149 
150  return _fe_var_data->get_fe_shape_data().n_basis();
151  }
152 
160  inline void
161  compute(ContextType& c,
162  vector_t& res,
163  matrix_t* jac = nullptr) const {
164 
165  Assert0(_fe_var_data, "FE data not initialized.");
166  Assert0(Dim==3 || _section, "Section property not initialized");
167  Assert0(_load, "Source load not initialized");
168 
169  const typename FEVarType::fe_shape_deriv_t
170  &fe = _fe_var_data->get_fe_shape_data();
171 
172  for (uint_t i=0; i<fe.n_q_points(); i++) {
173 
174  c.qp = i;
176  SectionAreaType,
177  SourceLoadFieldType,
178  ContextType,
179  Dim>(_load, _section, c);
180 
181  for (uint_t k=0; k<fe.n_basis(); k++)
182  res(k) -= fe.detJxW(i) * fe.phi(i, k) * p;
183  }
184  }
185 
186 
194  template <typename ScalarFieldType>
195  inline void derivative(ContextType& c,
196  const ScalarFieldType& f,
197  vector_t& res,
198  matrix_t* jac = nullptr) const {
199 
200  Assert0(_fe_var_data, "FE data not initialized.");
201  Assert0(Dim==3 || _section, "Section property not initialized");
202  Assert0(_load, "Source load not initialized");
203 
204  const typename FEVarType::fe_shape_deriv_t
205  &fe = _fe_var_data->get_fe_shape_data();
206 
207  for (uint_t i=0; i<fe.n_q_points(); i++) {
208 
209  c.qp = i;
210  scalar_t p =
211  source_load_derivative_multiplier<scalar_t,
212  SectionAreaType,
213  SourceLoadFieldType,
214  ContextType,
215  ScalarFieldType,
216  Dim>(_load, _section, c, f);
217 
218  for (uint_t k=0; k<fe.n_basis(); k++)
219  res(k) -= fe.detJxW(i) * fe.phi(i, k) * p;
220  }
221  }
222 
223 private:
224 
225  const SectionAreaType *_section;
226  const SourceLoadFieldType *_load;
227  const FEVarType *_fe_var_data;
228 };
229 
230 } // namespace Conduction
231 } // namespace Physics
232 } // namespace MAST
233 
234 
235 #endif // __mast_conduction_source_load_h__
This class implements the discrete evaluation of the conduction (Laplace operator) kernel defined as ...
typename FEVarType::scalar_t scalar_t
typename Eigen::Matrix< scalar_t, Eigen::Dynamic, Eigen::Dynamic > matrix_t
#define Assert1(cond, v1, msg)
Definition: exceptions.hpp:143
typename Eigen::Matrix< scalar_t, Eigen::Dynamic, 1 > vector_t
typename FEVarType::fe_shape_deriv_t::scalar_t basis_scalar_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...
#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.
void set_fe_var_data(const FEVarType &fe)
void set_section_area(const SectionAreaType &s)
Provides the section area through the object s.
void set_source(const SourceLoadFieldType &s)
std::enable_if< Dim==3, ScalarType >::type source_load_multiplier(const SourceLoadFieldType *f, const SectionAreaType *s, ContextType &c)