SmallMatrix_impl.h
Go to the documentation of this file.
1 // This file is part of the MercuryDPM project (https://www.mercurydpm.org).
2 // Copyright (c), The MercuryDPM Developers Team. All rights reserved.
3 // License: BSD 3-Clause License; see the LICENSE file in the root directory.
4 
5 /*
6  This file forms part of hpGEM. This package has been developed over a number of years by various people at the University of Twente and a full list of contributors can be found at
7  http://hpgem.org/about-the-code/team
8 
9  This code is distributed using BSD 3-Clause License. A copy of which can found below.
10 
11 
12  Copyright (c) 2014, University of Twente
13  All rights reserved.
14 
15  Redistribution and use in source and binary forms, with or without modification, are permitted provided that the following conditions are met:
16 
17  1. Redistributions of source code must retain the above copyright notice, this list of conditions and the following disclaimer.
18 
19  2. Redistributions in binary form must reproduce the above copyright notice, this list of conditions and the following disclaimer in the documentation and/or other materials provided with the distribution.
20 
21  3. Neither the name of the copyright holder nor the names of its contributors may be used to endorse or promote products derived from this software without specific prior written permission.
22 
23  THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
24  */
25 
26 
27 //Note: This code is copied and adapted from hpGEM (see license above), version 22th of January 2016. It has been
28 //integrated into MercuryDPM at 16th of March 2017.
29 
30 
31 #include "SmallMatrix.h"
32 
33 
34 extern "C"
35 {
36 
38 void dgemv_(const char* trans, int* m, int* n, double* alpha, double* A, int* LDA, double* x, int* incx, double* beta,
39  double* y, int* incy);
40 
42 int
43 dgemm_(const char* transA, const char* transB, int* M, int* N, int* k, double* alpha, double* A, int* LDA, double* B,
44  int* LDB, double* beta, double* C, int* LDC);
45 
47 int daxpy_(int* N, double* DA, double* DX, int* INCX, double* DY, int* INCY);
48 
50 void dgetrf_(int* M, int* N, double* A, int* lda, int* IPIV, int* INFO);
51 
53 void dgetri_(int* N, double* A, int* lda, int* IPIV, double* WORK, int* lwork, int* INFO);
54 
56 void dgesv_(int* N, int* NRHS, double* A, int* lda, int* IPIV, double* B, int* LDB, int* INFO);
57 }
58 
59 
60 template<unsigned int numberOfRows, unsigned int numberOfColumns>
62 {
63  if (numberOfRows == 0)
64  {
65  logger(WARN, "Trying to multiply a vector with a matrix without any rows.");
67  }
68  if (numberOfColumns == 0)
69  {
70  logger(WARN, "Trying to multiply a vector with a matrix without any columns.");
72  }
73  int nr = numberOfRows;
74  int nc = numberOfColumns;
75 
76  int i_one = 1;
77  double d_one = 1.0;
78  double d_zero = 0.0;
79 
81 
82  logger(DEBUG, "Matrix size: % x % \n Vector size: %", nr, nc, right.size());
83 
84  dgemv_("N", &nr, &nc, &d_one, this->data(), &nr, right.data(), &i_one, &d_zero, result.data(), &i_one);
85  return result;
86 }
87 
88 template<unsigned int numberOfRows, unsigned int numberOfColumns>
91 {
92  if (numberOfRows == 0)
93  {
94  logger(WARN, "Trying to multiply a vector with a matrix without any rows.");
96  }
97  if (numberOfColumns == 0)
98  {
99  logger(WARN, "Trying to multiply a vector with a matrix without any columns.");
100  return SmallVector<numberOfRows>();
101  }
102  int nr = numberOfRows;
103  int nc = numberOfColumns;
104 
105  int i_one = 1;
106  double d_one = 1.0;
107  double d_zero = 0.0;
108 
110 
111  logger(DEBUG, "Matrix size: % x % \n Vector size: %", nr, nc, right.size());
112 
113  dgemv_("N", &nr, &nc, &d_one, (const_cast<double*>(this->data())), &nr, right.data(), &i_one, &d_zero,
114  result.data(), &i_one);
115  return result;
116 }
117 
118 template<unsigned int numberOfRows, unsigned int numberOfColumns>
119 template<unsigned int K>
122 {
123  int i = numberOfRows;
124  int j = numberOfColumns;
125  int k = K;
126 
127  if (numberOfColumns == 0)
128  {
129  logger(WARN, "Trying to multiply a matrix with a matrix without any columns.");
131  }
132  //The result of the matrix is left.numberOfRows, right.numberOfColumns()
134 
135  double d_one = 1.0;
136  double d_zero = 0.0;
137 
138  //Let the actual multiplication be done by Fortran
139  dgemm_("N", "N", &i, &k, &j, &d_one, this->data(), &i, const_cast<double*>(other.data()), &j, &d_zero, C.data(),
140  &i);
141 
142  return C;
143 }
144 
145 template<unsigned int numberOfRows, unsigned int numberOfColumns>
146 template<unsigned int K>
149 {
150  int i = numberOfRows;
151  int j = numberOfColumns;
152  int k = K;
153 
154  if (numberOfColumns == 0)
155  {
156  logger(WARN, "Trying to multiply a matrix with a matrix without any columns.");
158  }
159  //The result of the matrix is left.Nrows, right.NCols()
161 
162  double d_one = 1.0;
163  double d_zero = 0.0;
164 
165  //Let the actual multiplication be done by Fortran
166  dgemm_("N", "N", &i, &k, &j, &d_one, const_cast<double*>(this->data()), &i, const_cast<double*>(other.data()), &j,
167  &d_zero, C.data(), &i);
168 
169  return C;
170 }
171 
172 template<unsigned int numberOfRows, unsigned int numberOfColumns>
175 {
176  //blas does not support in-place multiply
177  return (*this) = (*this) * other;
178 }
179 
180 template<unsigned int numberOfRows, unsigned int numberOfColumns>
182 {
183  //copied from MiddleSizeMatrix to prevent constructing a temporary MiddleSizeMatrix
184  logger.assert_debug(numberOfColumns == numberOfRows - 1,
185  "Matrix has wrong dimensions to construct the wedge stuff vector");
187 
188  switch (numberOfRows)
189  {
190  case 2:
191  result[0] = -(*this)(1, 0);
192  result[1] = +(*this)(0, 0);
193  break;
194  case 3:
195  result[0] = (*this)(1, 0) * (*this)(2, 1) - (*this)(2, 0) * (*this)(1, 1);
196  result[1] = (*this)(0, 1) * (*this)(2, 0) - (*this)(0, 0) * (*this)(2, 1); // includes minus sign already!
197  result[2] = (*this)(0, 0) * (*this)(1, 1) - (*this)(1, 0) * (*this)(0, 1);
198  break;
199  case 4:
200  result[0] = (*this)(1, 0) * (-(*this)(2, 1) * (*this)(3, 2) + (*this)(3, 1) * (*this)(2, 2)) +
201  (*this)(2, 0) * ((*this)(1, 1) * (*this)(3, 2) - (*this)(3, 1) * (*this)(1, 2)) +
202  (*this)(3, 0) * (-(*this)(1, 1) * (*this)(2, 2) + (*this)(2, 1) * (*this)(1, 2));
203 
204  result[1] = (*this)(0, 0) * ((*this)(2, 1) * (*this)(3, 2) - (*this)(3, 1) * (*this)(2, 2)) +
205  (*this)(2, 0) * (-(*this)(0, 1) * (*this)(3, 2) + (*this)(3, 1) * (*this)(0, 2)) +
206  (*this)(3, 0) * ((*this)(0, 1) * (*this)(2, 2) - (*this)(2, 1) * (*this)(0, 2));
207  result[2] = (*this)(0, 0) * (-(*this)(1, 1) * (*this)(3, 2) + (*this)(3, 1) * (*this)(1, 2)) +
208  (*this)(1, 0) * ((*this)(0, 1) * (*this)(3, 2) - (*this)(3, 1) * (*this)(0, 2)) +
209  (*this)(3, 0) * (-(*this)(0, 1) * (*this)(1, 2) + (*this)(1, 1) * (*this)(0, 2));
210  result[3] = (*this)(0, 0) * ((*this)(1, 1) * (*this)(2, 2) - (*this)(2, 1) * (*this)(1, 2)) +
211  (*this)(1, 0) * (-(*this)(0, 1) * (*this)(2, 2) + (*this)(2, 1) * (*this)(0, 2)) +
212  (*this)(2, 0) * ((*this)(0, 1) * (*this)(1, 2) - (*this)(1, 1) * (*this)(0, 2));
213  break;
214  default:
215  logger(ERROR, "Wedge product not implemented for this dimension");
216  } //end switch
217 
218  return (result);
219 }
220 
221 template<unsigned int numberOfRows, unsigned int numberOfColumns>
223 {
224  int nr = numberOfRows;
225  int nc = numberOfColumns;
226  int nPivot = std::min(numberOfRows, numberOfColumns);
227  int iPivot[nPivot];
228 
229  SmallMatrix result(*this);
230 
231  int info;
232 
233  dgetrf_(&nr, &nc, result.data(), &nr, iPivot, &info);
234 
235  return result;
236 }
237 
238 //class template specialization for this one function is a waste of code duplication
239 //just let the compiler figure out which case it needs
240 template<unsigned int numberOfRows, unsigned int numberOfColumns>
242 {
243  logger.assert_debug(numberOfRows == numberOfColumns, "Matrix should be square to have a determinant!");
244 
245  switch (numberOfRows)
246  {
247  case 0:
248  return 1;
249  case 1:
250  return (*this)(0, 0);
251  case 2:
252  return (*this)(0, 0) * (*this)(1, 1) - (*this)(0, 1) * (*this)(1, 0);
253 
254  case 3:
255  return (*this)(0, 0) * ((*this)(1, 1) * (*this)(2, 2) - (*this)(1, 2) * (*this)(2, 1)) -
256  (*this)(0, 1) * ((*this)(1, 0) * (*this)(2, 2) - (*this)(2, 0) * (*this)(1, 2)) +
257  (*this)(0, 2) * ((*this)(1, 0) * (*this)(2, 1) - (*this)(2, 0) * (*this)(1, 1));
258 
259  case 4:
260  return ((*this)(3, 0) * (*this)(2, 1) * (*this)(0, 3) - (*this)(2, 0) * (*this)(3, 1) * (*this)(0, 3)) *
261  (*this)(1, 2) +
262  (-(*this)(3, 0) * (*this)(0, 3) * (*this)(2, 2) + (*this)(2, 0) * (*this)(0, 3) * (*this)(3, 2)) *
263  (*this)(1, 1) +
264  ((*this)(3, 1) * (*this)(0, 3) * (*this)(2, 2) - (*this)(2, 1) * (*this)(0, 3) * (*this)(3, 2)) *
265  (*this)(1, 0) +
266  (-(*this)(3, 0) * (*this)(2, 1) * (*this)(1, 3) + (*this)(2, 0) * (*this)(3, 1) * (*this)(1, 3) +
267  (-(*this)(2, 0) * (*this)(3, 3) + (*this)(3, 0) * (*this)(2, 3)) * (*this)(1, 1) +
268  ((*this)(2, 1) * (*this)(3, 3) - (*this)(3, 1) * (*this)(2, 3)) * (*this)(1, 0)) * (*this)(0, 2) +
269  ((*this)(3, 0) * (*this)(1, 3) * (*this)(2, 2) - (*this)(2, 0) * (*this)(1, 3) * (*this)(3, 2) +
270  ((*this)(2, 0) * (*this)(3, 3) - (*this)(3, 0) * (*this)(2, 3)) * (*this)(1, 2) +
271  (-(*this)(2, 2) * (*this)(3, 3) + (*this)(2, 3) * (*this)(3, 2)) * (*this)(1, 0)) * (*this)(0, 1) +
272  (-(*this)(3, 1) * (*this)(1, 3) * (*this)(2, 2) + (*this)(2, 1) * (*this)(1, 3) * (*this)(3, 2) +
273  ((*this)(3, 1) * (*this)(2, 3) - (*this)(2, 1) * (*this)(3, 3)) * (*this)(1, 2) +
274  (*this)(1, 1) * ((*this)(2, 2) * (*this)(3, 3) - (*this)(2, 3) * (*this)(3, 2))) * (*this)(0, 0);
275  // ... says Maple; this can possibly be done more efficiently,
276  // maybe even with LU (with pivoting, though...)
277  default:
278  logger(ERROR, "Computing the Determinant for size % is not implemented", numberOfRows);
279  break;
280  }
281  return 0;
282 }
283 
284 template<unsigned int numberOfRows, unsigned int numberOfColumns>
286 {
287  logger.assert_debug(numberOfRows == numberOfColumns, "Cannot invert a non-square matrix");
289 
290  int nr = numberOfRows;
291  int nc = numberOfColumns;
292 
293  int nPivot = numberOfRows;
294  int iPivot[nPivot];
295 
296  int info = 0;
297 
298  dgetrf_(&nr, &nc, result.data(), &nr, iPivot, &info);
299 
300  int lwork = numberOfRows * numberOfColumns;
301  SmallMatrix work;
302  dgetri_(&nc, result.data(), &nc, iPivot, work.data(), &lwork, &info);
303 
304  return result;
305 }
306 
307 template<unsigned int numberOfRows, unsigned int numberOfColumns>
308 template<unsigned int numberOfRightHandSideColumns>
310 {
311  logger.assert_debug(numberOfRows == numberOfColumns, "can only solve for square matrixes");
312 
313  int n = numberOfRows;
314  int nrhs = numberOfRightHandSideColumns;
315  int info;
316 
317  int IPIV[numberOfRows];
319  dgesv_(&n, &nrhs, matThis.data(), &n, IPIV, B.data(), &n, &info);
320 }
321 
322 template<unsigned int numberOfRows, unsigned int numberOfColumns>
324 {
325  logger.assert_debug(numberOfRows == numberOfColumns, "can only solve for square matrixes");
326 
327  int n = numberOfRows;
328  int nrhs = 1;
329  int info;
330 
331  int IPIV[numberOfRows];
332  SmallMatrix matThis = *this;
333  dgesv_(&n, &nrhs, matThis.data(), &n, IPIV, b.data(), &n, &info);
334 }
335 
336 template<unsigned int numberOfRows, unsigned int numberOfColumns>
338 {
339  if (numberOfColumns == 0)
340  {
341  logger(WARN, "Trying to multiply a vector with a matrix without any columns.");
343  }
344  int nr = numberOfRows;
345  int nc = numberOfColumns;
346 
347  int i_one = 1;
348  double d_one = 1.0;
349  double d_zero = 0.0;
350 
352 
353  logger(DEBUG, "Matrix size: % x % \n Vector size: %", nr, nc, vec.size());
354 
355  dgemv_("T", &nr, &nc, &d_one, mat.data(), &nr, vec.data(), &i_one, &d_zero, result.data(), &i_one);
356  return result;
357 }
358 
int i
Definition: BiCGSTAB_step_by_step.cpp:9
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
Logger< MERCURYDPM_LOGLEVEL > logger("MercuryKernel")
Definition of different loggers with certain modules. A user can define its own custom logger here.
LL< Log::ERROR > ERROR
Error log level.
Definition: Logger.cc:32
LL< Log::WARN > WARN
Warning log level.
Definition: Logger.cc:33
int data[]
Definition: Map_placement_new.cpp:1
int daxpy_(int *N, double *DA, double *DX, int *INCX, double *DY, int *INCY)
This is the gerneral scalar times vector + vector from blas, hence from blas level 1....
void dgemv_(const char *trans, int *m, int *n, double *alpha, double *A, int *LDA, double *x, int *incx, double *beta, double *y, int *incy)
This does matrix times vector and is from blas level 2.
void dgetrf_(int *M, int *N, double *A, int *lda, int *IPIV, int *INFO)
This is LU factorisation of the matrix A. This has been taken from LAPACK.
void dgesv_(int *N, int *NRHS, double *A, int *lda, int *IPIV, double *B, int *LDB, int *INFO)
This is used for solve Ax=B for x. Again this is from LAPACK.
int dgemm_(const char *transA, const char *transB, int *M, int *N, int *k, double *alpha, double *A, int *LDA, double *B, int *LDB, double *beta, double *C, int *LDC)
This is the gernal matrix multiplication from blas level 3.
void dgetri_(int *N, double *A, int *lda, int *IPIV, double *WORK, int *lwork, int *INFO)
This is the inverse calulation also from LAPACK. Calculates inverse if you pass it the LU factorisati...
Scalar * b
Definition: benchVecAdd.cpp:17
Matrix< SCALARA, Dynamic, Dynamic, opt_A > A
Definition: bench_gemm.cpp:47
Matrix< Scalar, Dynamic, Dynamic > C
Definition: bench_gemm.cpp:49
The matrix class, also used for vectors and row-vectors.
Definition: Eigen/Eigen/src/Core/Matrix.h:186
constexpr Storage & data()
Definition: SparseMatrix.h:205
operator*(const MatrixBase< Derived > &matrix, const UniformScaling< Scalar > &s)
Definition: Eigen/src/Geometry/Scaling.h:133
Data type for small dense matrix.
Definition: SmallMatrix.h:48
Mdouble determinant() const
Definition: SmallMatrix_impl.h:241
Mdouble * data()
Definition: SmallMatrix.h:331
void solve(SmallMatrix< numberOfRows, numberOfRightHandSideColumns > &B) const
solves Ax=B where A is the current matrix and B is passed in. The result is returned in B.
Definition: SmallMatrix_impl.h:309
SmallVector< numberOfRows > computeWedgeStuffVector() const
computeWedgeStuffVector.
Definition: SmallMatrix_impl.h:181
SmallMatrix & operator*=(const Mdouble &scalar)
Does matrix A_ij=scalar*A_ij.
Definition: SmallMatrix.h:200
SmallMatrix LUfactorisation() const
Return the LUfactorisation of the matrix.
Definition: SmallMatrix_impl.h:222
SmallVector< numberOfRows > operator*(SmallVector< numberOfColumns > &right)
Defines Matrix A times vector B and return vector C i.e. C_,j= A_ij B_,j.
Definition: SmallMatrix_impl.h:61
SmallMatrix inverse() const
return the inverse in the vector result. The size of result matches the matrix.
Definition: SmallMatrix_impl.h:285
Definition: SmallVector.h:42
unsigned int size() const
Definition: SmallVector.h:213
const Mdouble * data() const
Definition: SmallVector.h:218
Definition: matrices.h:74
@ N
Definition: constructor.cpp:22
#define min(a, b)
Definition: datatypes.h:22
Scalar * y
Definition: level1_cplx_impl.h:128
int RealScalar int RealScalar int * incy
Definition: level1_cplx_impl.h:124
RealScalar alpha
Definition: level1_cplx_impl.h:151
RealScalar RealScalar int * incx
Definition: level1_cplx_impl.h:27
const char const int const RealScalar const RealScalar const int * lda
Definition: level2_cplx_impl.h:20
int * m
Definition: level2_cplx_impl.h:294
int info
Definition: level2_cplx_impl.h:39
Scalar beta
Definition: level2_cplx_impl.h:36
char * trans
Definition: level2_impl.h:240
char char char int int * k
Definition: level2_impl.h:374
#define DEBUG
Definition: main.h:181
#define INFO(i)
Definition: mumps_solver.h:54
double K
Wave number.
Definition: sphere_scattering.cc:115
@ LDB
Definition: octree.h:49
list x
Definition: plotDoE.py:28
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2