MAST
basis_matrix.h
Go to the documentation of this file.
1 /*
2  * MAST: Multidisciplinary-design Adaptation and Sensitivity Toolkit
3  * Copyright (C) 2013-2019 Manav Bhatia
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 
21 #ifndef __mast__structural_basis_matrix_h__
22 #define __mast__structural_basis_matrix_h__
23 
24 // C++ includes
25 #include <vector>
26 
27 // libmesh includes
28 #include "libmesh/shell_matrix.h"
29 #include "libmesh/numeric_vector.h"
30 
31 
32 namespace MAST {
33 
34  template <typename T>
35  class BasisMatrix:
36  public libMesh::ShellMatrix<T>
37  {
38  public:
39  BasisMatrix(const libMesh::Parallel::Communicator &comm_in);
40 
41  virtual ~BasisMatrix();
42 
43 
48  virtual libMesh::numeric_index_type m () const {
49 
50  libmesh_assert(modes.size() > 0);
51  return modes[0]->size();
52  }
53 
58  virtual libMesh::numeric_index_type n () const {
59 
60  libmesh_assert(modes.size() > 0);
61  return (libMesh::numeric_index_type) modes.size();
62  }
63 
64 
69  template <typename VecType>
70  void vector_mult (libMesh::NumericVector<T>& dest,
71  const VecType& arg) const {
72 
73  libmesh_assert(modes.size() > 0);
74  libmesh_assert_equal_to(m(), dest.size());
75  libmesh_assert_equal_to(n(), arg.size());
76 
77  dest.zero();
78  for (unsigned int i=0; i<n(); i++)
79  dest.add(arg(i), *(modes[i]));
80  }
81 
86  virtual void
87  vector_mult (libMesh::NumericVector<T>& dest,
88  const libMesh::NumericVector<T>& arg) const {
89 
90  // not defined for multiplcation with libMesh::NumericVector
91  libmesh_assert(false);
92  }
93 
94 
99  template <typename VecType>
100  void
101  vector_mult_transpose (VecType& dest,
102  const libMesh::NumericVector<T>& arg) const {
103 
104  libmesh_assert(modes.size() > 0);
105  libmesh_assert_equal_to(m(), arg.size());
106  libmesh_assert_equal_to(n(), dest.size());
107 
108  dest.setZero(n());
109  for (unsigned int i=0; i<n(); i++)
110  dest(i) = arg.dot(*(modes[i]));
111  }
112 
113 
118  virtual void
119  vector_mult_transpose (libMesh::NumericVector<T>& dest,
120  const libMesh::NumericVector<T>& arg) const {
121 
122  // not defined for multiplcation with libMesh::NumericVector
123  libmesh_assert(false);
124  }
125 
129  virtual void
130  vector_mult_add (libMesh::NumericVector<T>& dest,
131  const libMesh::NumericVector<T>& arg) const {
132 
133  // not defined for multiplcation with libMesh::NumericVector
134  libmesh_assert(false);
135  }
136 
137 
141  virtual void
142  get_diagonal (libMesh::NumericVector<T>& dest) const {
143 
144  // not defined for multiplcation with libMesh::NumericVector
145  libmesh_assert(false);
146  }
147 
148 
153  virtual libMesh::NumericVector<T>&
154  basis(unsigned int i) {
155 
156  // not defined for multiplcation with libMesh::NumericVector
157  libmesh_assert(modes.size() > 0);
158  libmesh_assert_less(i, modes.size());
159 
160  return *(modes[i]);
161  }
162 
166  std::vector<libMesh::NumericVector<T>*> modes;
167  };
168 
169 }
170 
171 #endif // __mast__structural_basis_matrix_h__
virtual void vector_mult(libMesh::NumericVector< T > &dest, const libMesh::NumericVector< T > &arg) const
Multiplies the matrix with arg and stores the result in dest.
Definition: basis_matrix.h:87
virtual libMesh::numeric_index_type m() const
Definition: basis_matrix.h:48
BasisMatrix(const libMesh::Parallel::Communicator &comm_in)
void vector_mult_transpose(VecType &dest, const libMesh::NumericVector< T > &arg) const
Multiplies the matrix with arg and stores the result in dest.
Definition: basis_matrix.h:101
virtual void vector_mult_transpose(libMesh::NumericVector< T > &dest, const libMesh::NumericVector< T > &arg) const
Multiplies the transpose of matrix with arg and stores the result in dest.
Definition: basis_matrix.h:119
virtual void vector_mult_add(libMesh::NumericVector< T > &dest, const libMesh::NumericVector< T > &arg) const
Multiplies the matrix with arg and adds the result to dest.
Definition: basis_matrix.h:130
virtual libMesh::NumericVector< T > & basis(unsigned int i)
Returns the vector that defines the i^th basis vector.
Definition: basis_matrix.h:154
virtual libMesh::numeric_index_type n() const
Definition: basis_matrix.h:58
void vector_mult(libMesh::NumericVector< T > &dest, const VecType &arg) const
Multiplies the matrix with arg and stores the result in dest.
Definition: basis_matrix.h:70
virtual void get_diagonal(libMesh::NumericVector< T > &dest) const
Copies the diagonal part of the matrix into dest.
Definition: basis_matrix.h:142
virtual ~BasisMatrix()
std::vector< libMesh::NumericVector< T > * > modes
vector of modes
Definition: basis_matrix.h:166