MAST3
Multidisciplinary-design Adaptation and Sensitivity Toolkit (MAST)
gradient_operator.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_gradient_operator_h__
21 #define __mast_conduction_gradient_operator_h__
22 
23 // MAST includes
25 #include <mast/base/exceptions.hpp>
27 
28 
29 namespace MAST {
30 namespace Physics {
31 namespace Conduction {
32 namespace GradientOperator {
33 
34 
35 template <typename NodalScalarType,
36  typename VarScalarType,
37  typename FEVarType,
38  uint_t Dim>
39 inline void
40 gradient_operator(const FEVarType &fe_var,
41  const uint_t qp,
42  typename Eigen::Matrix<VarScalarType, Dim, 1> &epsilon,
44 
45  epsilon.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(epsilon.size() == Dim,
52  epsilon.size(),
53  "Invalid gradient dimension");
54  Assert1(Bmat.m() == Dim,
55  Bmat.m(),
56  "Invalid gradient dimension");
57  Assert2(Bmat.n() == fe.n_basis(),
58  Bmat.n(), fe.n_basis(),
59  "Incompatible Operator size.");
60 
61 
62  for (uint_t i=0; i<Dim; i++) {
63 
64  Bmat.set_shape_function(i, 0, fe.dphi_dx(qp, i)); // du/dxi
65  epsilon(i) = fe_var.du_dx(qp, 0, i); // du/dxi
66  }
67 }
68 
69 } // namespace GradientOperator
70 } // namespace Conduction
71 } // namespace Physics
72 } // namespace MAST
73 
74 
75 #endif // __mast_linear_elastic_strain_operator_h__
#define Assert1(cond, v1, msg)
Definition: exceptions.hpp:143
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
void gradient_operator(const FEVarType &fe_var, const uint_t qp, typename Eigen::Matrix< VarScalarType, Dim, 1 > &epsilon, MAST::Numerics::FEMOperatorMatrix< NodalScalarType > &Bmat)