MAST
utility.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 #ifndef __mast_utility__
21 #define __mast_utility__
22 
23 
24 // MAST includes
25 #include "base/mast_data_types.h"
26 
27 // libMesh includes
28 #include "libmesh/parallel.h"
29 
30 
31 namespace MAST {
32 
33  template <typename ValType>
34  void transform_to_elem_vector(libMesh::DenseVector<ValType>& v,
35  const DenseRealVector& v_real);
36 
37 
38  template <>
40  const DenseRealVector& v_real) {
41  // make sure that the real vector is twice the size of the dense vector
42  const unsigned int n = v.size();
43  libmesh_assert_equal_to(v_real.size(), n);
44  v = v_real;
45  }
46 
47 
48  template <>
50  const DenseRealVector& v_real) {
51  // make sure that the real vector is twice the size of the dense vector
52  const unsigned int n = v.size();
53  libmesh_assert_equal_to(v_real.size(), 2*n);
54 
55  for (unsigned int i=0; i<n; i++)
56  v(i) = Complex(v_real(i), v_real(i+n));
57  }
58 
59 
60 
61  template <typename ValType>
62  void transform_to_elem_matrix(libMesh::DenseMatrix<ValType>& m,
63  const DenseRealMatrix& m_real);
64 
65 
66  template <>
68  const DenseRealMatrix& m_real) {
69  // make sure that the real vector is twice the size of the dense vector
70  const unsigned int mm = m.m(), nn = m.n();
71  libmesh_assert_equal_to(m_real.m(), mm);
72  libmesh_assert_equal_to(m_real.n(), nn);
73  m = m_real;
74  }
75 
76 
77  template <>
79  const DenseRealMatrix& m_real) {
80  // make sure that the real vector is twice the size of the dense vector
81  const unsigned int mm = m.m(), nn = m.n();
82  libmesh_assert_equal_to(m_real.m(), 2*mm);
83  libmesh_assert_equal_to(m_real.n(), 2*nn);
84 
85  for (unsigned int i=0; i<mm; i++)
86  for (unsigned int j=0; j<nn; j++)
87  m(i,j) = Complex(m_real(i,j), -m_real(i,j+nn));
88  }
89 
90 
97  template <typename ValType>
98  void add_to_assembled_vector(RealVectorX& assembled_vec,
99  const ValType& elem_vec);
100 
101 
108  template <typename ValType>
109  inline void add_to_assembled_matrix(RealMatrixX& assembled_mat,
110  const ValType& elem_mat);
111 
112 
113  template <>
114  inline void
116  const RealMatrixX& elem_mat) {
117  assembled_mat += elem_mat;
118  }
119 
120 
121  template <>
122  inline void
124  const RealVectorX& elem_vec) {
125  assembled_vec += elem_vec;
126  }
127 
128 
129  template <>
130  inline void
132  const ComplexMatrixX& elem_mat) {
133 
134  // make sure the the assembled mat is twice the size of the elem mat
135  const unsigned int
136  m = (unsigned int)elem_mat.rows(),
137  n = (unsigned int)elem_mat.cols();
138  libmesh_assert_equal_to(assembled_mat.rows(), m*2);
139  libmesh_assert_equal_to(assembled_mat.cols(), n*2);
140  for (unsigned int i=0; i<m; i++)
141  for (unsigned int j=0; j<n; j++) {
142  assembled_mat(i,j) += std::real(elem_mat(i,j));
143  assembled_mat(i+m,j+n) += std::real(elem_mat(i,j));
144  assembled_mat(i,j+n) += -std::imag(elem_mat(i,j));
145  assembled_mat(i+m,j) += std::imag(elem_mat(i,j));
146  }
147  }
148 
149 
150  template <>
151  inline void
153  const ComplexVectorX& elem_vec) {
154  // make sure the the assembled mat is twice the size of the elem mat
155  const unsigned int n = (unsigned int)elem_vec.size();
156  libmesh_assert_equal_to(assembled_vec.size(), n*2);
157 
158  for (unsigned int i=0; i<n; i++) {
159  assembled_vec(i) += std::real(elem_vec(i));
160  assembled_vec(i+n) += std::imag(elem_vec(i));
161  }
162  }
163 
164 
165 
166  inline void
167  copy (DenseRealMatrix& m1, const RealMatrixX& m2) {
168 
169  const unsigned int m=(unsigned int)m2.rows(), n=(unsigned int)m2.cols();
170  m1.resize(m, n);
171  for (unsigned int i=0; i<m; i++)
172  for (unsigned int j=0; j<n; j++)
173  m1(i,j) = m2(i,j);
174  }
175 
176 
177 
178  inline void
179  copy (RealMatrixX& m2, const DenseRealMatrix& m1) {
180 
181  const unsigned int m=(unsigned int)m1.m(), n=(unsigned int)m1.n();
182  m2.setZero(m, n);
183  for (unsigned int i=0; i<m; i++)
184  for (unsigned int j=0; j<n; j++)
185  m2(i,j) = m1(i,j);
186  }
187 
188 
189 
190  inline void
191  copy (DenseRealVector& v1, const RealVectorX& v2) {
192 
193  const unsigned int m=(unsigned int)v2.rows();
194  v1.resize(m);
195  for (unsigned int i=0; i<m; i++)
196  v1(i) = v2(i);
197  }
198 
199 
200  inline void
201  copy (RealVectorX& v1, const DenseRealVector& v2) {
202 
203  const unsigned int m=(unsigned int)v2.size();
204 
205  v1 = RealVectorX::Zero(m);
206 
207  for (unsigned int i=0; i<m; i++)
208  v1(i) = v2(i);
209  }
210 
211 
212 
213  inline void
214  parallel_sum (const libMesh::Parallel::Communicator& c,
215  RealMatrixX& mat) {
216 
217  const unsigned int
218  m = (unsigned int)mat.rows(),
219  n = (unsigned int)mat.cols();
220 
221  std::vector<Real> vals(m*n);
222  for (unsigned int i=0; i<m; i++)
223  for (unsigned int j=0; j<n; j++)
224  vals[m*j+i] = mat(i,j);
225 
226  c.sum(vals);
227 
228  for (unsigned int i=0; i<m; i++)
229  for (unsigned int j=0; j<n; j++)
230  mat(i,j) = vals[m*j+i];
231  }
232 
233 
234 
235  inline void
236  parallel_sum (const libMesh::Parallel::Communicator& c,
237  ComplexMatrixX& mat) {
238 
239  const unsigned int
240  m = (unsigned int)mat.rows(),
241  n = (unsigned int)mat.cols();
242 
243  std::vector<Complex> vals(m*n);
244  for (unsigned int i=0; i<m; i++)
245  for (unsigned int j=0; j<n; j++)
246  vals[m*j+i] = mat(i,j);
247 
248  c.sum(vals);
249 
250  for (unsigned int i=0; i<m; i++)
251  for (unsigned int j=0; j<n; j++)
252  mat(i,j) = vals[m*j+i];
253 
254  }
255 }
256 
257 
258 
259 
260 #endif // __mast_utility__
261 
libMesh::DenseMatrix< Real > DenseRealMatrix
void add_to_assembled_matrix(RealMatrixX &assembled_mat, const ValType &elem_mat)
All calculations in MAST are done using Real numbers.
void transform_to_elem_vector(libMesh::DenseVector< ValType > &v, const DenseRealVector &v_real)
Matrix< Complex, Dynamic, 1 > ComplexVectorX
libMesh::Complex Complex
libMesh::DenseVector< Real > DenseRealVector
void transform_to_elem_matrix(libMesh::DenseMatrix< ValType > &m, const DenseRealMatrix &m_real)
Matrix< Real, Dynamic, Dynamic > RealMatrixX
void parallel_sum(const libMesh::Parallel::Communicator &c, RealMatrixX &mat)
Definition: utility.h:214
Matrix< Complex, Dynamic, Dynamic > ComplexMatrixX
void copy(DenseRealMatrix &m1, const RealMatrixX &m2)
Definition: utility.h:167
libMesh::DenseMatrix< Complex > DenseComplexMatrix
Matrix< Real, Dynamic, 1 > RealVectorX
libMesh::DenseVector< Complex > DenseComplexVector
void add_to_assembled_vector(RealVectorX &assembled_vec, const ValType &elem_vec)
All calculations in MAST are done using Real numbers.