MAST3
Multidisciplinary-design Adaptation and Sensitivity Toolkit (MAST)
residual_sensitivity.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_libmesh_residual_sensitivity_h__
21 #define __mast_libmesh_residual_sensitivity_h__
22 
23 // MAST includes
25 #include <mast/base/exceptions.hpp>
29 
30 // libMesh includes
31 #include <libmesh/nonlinear_implicit_system.h>
32 #include <libmesh/dof_map.h>
33 
34 
35 namespace MAST {
36 namespace Base {
37 namespace Assembly {
38 namespace libMeshWrapper {
39 
40 template <typename ScalarType,
41  typename ElemOpsType>
43 
44 public:
45 
47  "Scalar type of assembly and element operations must be same");
48 
50  _e_ops (nullptr)
51  { }
52 
53  virtual ~ResidualSensitivity() { }
54 
55  inline void set_elem_ops(ElemOpsType& e_ops) { _e_ops = &e_ops; }
56 
57  template <typename VecType, typename MatType, typename ContextType, typename ScalarFieldType>
58  inline void assemble(ContextType &c,
59  const ScalarFieldType& f,
60  const VecType &X,
61  VecType *R,
62  MatType *J) {
63 
64  Assert0( R || J, "Atleast one assembled quantity should be specified.");
65 
68 
69  // iterate over each element, initialize it and get the relevant
70  // analysis quantities
72  sol_accessor(*c.sys, X);
73 
74  using elem_vector_t = typename ElemOpsType::vector_t;
75  using elem_matrix_t = typename ElemOpsType::matrix_t;
76 
77  elem_vector_t res_e;
78  elem_matrix_t jac_e;
79 
80 
81  libMesh::MeshBase::const_element_iterator
82  el = c.mesh->active_local_elements_begin(),
83  end_el = c.mesh->active_local_elements_end();
84 
85  for ( ; el != end_el; ++el) {
86 
87  // set element in the context, which will be used for the initialization routines
88  c.elem = *el;
89 
90  sol_accessor.init(*c.elem);
91 
92  res_e.setZero(sol_accessor.n_dofs());
93  if (J) jac_e.setZero(sol_accessor.n_dofs(), sol_accessor.n_dofs());
94 
95  // perform the element level calculations
96  _e_ops->derivative(c, f, sol_accessor, res_e, J?&jac_e:nullptr);
97 
98  // constrain the quantities to account for hanging dofs,
99  // Dirichlet constraints, etc.
100  if (R && J)
102  <ScalarType, VecType, MatType, elem_vector_t, elem_matrix_t>
103  (*R, *J, c.sys->get_dof_map(), sol_accessor.dof_indices(), res_e, jac_e);
104  else if (R)
106  <ScalarType, VecType, elem_vector_t>
107  (*R, c.sys->get_dof_map(), sol_accessor.dof_indices(), res_e);
108  else
110  <ScalarType, MatType, elem_matrix_t>
111  (*J, c.sys->get_dof_map(), sol_accessor.dof_indices(), jac_e);
112  }
113 
114  // parallel matrix/vector require finalization of communication
117  }
118 
119 private:
120 
122  ElemOpsType *_e_ops;
123 };
124 
125 } // namespace libMeshWrapper
126 } // namespace Assembly
127 } // namespace Base
128 } // namespace MAST
129 
130 #endif // __mast_libmesh_residual_sensitivity_h__
void init(const libMesh::Elem &e)
Definition: accessor.hpp:69
void setZero(ValType &m)
Definition: utility.hpp:44
void assemble(ContextType &c, const ScalarFieldType &f, const VecType &X, VecType *R, MatType *J)
void finalize(ValType &m)
Definition: utility.hpp:252
const std::vector< libMesh::dof_id_type > & dof_indices() const
Definition: accessor.hpp:55
std::enable_if< Dim< 3, ScalarType >::typesource_load_multiplier(const SourceLoadFieldType *f, const SectionAreaType *s, ContextType &c) { Assert0(f, "Invalid pointer");Assert0(s, "Invalid pointer");return f-> value(c) *s -> value(c)
std::enable_if< std::is_same< ScalarType, real_t >::value, void >::type constrain_and_add_matrix_and_vector(VecType &v, MatType &m, const libMesh::DofMap &dof_map, std::vector< libMesh::dof_id_type > &dof_indices, SubVecType &v_sub, SubMatType &m_sub)
Definition: utility.hpp:87
std::enable_if< std::is_same< ScalarType, real_t >::value, void >::type constrain_and_add_vector(VecType &v, const libMesh::DofMap &dof_map, std::vector< libMesh::dof_id_type > &dof_indices, SubVecType &v_sub)
Definition: utility.hpp:113
#define Assert0(cond, msg)
Definition: exceptions.hpp:134
std::enable_if< std::is_same< ScalarType, real_t >::value, void >::type constrain_and_add_matrix(MatType &m, const libMesh::DofMap &dof_map, std::vector< libMesh::dof_id_type > &dof_indices, SubMatType &m_sub)
Definition: utility.hpp:132