dense_matrix_impl.h
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2018 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 #ifndef LIBMESH_DENSE_MATRIX_IMPL_H
20 #define LIBMESH_DENSE_MATRIX_IMPL_H
21 
22 // C++ Includes
23 #include <cstdlib> // *must* precede <cmath> for proper std:abs() on PGI, Sun Studio CC
24 #include <cmath> // for sqrt
25 
26 // Local Includes
27 #include "libmesh/dense_matrix.h"
28 #include "libmesh/dense_vector.h"
29 #include "libmesh/libmesh.h"
30 
31 namespace libMesh
32 {
33 
34 
35 
36 // ------------------------------------------------------------
37 // Dense Matrix member functions
38 
39 template<typename T>
41 {
42  if (this->use_blas_lapack)
43  this->_multiply_blas(M2, LEFT_MULTIPLY);
44  else
45  {
46  // (*this) <- M2 * (*this)
47  // Where:
48  // (*this) = (m x n),
49  // M2 = (m x p),
50  // M3 = (p x n)
51 
52  // M3 is a copy of *this before it gets resize()d
53  DenseMatrix<T> M3(*this);
54 
55  // Resize *this so that the result can fit
56  this->resize (M2.m(), M3.n());
57 
58  // Call the multiply function in the base class
59  this->multiply(*this, M2, M3);
60  }
61 }
62 
63 
64 
65 template<typename T>
66 template<typename T2>
68 {
69  // (*this) <- M2 * (*this)
70  // Where:
71  // (*this) = (m x n),
72  // M2 = (m x p),
73  // M3 = (p x n)
74 
75  // M3 is a copy of *this before it gets resize()d
76  DenseMatrix<T> M3(*this);
77 
78  // Resize *this so that the result can fit
79  this->resize (M2.m(), M3.n());
80 
81  // Call the multiply function in the base class
82  this->multiply(*this, M2, M3);
83 }
84 
85 
86 
87 template<typename T>
89 {
90  if (this->use_blas_lapack)
91  this->_multiply_blas(A, LEFT_MULTIPLY_TRANSPOSE);
92  else
93  {
94  //Check to see if we are doing (A^T)*A
95  if (this == &A)
96  {
97  //libmesh_here();
98  DenseMatrix<T> B(*this);
99 
100  // Simple but inefficient way
101  // return this->left_multiply_transpose(B);
102 
103  // More efficient, but more code way
104  // If A is mxn, the result will be a square matrix of Size n x n.
105  const unsigned int n_rows = A.m();
106  const unsigned int n_cols = A.n();
107 
108  // resize() *this and also zero out all entries.
109  this->resize(n_cols,n_cols);
110 
111  // Compute the lower-triangular part
112  for (unsigned int i=0; i<n_cols; ++i)
113  for (unsigned int j=0; j<=i; ++j)
114  for (unsigned int k=0; k<n_rows; ++k) // inner products are over n_rows
115  (*this)(i,j) += B(k,i)*B(k,j);
116 
117  // Copy lower-triangular part into upper-triangular part
118  for (unsigned int i=0; i<n_cols; ++i)
119  for (unsigned int j=i+1; j<n_cols; ++j)
120  (*this)(i,j) = (*this)(j,i);
121  }
122 
123  else
124  {
125  DenseMatrix<T> B(*this);
126 
127  this->resize (A.n(), B.n());
128 
129  libmesh_assert_equal_to (A.m(), B.m());
130  libmesh_assert_equal_to (this->m(), A.n());
131  libmesh_assert_equal_to (this->n(), B.n());
132 
133  const unsigned int m_s = A.n();
134  const unsigned int p_s = A.m();
135  const unsigned int n_s = this->n();
136 
137  // Do it this way because there is a
138  // decent chance (at least for constraint matrices)
139  // that A.transpose(i,k) = 0.
140  for (unsigned int i=0; i<m_s; i++)
141  for (unsigned int k=0; k<p_s; k++)
142  if (A.transpose(i,k) != 0.)
143  for (unsigned int j=0; j<n_s; j++)
144  (*this)(i,j) += A.transpose(i,k)*B(k,j);
145  }
146  }
147 
148 }
149 
150 
151 
152 template<typename T>
153 template<typename T2>
155 {
156  //Check to see if we are doing (A^T)*A
157  if (this == &A)
158  {
159  //libmesh_here();
160  DenseMatrix<T> B(*this);
161 
162  // Simple but inefficient way
163  // return this->left_multiply_transpose(B);
164 
165  // More efficient, but more code way
166  // If A is mxn, the result will be a square matrix of Size n x n.
167  const unsigned int n_rows = A.m();
168  const unsigned int n_cols = A.n();
169 
170  // resize() *this and also zero out all entries.
171  this->resize(n_cols,n_cols);
172 
173  // Compute the lower-triangular part
174  for (unsigned int i=0; i<n_cols; ++i)
175  for (unsigned int j=0; j<=i; ++j)
176  for (unsigned int k=0; k<n_rows; ++k) // inner products are over n_rows
177  (*this)(i,j) += B(k,i)*B(k,j);
178 
179  // Copy lower-triangular part into upper-triangular part
180  for (unsigned int i=0; i<n_cols; ++i)
181  for (unsigned int j=i+1; j<n_cols; ++j)
182  (*this)(i,j) = (*this)(j,i);
183  }
184 
185  else
186  {
187  DenseMatrix<T> B(*this);
188 
189  this->resize (A.n(), B.n());
190 
191  libmesh_assert_equal_to (A.m(), B.m());
192  libmesh_assert_equal_to (this->m(), A.n());
193  libmesh_assert_equal_to (this->n(), B.n());
194 
195  const unsigned int m_s = A.n();
196  const unsigned int p_s = A.m();
197  const unsigned int n_s = this->n();
198 
199  // Do it this way because there is a
200  // decent chance (at least for constraint matrices)
201  // that A.transpose(i,k) = 0.
202  for (unsigned int i=0; i<m_s; i++)
203  for (unsigned int k=0; k<p_s; k++)
204  if (A.transpose(i,k) != 0.)
205  for (unsigned int j=0; j<n_s; j++)
206  (*this)(i,j) += A.transpose(i,k)*B(k,j);
207  }
208 }
209 
210 
211 
212 template<typename T>
214 {
215  if (this->use_blas_lapack)
216  this->_multiply_blas(M3, RIGHT_MULTIPLY);
217  else
218  {
219  // (*this) <- M3 * (*this)
220  // Where:
221  // (*this) = (m x n),
222  // M2 = (m x p),
223  // M3 = (p x n)
224 
225  // M2 is a copy of *this before it gets resize()d
226  DenseMatrix<T> M2(*this);
227 
228  // Resize *this so that the result can fit
229  this->resize (M2.m(), M3.n());
230 
231  this->multiply(*this, M2, M3);
232  }
233 }
234 
235 
236 
237 template<typename T>
238 template<typename T2>
240 {
241  // (*this) <- M3 * (*this)
242  // Where:
243  // (*this) = (m x n),
244  // M2 = (m x p),
245  // M3 = (p x n)
246 
247  // M2 is a copy of *this before it gets resize()d
248  DenseMatrix<T> M2(*this);
249 
250  // Resize *this so that the result can fit
251  this->resize (M2.m(), M3.n());
252 
253  this->multiply(*this, M2, M3);
254 }
255 
256 
257 
258 
259 template<typename T>
261 {
262  if (this->use_blas_lapack)
263  this->_multiply_blas(B, RIGHT_MULTIPLY_TRANSPOSE);
264  else
265  {
266  //Check to see if we are doing B*(B^T)
267  if (this == &B)
268  {
269  //libmesh_here();
270  DenseMatrix<T> A(*this);
271 
272  // Simple but inefficient way
273  // return this->right_multiply_transpose(A);
274 
275  // More efficient, more code
276  // If B is mxn, the result will be a square matrix of Size m x m.
277  const unsigned int n_rows = B.m();
278  const unsigned int n_cols = B.n();
279 
280  // resize() *this and also zero out all entries.
281  this->resize(n_rows,n_rows);
282 
283  // Compute the lower-triangular part
284  for (unsigned int i=0; i<n_rows; ++i)
285  for (unsigned int j=0; j<=i; ++j)
286  for (unsigned int k=0; k<n_cols; ++k) // inner products are over n_cols
287  (*this)(i,j) += A(i,k)*A(j,k);
288 
289  // Copy lower-triangular part into upper-triangular part
290  for (unsigned int i=0; i<n_rows; ++i)
291  for (unsigned int j=i+1; j<n_rows; ++j)
292  (*this)(i,j) = (*this)(j,i);
293  }
294 
295  else
296  {
297  DenseMatrix<T> A(*this);
298 
299  this->resize (A.m(), B.m());
300 
301  libmesh_assert_equal_to (A.n(), B.n());
302  libmesh_assert_equal_to (this->m(), A.m());
303  libmesh_assert_equal_to (this->n(), B.m());
304 
305  const unsigned int m_s = A.m();
306  const unsigned int p_s = A.n();
307  const unsigned int n_s = this->n();
308 
309  // Do it this way because there is a
310  // decent chance (at least for constraint matrices)
311  // that B.transpose(k,j) = 0.
312  for (unsigned int j=0; j<n_s; j++)
313  for (unsigned int k=0; k<p_s; k++)
314  if (B.transpose(k,j) != 0.)
315  for (unsigned int i=0; i<m_s; i++)
316  (*this)(i,j) += A(i,k)*B.transpose(k,j);
317  }
318  }
319 }
320 
321 
322 
323 template<typename T>
324 template<typename T2>
326 {
327  //Check to see if we are doing B*(B^T)
328  if (this == &B)
329  {
330  //libmesh_here();
331  DenseMatrix<T> A(*this);
332 
333  // Simple but inefficient way
334  // return this->right_multiply_transpose(A);
335 
336  // More efficient, more code
337  // If B is mxn, the result will be a square matrix of Size m x m.
338  const unsigned int n_rows = B.m();
339  const unsigned int n_cols = B.n();
340 
341  // resize() *this and also zero out all entries.
342  this->resize(n_rows,n_rows);
343 
344  // Compute the lower-triangular part
345  for (unsigned int i=0; i<n_rows; ++i)
346  for (unsigned int j=0; j<=i; ++j)
347  for (unsigned int k=0; k<n_cols; ++k) // inner products are over n_cols
348  (*this)(i,j) += A(i,k)*A(j,k);
349 
350  // Copy lower-triangular part into upper-triangular part
351  for (unsigned int i=0; i<n_rows; ++i)
352  for (unsigned int j=i+1; j<n_rows; ++j)
353  (*this)(i,j) = (*this)(j,i);
354  }
355 
356  else
357  {
358  DenseMatrix<T> A(*this);
359 
360  this->resize (A.m(), B.m());
361 
362  libmesh_assert_equal_to (A.n(), B.n());
363  libmesh_assert_equal_to (this->m(), A.m());
364  libmesh_assert_equal_to (this->n(), B.m());
365 
366  const unsigned int m_s = A.m();
367  const unsigned int p_s = A.n();
368  const unsigned int n_s = this->n();
369 
370  // Do it this way because there is a
371  // decent chance (at least for constraint matrices)
372  // that B.transpose(k,j) = 0.
373  for (unsigned int j=0; j<n_s; j++)
374  for (unsigned int k=0; k<p_s; k++)
375  if (B.transpose(k,j) != 0.)
376  for (unsigned int i=0; i<m_s; i++)
377  (*this)(i,j) += A(i,k)*B.transpose(k,j);
378  }
379 }
380 
381 
382 
383 
384 template<typename T>
386  const DenseVector<T> & arg) const
387 {
388  // Make sure the input sizes are compatible
389  libmesh_assert_equal_to (this->n(), arg.size());
390 
391  // Resize and clear dest.
392  // Note: DenseVector::resize() also zeros the vector.
393  dest.resize(this->m());
394 
395  // Short-circuit if the matrix is empty
396  if(this->m() == 0 || this->n() == 0)
397  return;
398 
399  if (this->use_blas_lapack)
400  this->_matvec_blas(1., 0., dest, arg);
401  else
402  {
403  const unsigned int n_rows = this->m();
404  const unsigned int n_cols = this->n();
405 
406  for (unsigned int i=0; i<n_rows; i++)
407  for (unsigned int j=0; j<n_cols; j++)
408  dest(i) += (*this)(i,j)*arg(j);
409  }
410 }
411 
412 
413 
414 template<typename T>
415 template<typename T2>
417  const DenseVector<T2> & arg) const
418 {
419  // Make sure the input sizes are compatible
420  libmesh_assert_equal_to (this->n(), arg.size());
421 
422  // Resize and clear dest.
423  // Note: DenseVector::resize() also zeros the vector.
424  dest.resize(this->m());
425 
426  // Short-circuit if the matrix is empty
427  if (this->m() == 0 || this->n() == 0)
428  return;
429 
430  const unsigned int n_rows = this->m();
431  const unsigned int n_cols = this->n();
432 
433  for (unsigned int i=0; i<n_rows; i++)
434  for (unsigned int j=0; j<n_cols; j++)
435  dest(i) += (*this)(i,j)*arg(j);
436 }
437 
438 
439 
440 template<typename T>
442  const DenseVector<T> & arg) const
443 {
444  // Make sure the input sizes are compatible
445  libmesh_assert_equal_to (this->m(), arg.size());
446 
447  // Resize and clear dest.
448  // Note: DenseVector::resize() also zeros the vector.
449  dest.resize(this->n());
450 
451  // Short-circuit if the matrix is empty
452  if (this->m() == 0)
453  return;
454 
455  if (this->use_blas_lapack)
456  {
457  this->_matvec_blas(1., 0., dest, arg, /*trans=*/true);
458  }
459  else
460  {
461  const unsigned int n_rows = this->m();
462  const unsigned int n_cols = this->n();
463 
464  // WORKS
465  // for (unsigned int j=0; j<n_cols; j++)
466  // for (unsigned int i=0; i<n_rows; i++)
467  // dest(j) += (*this)(i,j)*arg(i);
468 
469  // ALSO WORKS, (i,j) just swapped
470  for (unsigned int i=0; i<n_cols; i++)
471  for (unsigned int j=0; j<n_rows; j++)
472  dest(i) += (*this)(j,i)*arg(j);
473  }
474 }
475 
476 
477 
478 template<typename T>
479 template<typename T2>
481  const DenseVector<T2> & arg) const
482 {
483  // Make sure the input sizes are compatible
484  libmesh_assert_equal_to (this->m(), arg.size());
485 
486  // Resize and clear dest.
487  // Note: DenseVector::resize() also zeros the vector.
488  dest.resize(this->n());
489 
490  // Short-circuit if the matrix is empty
491  if (this->m() == 0)
492  return;
493 
494  const unsigned int n_rows = this->m();
495  const unsigned int n_cols = this->n();
496 
497  // WORKS
498  // for (unsigned int j=0; j<n_cols; j++)
499  // for (unsigned int i=0; i<n_rows; i++)
500  // dest(j) += (*this)(i,j)*arg(i);
501 
502  // ALSO WORKS, (i,j) just swapped
503  for (unsigned int i=0; i<n_cols; i++)
504  for (unsigned int j=0; j<n_rows; j++)
505  dest(i) += (*this)(j,i)*arg(j);
506 }
507 
508 
509 
510 template<typename T>
512  const T factor,
513  const DenseVector<T> & arg) const
514 {
515  // Short-circuit if the matrix is empty
516  if (this->m() == 0)
517  {
518  dest.resize(0);
519  return;
520  }
521 
522  if (this->use_blas_lapack)
523  this->_matvec_blas(factor, 1., dest, arg);
524  else
525  {
526  DenseVector<T> temp(arg.size());
527  this->vector_mult(temp, arg);
528  dest.add(factor, temp);
529  }
530 }
531 
532 
533 
534 template<typename T>
535 template<typename T2, typename T3>
537  const T2 factor,
538  const DenseVector<T3> & arg) const
539 {
540  // Short-circuit if the matrix is empty
541  if (this->m() == 0)
542  {
543  dest.resize(0);
544  return;
545  }
546 
548  temp(arg.size());
549  this->vector_mult(temp, arg);
550  dest.add(factor, temp);
551 }
552 
553 
554 
555 template<typename T>
557  unsigned int sub_n,
558  DenseMatrix<T> & dest) const
559 {
560  libmesh_assert( (sub_m <= this->m()) && (sub_n <= this->n()) );
561 
562  dest.resize(sub_m, sub_n);
563  for (unsigned int i=0; i<sub_m; i++)
564  for (unsigned int j=0; j<sub_n; j++)
565  dest(i,j) = (*this)(i,j);
566 }
567 
568 
569 
570 template<typename T>
572  const DenseVector<T> & b)
573 {
574  const unsigned int m = a.size();
575  const unsigned int n = b.size();
576 
577  this->resize(m, n);
578  for (unsigned int i = 0; i < m; ++i)
579  for (unsigned int j = 0; j < n; ++j)
580  (*this)(i,j) = a(i) * libmesh_conj(b(j));
581 }
582 
583 
584 
585 template<typename T>
586 void DenseMatrix<T>::get_principal_submatrix (unsigned int sub_m, DenseMatrix<T> & dest) const
587 {
588  get_principal_submatrix(sub_m, sub_m, dest);
589 }
590 
591 
592 
593 template<typename T>
595 {
596  dest.resize(this->n(), this->m());
597 
598  for (unsigned int i=0; i<dest.m(); i++)
599  for (unsigned int j=0; j<dest.n(); j++)
600  dest(i,j) = (*this)(j,i);
601 }
602 
603 
604 
605 
606 template<typename T>
608  DenseVector<T> & x)
609 {
610  // Check to be sure that the matrix is square before attempting
611  // an LU-solve. In general, one can compute the LU factorization of
612  // a non-square matrix, but:
613  //
614  // Overdetermined systems (m>n) have a solution only if enough of
615  // the equations are linearly-dependent.
616  //
617  // Underdetermined systems (m<n) typically have infinitely many
618  // solutions.
619  //
620  // We don't want to deal with either of these ambiguous cases here...
621  libmesh_assert_equal_to (this->m(), this->n());
622 
623  switch(this->_decomposition_type)
624  {
625  case NONE:
626  {
627  if (this->use_blas_lapack)
628  this->_lu_decompose_lapack();
629  else
630  this->_lu_decompose ();
631  break;
632  }
633 
634  case LU_BLAS_LAPACK:
635  {
636  // Already factored, just need to call back_substitute.
637  if (this->use_blas_lapack)
638  break;
639  }
640  libmesh_fallthrough();
641 
642  case LU:
643  {
644  // Already factored, just need to call back_substitute.
645  if (!(this->use_blas_lapack))
646  break;
647  }
648  libmesh_fallthrough();
649 
650  default:
651  libmesh_error_msg("Error! This matrix already has a different decomposition...");
652  }
653 
654  if (this->use_blas_lapack)
655  this->_lu_back_substitute_lapack (b, x);
656  else
657  this->_lu_back_substitute (b, x);
658 }
659 
660 
661 
662 
663 
664 
665 template<typename T>
667  DenseVector<T> & x ) const
668 {
669  const unsigned int
670  n_cols = this->n();
671 
672  libmesh_assert_equal_to (this->m(), n_cols);
673  libmesh_assert_equal_to (this->m(), b.size());
674 
675  x.resize (n_cols);
676 
677  // A convenient reference to *this
678  const DenseMatrix<T> & A = *this;
679 
680  // Temporary vector storage. We use this instead of
681  // modifying the RHS.
682  DenseVector<T> z = b;
683 
684  // Lower-triangular "top to bottom" solve step, taking into account pivots
685  for (unsigned int i=0; i<n_cols; ++i)
686  {
687  // Swap
688  if (_pivots[i] != static_cast<pivot_index_t>(i))
689  std::swap( z(i), z(_pivots[i]) );
690 
691  x(i) = z(i);
692 
693  for (unsigned int j=0; j<i; ++j)
694  x(i) -= A(i,j)*x(j);
695 
696  x(i) /= A(i,i);
697  }
698 
699  // Upper-triangular "bottom to top" solve step
700  const unsigned int last_row = n_cols-1;
701 
702  for (int i=last_row; i>=0; --i)
703  {
704  for (int j=i+1; j<static_cast<int>(n_cols); ++j)
705  x(i) -= A(i,j)*x(j);
706  }
707 }
708 
709 
710 
711 
712 
713 
714 
715 
716 template<typename T>
718 {
719  // If this function was called, there better not be any
720  // previous decomposition of the matrix.
721  libmesh_assert_equal_to (this->_decomposition_type, NONE);
722 
723  // Get the matrix size and make sure it is square
724  const unsigned int
725  n_rows = this->m();
726 
727  // A convenient reference to *this
728  DenseMatrix<T> & A = *this;
729 
730  _pivots.resize(n_rows);
731 
732  for (unsigned int i=0; i<n_rows; ++i)
733  {
734  // Find the pivot row by searching down the i'th column
735  _pivots[i] = i;
736 
737  // std::abs(complex) must return a Real!
738  Real the_max = std::abs( A(i,i) );
739  for (unsigned int j=i+1; j<n_rows; ++j)
740  {
741  Real candidate_max = std::abs( A(j,i) );
742  if (the_max < candidate_max)
743  {
744  the_max = candidate_max;
745  _pivots[i] = j;
746  }
747  }
748 
749  // libMesh::out << "the_max=" << the_max << " found at row " << _pivots[i] << std::endl;
750 
751  // If the max was found in a different row, interchange rows.
752  // Here we interchange the *entire* row, in Gaussian elimination
753  // you would only interchange the subrows A(i,j) and A(p(i),j), for j>i
754  if (_pivots[i] != static_cast<pivot_index_t>(i))
755  {
756  for (unsigned int j=0; j<n_rows; ++j)
757  std::swap( A(i,j), A(_pivots[i], j) );
758  }
759 
760 
761  // If the max abs entry found is zero, the matrix is singular
762  if (A(i,i) == libMesh::zero)
763  libmesh_error_msg("Matrix A is singular!");
764 
765  // Scale upper triangle entries of row i by the diagonal entry
766  // Note: don't scale the diagonal entry itself!
767  const T diag_inv = 1. / A(i,i);
768  for (unsigned int j=i+1; j<n_rows; ++j)
769  A(i,j) *= diag_inv;
770 
771  // Update the remaining sub-matrix A[i+1:m][i+1:m]
772  // by subtracting off (the diagonal-scaled)
773  // upper-triangular part of row i, scaled by the
774  // i'th column entry of each row. In terms of
775  // row operations, this is:
776  // for each r > i
777  // SubRow(r) = SubRow(r) - A(r,i)*SubRow(i)
778  //
779  // If we were scaling the i'th column as well, like
780  // in Gaussian elimination, this would 'zero' the
781  // entry in the i'th column.
782  for (unsigned int row=i+1; row<n_rows; ++row)
783  for (unsigned int col=i+1; col<n_rows; ++col)
784  A(row,col) -= A(row,i) * A(i,col);
785 
786  } // end i loop
787 
788  // Set the flag for LU decomposition
789  this->_decomposition_type = LU;
790 }
791 
792 
793 
794 template<typename T>
796 {
797  // We use the LAPACK svd implementation
798  _svd_lapack(sigma);
799 }
800 
801 
802 template<typename T>
805  DenseMatrix<Number> & VT)
806 {
807  // We use the LAPACK svd implementation
808  _svd_lapack(sigma, U, VT);
809 }
810 
811 
812 
813 template<typename T>
815  DenseVector<T> & x,
816  Real rcond) const
817 {
818  _svd_solve_lapack(rhs, x, rcond);
819 }
820 
821 
822 
823 template<typename T>
825  DenseVector<T> & lambda_imag)
826 {
827  // We use the LAPACK eigenvalue problem implementation
828  _evd_lapack(lambda_real, lambda_imag);
829 }
830 
831 
832 
833 template<typename T>
835  DenseVector<T> & lambda_imag,
836  DenseMatrix<T> & VL)
837 {
838  // We use the LAPACK eigenvalue problem implementation
839  _evd_lapack(lambda_real, lambda_imag, &VL, nullptr);
840 }
841 
842 
843 
844 template<typename T>
846  DenseVector<T> & lambda_imag,
847  DenseMatrix<T> & VR)
848 {
849  // We use the LAPACK eigenvalue problem implementation
850  _evd_lapack(lambda_real, lambda_imag, nullptr, &VR);
851 }
852 
853 
854 
855 template<typename T>
857  DenseVector<T> & lambda_imag,
858  DenseMatrix<T> & VL,
859  DenseMatrix<T> & VR)
860 {
861  // We use the LAPACK eigenvalue problem implementation
862  _evd_lapack(lambda_real, lambda_imag, &VL, &VR);
863 }
864 
865 
866 
867 template<typename T>
869 {
870  switch(this->_decomposition_type)
871  {
872  case NONE:
873  {
874  // First LU decompose the matrix.
875  // Note that the lu_decompose routine will check to see if the
876  // matrix is square so we don't worry about it.
877  if (this->use_blas_lapack)
878  this->_lu_decompose_lapack();
879  else
880  this->_lu_decompose ();
881  }
882  case LU:
883  case LU_BLAS_LAPACK:
884  {
885  // Already decomposed, don't do anything
886  break;
887  }
888  default:
889  libmesh_error_msg("Error! Can't compute the determinant under the current decomposition.");
890  }
891 
892  // A variable to keep track of the running product of diagonal terms.
893  T determinant = 1.;
894 
895  // Loop over diagonal terms, computing the product. In practice,
896  // be careful because this value could easily become too large to
897  // fit in a double or float. To be safe, one should keep track of
898  // the power (of 10) of the determinant in a separate variable
899  // and maintain an order 1 value for the determinant itself.
900  unsigned int n_interchanges = 0;
901  for (unsigned int i=0; i<this->m(); i++)
902  {
903  if (this->_decomposition_type==LU)
904  if (_pivots[i] != static_cast<pivot_index_t>(i))
905  n_interchanges++;
906 
907  // Lapack pivots are 1-based!
908  if (this->_decomposition_type==LU_BLAS_LAPACK)
909  if (_pivots[i] != static_cast<pivot_index_t>(i+1))
910  n_interchanges++;
911 
912  determinant *= (*this)(i,i);
913  }
914 
915  // Compute sign of determinant, depends on number of row interchanges!
916  // The sign should be (-1)^{n}, where n is the number of interchanges.
917  Real sign = n_interchanges % 2 == 0 ? 1. : -1.;
918 
919  return sign*determinant;
920 }
921 
922 
923 
924 // The cholesky solve function first decomposes the matrix
925 // with cholesky_decompose and then uses the cholesky_back_substitute
926 // routine to find the solution x.
927 template <typename T>
928 template <typename T2>
930  DenseVector<T2> & x)
931 {
932  // Check for a previous decomposition
933  switch(this->_decomposition_type)
934  {
935  case NONE:
936  {
937  this->_cholesky_decompose ();
938  break;
939  }
940 
941  case CHOLESKY:
942  {
943  // Already factored, just need to call back_substitute.
944  break;
945  }
946 
947  default:
948  libmesh_error_msg("Error! This matrix already has a different decomposition...");
949  }
950 
951  // Perform back substitution
952  this->_cholesky_back_substitute (b, x);
953 }
954 
955 
956 
957 
958 // This algorithm is based on the Cholesky decomposition in
959 // the Numerical Recipes in C book.
960 template<typename T>
962 {
963  // If we called this function, there better not be any
964  // previous decomposition of the matrix.
965  libmesh_assert_equal_to (this->_decomposition_type, NONE);
966 
967  // Shorthand notation for number of rows and columns.
968  const unsigned int
969  n_rows = this->m(),
970  n_cols = this->n();
971 
972  // Just to be really sure...
973  libmesh_assert_equal_to (n_rows, n_cols);
974 
975  // A convenient reference to *this
976  DenseMatrix<T> & A = *this;
977 
978  for (unsigned int i=0; i<n_rows; ++i)
979  {
980  for (unsigned int j=i; j<n_cols; ++j)
981  {
982  for (unsigned int k=0; k<i; ++k)
983  A(i,j) -= A(i,k) * A(j,k);
984 
985  if (i == j)
986  {
987 #ifndef LIBMESH_USE_COMPLEX_NUMBERS
988  if (A(i,j) <= 0.0)
989  libmesh_error_msg("Error! Can only use Cholesky decomposition with symmetric positive definite matrices.");
990 #endif
991 
992  A(i,i) = std::sqrt(A(i,j));
993  }
994  else
995  A(j,i) = A(i,j) / A(i,i);
996  }
997  }
998 
999  // Set the flag for CHOLESKY decomposition
1000  this->_decomposition_type = CHOLESKY;
1001 }
1002 
1003 
1004 
1005 template <typename T>
1006 template <typename T2>
1008  DenseVector<T2> & x) const
1009 {
1010  // Shorthand notation for number of rows and columns.
1011  const unsigned int
1012  n_rows = this->m(),
1013  n_cols = this->n();
1014 
1015  // Just to be really sure...
1016  libmesh_assert_equal_to (n_rows, n_cols);
1017 
1018  // A convenient reference to *this
1019  const DenseMatrix<T> & A = *this;
1020 
1021  // Now compute the solution to Ax =b using the factorization.
1022  x.resize(n_rows);
1023 
1024  // Solve for Ly=b
1025  for (unsigned int i=0; i<n_cols; ++i)
1026  {
1027  T2 temp = b(i);
1028 
1029  for (unsigned int k=0; k<i; ++k)
1030  temp -= A(i,k)*x(k);
1031 
1032  x(i) = temp / A(i,i);
1033  }
1034 
1035  // Solve for L^T x = y
1036  for (unsigned int i=0; i<n_cols; ++i)
1037  {
1038  const unsigned int ib = (n_cols-1)-i;
1039 
1040  for (unsigned int k=(ib+1); k<n_cols; ++k)
1041  x(ib) -= A(k,ib) * x(k);
1042 
1043  x(ib) /= A(ib,ib);
1044  }
1045 }
1046 
1047 
1048 
1049 
1050 
1051 
1052 
1053 
1054 // This routine is commented out since it is not really a memory
1055 // efficient implementation. Also, you don't *need* the inverse
1056 // for anything, instead just use lu_solve to solve Ax=b.
1057 // template<typename T>
1058 // void DenseMatrix<T>::inverse ()
1059 // {
1060 // // First LU decompose the matrix
1061 // // Note that the lu_decompose routine will check to see if the
1062 // // matrix is square so we don't worry about it.
1063 // if (!this->_lu_decomposed)
1064 // this->_lu_decompose();
1065 
1066 // // A unit vector which will be used as a rhs
1067 // // to pick off a single value each time.
1068 // DenseVector<T> e;
1069 // e.resize(this->m());
1070 
1071 // // An empty vector which will be used to hold the solution
1072 // // to the back substitutions.
1073 // DenseVector<T> x;
1074 // x.resize(this->m());
1075 
1076 // // An empty dense matrix to store the resulting inverse
1077 // // temporarily until we can overwrite A.
1078 // DenseMatrix<T> inv;
1079 // inv.resize(this->m(), this->n());
1080 
1081 // // Resize the passed in matrix to hold the inverse
1082 // inv.resize(this->m(), this->n());
1083 
1084 // for (unsigned int j=0; j<this->n(); ++j)
1085 // {
1086 // e.zero();
1087 // e(j) = 1.;
1088 // this->_lu_back_substitute(e, x, false);
1089 // for (unsigned int i=0; i<this->n(); ++i)
1090 // inv(i,j) = x(i);
1091 // }
1092 
1093 // // Now overwrite all the entries
1094 // *this = inv;
1095 // }
1096 
1097 
1098 } // namespace libMesh
1099 
1100 #endif // LIBMESH_DENSE_MATRIX_IMPL_H
virtual unsigned int size() const override
Definition: dense_vector.h:92
double abs(double a)
void _lu_back_substitute(const DenseVector< T > &b, DenseVector< T > &x) const
T transpose(const unsigned int i, const unsigned int j) const
boostcopy::enable_if_c< ScalarTraits< T2 >::value, void >::type add(const T2 factor, const DenseVector< T3 > &vec)
Definition: dense_vector.h:436
T libmesh_conj(T a)
void resize(const unsigned int n)
Definition: dense_vector.h:355
void evd_left(DenseVector< T > &lambda_real, DenseVector< T > &lambda_imag, DenseMatrix< T > &VL)
void evd_left_and_right(DenseVector< T > &lambda_real, DenseVector< T > &lambda_imag, DenseMatrix< T > &VL, DenseMatrix< T > &VR)
unsigned int m() const
const Number zero
Definition: libmesh.h:239
void outer_product(const DenseVector< T > &a, const DenseVector< T > &b)
void vector_mult_transpose(DenseVector< T > &dest, const DenseVector< T > &arg) const
void vector_mult_add(DenseVector< T > &dest, const T factor, const DenseVector< T > &arg) const
void svd(DenseVector< Real > &sigma)
void get_transpose(DenseMatrix< T > &dest) const
virtual void right_multiply(const DenseMatrixBase< T > &M2) override
void evd_right(DenseVector< T > &lambda_real, DenseVector< T > &lambda_imag, DenseMatrix< T > &VR)
void _cholesky_back_substitute(const DenseVector< T2 > &b, DenseVector< T2 > &x) const
void get_principal_submatrix(unsigned int sub_m, unsigned int sub_n, DenseMatrix< T > &dest) const
void right_multiply_transpose(const DenseMatrix< T > &A)
void lu_solve(const DenseVector< T > &b, DenseVector< T > &x)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void evd(DenseVector< T > &lambda_real, DenseVector< T > &lambda_imag)
void left_multiply_transpose(const DenseMatrix< T > &A)
static PetscErrorCode Mat * A
void swap(Iterator &lhs, Iterator &rhs)
void resize(const unsigned int new_m, const unsigned int new_n)
Definition: dense_matrix.h:792
virtual void left_multiply(const DenseMatrixBase< T > &M2) override
void cholesky_solve(const DenseVector< T2 > &b, DenseVector< T2 > &x)
unsigned int n() const
A matrix object used for finite element assembly and numerics.
Definition: dense_matrix.h:54
void svd_solve(const DenseVector< T > &rhs, DenseVector< T > &x, Real rcond=std::numeric_limits< Real >::epsilon()) const
void vector_mult(DenseVector< T > &dest, const DenseVector< T > &arg) const