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