deal.II version 9.7.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
lapack_full_matrix.cc
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2005 - 2025 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
16
25#include <deal.II/lac/vector.h>
26
27#include <iomanip>
28#include <iostream>
29
31
32using namespace LAPACKSupport;
33
34namespace internal
35{
37 {
38 // ZGEEV/CGEEV and DGEEV/SGEEV need different work arrays and different
39 // output arrays for eigenvalues. This makes working with generic scalar
40 // types a bit difficult. To get around this, geev_helper has the same
41 // signature for real and complex arguments, but it ignores some
42 // parameters when called with a real type and ignores different
43 // parameters when called with a complex type.
44 template <typename T>
45 void
46 geev_helper(const char vl,
47 const char vr,
49 const types::blas_int n_rows,
50 std::vector<T> &real_part_eigenvalues,
51 std::vector<T> &imag_part_eigenvalues,
52 std::vector<T> &left_eigenvectors,
53 std::vector<T> &right_eigenvectors,
54 std::vector<T> &real_work,
55 std::vector<T> & /*complex_work*/,
56 const types::blas_int work_flag,
57 types::blas_int &info)
58 {
59 static_assert(std::is_same_v<T, double> || std::is_same_v<T, float>,
60 "Only implemented for double and float");
61 Assert(matrix.size() == static_cast<std::size_t>(n_rows * n_rows),
63 Assert(static_cast<std::size_t>(n_rows) <= real_part_eigenvalues.size(),
65 Assert(static_cast<std::size_t>(n_rows) <= imag_part_eigenvalues.size(),
67 if (vl == 'V')
68 Assert(static_cast<std::size_t>(n_rows * n_rows) <=
69 left_eigenvectors.size(),
71 if (vr == 'V')
72 Assert(static_cast<std::size_t>(n_rows * n_rows) <=
73 right_eigenvectors.size(),
75 Assert(work_flag == -1 ||
76 static_cast<std::size_t>(2 * n_rows) <= real_work.size(),
78 Assert(work_flag == -1 || std::max<long int>(1, 3 * n_rows) <= work_flag,
80 geev(&vl,
81 &vr,
82 &n_rows,
83 matrix.data(),
84 &n_rows,
85 real_part_eigenvalues.data(),
86 imag_part_eigenvalues.data(),
87 left_eigenvectors.data(),
88 &n_rows,
89 right_eigenvectors.data(),
90 &n_rows,
91 real_work.data(),
92 &work_flag,
93 &info);
94 }
95
96
97
98 template <typename T>
99 void
100 geev_helper(const char vl,
101 const char vr,
102 AlignedVector<std::complex<T>> &matrix,
103 const types::blas_int n_rows,
104 std::vector<T> & /*real_part_eigenvalues*/,
105 std::vector<std::complex<T>> &eigenvalues,
106 std::vector<std::complex<T>> &left_eigenvectors,
107 std::vector<std::complex<T>> &right_eigenvectors,
108 std::vector<std::complex<T>> &complex_work,
109 std::vector<T> &real_work,
110 const types::blas_int work_flag,
111 types::blas_int &info)
112 {
113 static_assert(
114 std::is_same_v<T, double> || std::is_same_v<T, float>,
115 "Only implemented for std::complex<double> and std::complex<float>");
116 Assert(matrix.size() == static_cast<std::size_t>(n_rows * n_rows),
118 Assert(static_cast<std::size_t>(n_rows) <= eigenvalues.size(),
120 if (vl == 'V')
121 Assert(static_cast<std::size_t>(n_rows * n_rows) <=
122 left_eigenvectors.size(),
124 if (vr == 'V')
125 Assert(static_cast<std::size_t>(n_rows * n_rows) <=
126 right_eigenvectors.size(),
128 Assert(work_flag == -1 ||
129 std::max<std::size_t>(1, work_flag) <= real_work.size(),
131 Assert(work_flag == -1 || std::max<long int>(1, 2 * n_rows) <= work_flag,
133
134 geev(&vl,
135 &vr,
136 &n_rows,
137 matrix.data(),
139 eigenvalues.data(),
140 left_eigenvectors.data(),
141 &n_rows,
142 right_eigenvectors.data(),
143 &n_rows,
144 complex_work.data(),
145 &work_flag,
146 real_work.data(),
147 &info);
148 }
149
150
151
152 template <typename T>
153 void
154 gesdd_helper(const char job,
155 const types::blas_int n_rows,
156 const types::blas_int n_cols,
158 std::vector<T> &singular_values,
159 AlignedVector<T> &left_vectors,
160 AlignedVector<T> &right_vectors,
161 std::vector<T> &real_work,
162 std::vector<T> & /*complex work*/,
163 std::vector<types::blas_int> &integer_work,
164 const types::blas_int work_flag,
165 types::blas_int &info)
166 {
167 Assert(job == 'A' || job == 'S' || job == 'O' || job == 'N',
169 Assert(static_cast<std::size_t>(n_rows * n_cols) == matrix.size(),
171 Assert(std::min<std::size_t>(n_rows, n_cols) <= singular_values.size(),
173 Assert(8 * std::min<std::size_t>(n_rows, n_cols) <= integer_work.size(),
175 Assert(work_flag == -1 ||
176 static_cast<std::size_t>(work_flag) <= real_work.size(),
178 gesdd(&job,
179 &n_rows,
180 &n_cols,
181 matrix.data(),
182 &n_rows,
183 singular_values.data(),
184 left_vectors.data(),
185 &n_rows,
186 right_vectors.data(),
187 &n_cols,
188 real_work.data(),
189 &work_flag,
190 integer_work.data(),
191 &info);
192 }
193
194
195
196 template <typename T>
197 void
198 gesdd_helper(const char job,
199 const types::blas_int n_rows,
200 const types::blas_int n_cols,
201 AlignedVector<std::complex<T>> &matrix,
202 std::vector<T> &singular_values,
203 AlignedVector<std::complex<T>> &left_vectors,
204 AlignedVector<std::complex<T>> &right_vectors,
205 std::vector<std::complex<T>> &work,
206 std::vector<T> &real_work,
207 std::vector<types::blas_int> &integer_work,
208 const types::blas_int &work_flag,
209 types::blas_int &info)
210 {
211 Assert(job == 'A' || job == 'S' || job == 'O' || job == 'N',
213 Assert(static_cast<std::size_t>(n_rows * n_cols) == matrix.size(),
215 Assert(static_cast<std::size_t>(std::min(n_rows, n_cols)) <=
216 singular_values.size(),
218 Assert(8 * std::min<std::size_t>(n_rows, n_cols) <= integer_work.size(),
220 Assert(work_flag == -1 ||
221 static_cast<std::size_t>(work_flag) <= real_work.size(),
223
224 gesdd(&job,
225 &n_rows,
226 &n_cols,
227 matrix.data(),
228 &n_rows,
229 singular_values.data(),
230 left_vectors.data(),
231 &n_rows,
232 right_vectors.data(),
233 &n_cols,
234 work.data(),
235 &work_flag,
236 real_work.data(),
237 integer_work.data(),
238 &info);
239 }
240 } // namespace LAPACKFullMatrixImplementation
241} // namespace internal
242
243
244
245template <typename number>
251
252
253
254template <typename number>
256 : TransposeTable<number>(m, n)
257 , state(matrix)
259{}
260
261
262
263template <typename number>
269
270
271
272template <typename number>
275{
277 state = M.state;
278 property = M.property;
279 return *this;
280}
281
282
283
284template <typename number>
285void
291
292
293
294template <typename number>
295void
297{
298 const size_type s = std::min(std::min(this->m(), n), this->n());
299 TransposeTable<number> copy(std::move(*this));
301 for (size_type i = 0; i < s; ++i)
302 for (size_type j = 0; j < s; ++j)
303 (*this)(i, j) = copy(i, j);
304}
305
306
307
308template <typename number>
309void
311 const std::array<number, 3> &csr,
312 const size_type i,
313 const size_type k,
314 const bool left)
315{
316 auto &A = *this;
317 // see Golub 2013 "Matrix computations", p241 5.1.9 Applying Givens
318 // Rotations but note the difference in notation, namely the sign of s: we
319 // have G * A, where G[1,1] = s
320 if (left)
321 {
322 for (size_type j = 0; j < A.n(); ++j)
323 {
324 const number t = A(i, j);
325 A(i, j) = csr[0] * A(i, j) + csr[1] * A(k, j);
326 A(k, j) = -csr[1] * t + csr[0] * A(k, j);
327 }
328 }
329 else
330 {
331 for (size_type j = 0; j < A.m(); ++j)
332 {
333 const number t = A(j, i);
334 A(j, i) = csr[0] * A(j, i) + csr[1] * A(j, k);
335 A(j, k) = -csr[1] * t + csr[0] * A(j, k);
336 }
337 }
338}
339
340
341
342template <typename number>
343void
345 const size_type col)
346{
347 AssertIndexRange(row, this->m());
348 AssertIndexRange(col, this->n());
349
350 const size_type nrows = this->m() - 1;
351 const size_type ncols = this->n() - 1;
352
353 TransposeTable<number> copy(std::move(*this));
354 this->TransposeTable<number>::reinit(nrows, ncols);
355
356 for (size_type j = 0; j < ncols; ++j)
357 {
358 const size_type jj = (j < col ? j : j + 1);
359 for (size_type i = 0; i < nrows; ++i)
360 {
361 const size_type ii = (i < row ? i : i + 1);
362 (*this)(i, j) = copy(ii, jj);
363 }
364 }
365}
366
367
368
369template <typename number>
370void
376
377
378
379template <typename number>
380template <typename number2>
383{
384 Assert(this->m() == M.m(), ExcDimensionMismatch(this->m(), M.m()));
385 Assert(this->n() == M.n(), ExcDimensionMismatch(this->n(), M.n()));
386 for (size_type i = 0; i < this->m(); ++i)
387 for (size_type j = 0; j < this->n(); ++j)
388 (*this)(i, j) = M(i, j);
389
391 property = LAPACKSupport::general;
392 return *this;
393}
394
395
396
397template <typename number>
398template <typename number2>
401{
402 Assert(this->m() == M.n(), ExcDimensionMismatch(this->m(), M.n()));
403 Assert(this->n() == M.m(), ExcDimensionMismatch(this->n(), M.m()));
404 for (size_type i = 0; i < this->m(); ++i)
405 for (size_type j = 0; j < this->n(); ++j)
406 (*this)(i, j) = M.el(i, j);
407
409 property = LAPACKSupport::general;
410 return *this;
411}
412
413
414
415template <typename number>
418{
420
421 if (this->n_elements() != 0)
422 this->reset_values();
423
425 return *this;
426}
427
429
430template <typename number>
433{
436 ExcState(state));
437
438 AssertIsFinite(factor);
439 const char type = 'G';
440 const number cfrom = 1.;
441 const types::blas_int m = this->m();
442 const types::blas_int n = this->n();
443 const types::blas_int lda = this->m();
444 types::blas_int info = 0;
445 // kl and ku will not be referenced for type = G (dense matrices).
446 const types::blas_int kl = 0;
447 number *values = this->values.data();
448
449 lascl(&type, &kl, &kl, &cfrom, &factor, &m, &n, values, &lda, &info);
450
451 // Negative return value implies a wrong argument. This should be internal.
452 Assert(info >= 0, ExcInternalError());
453
454 return *this;
455}
457
458
459template <typename number>
462{
465 ExcState(state));
466
467 AssertIsFinite(factor);
468 Assert(factor != number(0.), ExcZero());
469
470 const char type = 'G';
471 const number cto = 1.;
472 const types::blas_int m = this->m();
473 const types::blas_int n = this->n();
474 const types::blas_int lda = this->m();
475 types::blas_int info = 0;
476 // kl and ku will not be referenced for type = G (dense matrices).
477 const types::blas_int kl = 0;
478 number *values = this->values.data();
479
480 lascl(&type, &kl, &kl, &factor, &cto, &m, &n, values, &lda, &info);
481
482 // Negative return value implies a wrong argument. This should be internal.
483 Assert(info >= 0, ExcInternalError());
484
485 return *this;
486}
487
488
489
490template <typename number>
491void
493{
496 ExcState(state));
497
498 Assert(m() == A.m(), ExcDimensionMismatch(m(), A.m()));
499 Assert(n() == A.n(), ExcDimensionMismatch(n(), A.n()));
500
502
503 // BLAS does not offer functions to add matrices.
504 // LapackFullMatrix is stored in contiguous array
505 // ==> use BLAS 1 for adding vectors
506 const types::blas_int n = this->m() * this->n();
507 const types::blas_int inc = 1;
508 number *values = this->values.data();
509 const number *values_A = A.values.data();
510
511 axpy(&n, &a, values_A, &inc, values, &inc);
512}
513
514
515
516namespace
517{
518 template <typename number>
519 void
520 cholesky_rank1(LAPACKFullMatrix<number> &A,
521 const number a,
522 const Vector<number> &v)
523 {
524 const typename LAPACKFullMatrix<number>::size_type N = A.n();
525 Vector<number> z(v);
526 // Cholesky update / downdate, see
527 // 6.5.4 Cholesky Updating and Downdating, Golub 2013 Matrix computations
528 // Note that potrf() is called with LAPACKSupport::L , so the
529 // factorization is stored in lower triangular part.
530 // Also see discussion here
531 // http://icl.cs.utk.edu/lapack-forum/viewtopic.php?f=2&t=2646
532 if (a > 0.)
533 {
534 // simple update via a sequence of Givens rotations.
535 // Observe that
536 //
537 // | L^T |T | L^T |
538 // A <-- | | | | = L L^T + z z^T
539 // | z^T | | z^T |
540 //
541 // so we can get updated factor by doing a sequence of Givens
542 // rotations to make the matrix lower-triangular
543 // Also see LINPACK's dchud http://www.netlib.org/linpack/dchud.f
544 z *= std::sqrt(a);
545 for (typename LAPACKFullMatrix<number>::size_type k = 0; k < N; ++k)
546 {
547 const std::array<number, 3> csr =
549 A(k, k) = csr[2];
550 for (typename LAPACKFullMatrix<number>::size_type i = k + 1; i < N;
551 ++i)
552 {
553 const number t = A(i, k);
554 A(i, k) = csr[0] * A(i, k) + csr[1] * z(i);
555 z(i) = -csr[1] * t + csr[0] * z(i);
556 }
557 }
559 else
560 {
561 // downdating is not always possible as we may end up with
562 // negative definite matrix. If it's possible, then it boils
563 // down to application of hyperbolic rotations.
564 // Observe that
565 //
566 // | L^T |T | L^T |
567 // A <-- | | S | | = L L^T - z z^T
568 // | z^T | | z^T |
569 //
570 // |In 0 |
571 // S := | |
572 // |0 -1 |
573 //
574 // We are looking for H which is S-orthogonal (HSH^T=S) and
575 // can restore upper-triangular factor of the factorization of A above.
576 // We will use Hyperbolic rotations to do the job
577 //
578 // | c -s | | x1 | | r |
579 // | | = | | = | |, c^2 - s^2 = 1
580 // |-s c | | x2 | | 0 |
581 //
582 // which have real solution only if x2 <= x1.
583 // See also Linpack's http://www.netlib.org/linpack/dchdd.f and
584 // https://infoscience.epfl.ch/record/161468/files/cholupdate.pdf and
585 // "Analysis of a recursive Least Squares Hyperbolic Rotation Algorithm
586 // for Signal Processing", Alexander, Pan, Plemmons, 1988.
587 z *= std::sqrt(-a);
588 for (typename LAPACKFullMatrix<number>::size_type k = 0; k < N; ++k)
589 {
590 const std::array<number, 3> csr =
592 A(k, k) = csr[2];
593 for (typename LAPACKFullMatrix<number>::size_type i = k + 1; i < N;
594 ++i)
595 {
596 const number t = A(i, k);
597 A(i, k) = csr[0] * A(i, k) - csr[1] * z(i);
598 z(i) = -csr[1] * t + csr[0] * z(i);
599 }
600 }
601 }
602 }
604
605 template <typename number>
606 void
607 cholesky_rank1(LAPACKFullMatrix<std::complex<number>> & /*A*/,
608 const std::complex<number> /*a*/,
609 const Vector<std::complex<number>> & /*v*/)
610 {
612 }
613} // namespace
614
615
616
617template <typename number>
618void
620{
622
623 Assert(n() == m(), ExcInternalError());
624 Assert(m() == v.size(), ExcDimensionMismatch(m(), v.size()));
625
629 {
630 {
631 const types::blas_int N = this->m();
632 const char uplo = LAPACKSupport::U;
633 const types::blas_int lda = N;
634 const types::blas_int incx = 1;
635
636 syr(&uplo, &N, &a, v.begin(), &incx, this->values.begin(), &lda);
637 }
638
639 const size_type N = this->m();
640 // FIXME: we should really only update upper or lower triangular parts
641 // of symmetric matrices and make sure the interface is consistent,
642 // for example operator(i,j) gives correct results regardless of storage.
643 for (size_type i = 0; i < N; ++i)
644 for (size_type j = 0; j < i; ++j)
645 (*this)(i, j) = (*this)(j, i);
647 else if (state == LAPACKSupport::cholesky)
648 {
649 cholesky_rank1(*this, a, v);
650 }
651 else
652 AssertThrow(false, ExcState(state));
653}
655
656
657template <typename number>
658void
660 const Vector<number> &v,
661 const bool adding) const
662{
663 const types::blas_int mm = this->m();
664 const types::blas_int nn = this->n();
665 const number alpha = 1.;
666 const number beta = (adding ? 1. : 0.);
667 const number null = 0.;
668
669 // use trmv for triangular matrices
671 (mm == nn) && state == matrix)
672 {
673 Assert(adding == false, ExcNotImplemented());
674
675 AssertDimension(v.size(), this->n());
676 AssertDimension(w.size(), this->m());
677
678 const char diag = 'N';
679 const char trans = 'N';
680 const char uplo =
682
683 w = v;
684
685 const types::blas_int N = mm;
686 const types::blas_int lda = N;
687 const types::blas_int incx = 1;
688
689 trmv(
690 &uplo, &trans, &diag, &N, this->values.data(), &lda, w.data(), &incx);
691
692 return;
694
695 switch (state)
696 {
697 case matrix:
698 case inverse_matrix:
699 {
700 AssertDimension(v.size(), this->n());
701 AssertDimension(w.size(), this->m());
702
703 gemv("N",
704 &mm,
705 &nn,
706 &alpha,
707 this->values.data(),
708 &mm,
709 v.data(),
710 &one,
711 &beta,
712 w.data(),
713 &one);
714 break;
715 }
716 case svd:
717 {
718 std::lock_guard<std::mutex> lock(mutex);
719 AssertDimension(v.size(), this->n());
720 AssertDimension(w.size(), this->m());
721 // Compute V^T v
722 work.resize(std::max(mm, nn));
723 gemv("N",
724 &nn,
725 &nn,
726 &alpha,
727 svd_vt->values.data(),
728 &nn,
729 v.data(),
730 &one,
731 &null,
732 work.data(),
733 &one);
734 // Multiply by singular values
735 for (size_type i = 0; i < wr.size(); ++i)
736 work[i] *= wr[i];
737 // Multiply with U
738 gemv("N",
739 &mm,
740 &mm,
741 &alpha,
742 svd_u->values.data(),
743 &mm,
744 work.data(),
745 &one,
746 &beta,
747 w.data(),
748 &one);
749 break;
750 }
751 case inverse_svd:
752 {
753 std::lock_guard<std::mutex> lock(mutex);
754 AssertDimension(w.size(), this->n());
755 AssertDimension(v.size(), this->m());
756 // Compute U^T v
757 work.resize(std::max(mm, nn));
758 gemv("T",
759 &mm,
760 &mm,
761 &alpha,
762 svd_u->values.data(),
763 &mm,
764 v.data(),
765 &one,
766 &null,
767 work.data(),
768 &one);
769 // Multiply by singular values
770 for (size_type i = 0; i < wr.size(); ++i)
771 work[i] *= wr[i];
772 // Multiply with V
773 gemv("T",
774 &nn,
775 &nn,
776 &alpha,
777 svd_vt->values.data(),
778 &nn,
779 work.data(),
780 &one,
781 &beta,
782 w.data(),
783 &one);
784 break;
785 }
786 default:
787 Assert(false, ExcState(state));
788 }
789}
790
791
792
793template <typename number>
794void
796 const Vector<number> &v,
797 const bool adding) const
798{
799 const types::blas_int mm = this->m();
800 const types::blas_int nn = this->n();
801 const number alpha = 1.;
802 const number beta = (adding ? 1. : 0.);
803 const number null = 0.;
804
805 // use trmv for triangular matrices
807 (mm == nn) && state == matrix)
808 {
809 Assert(adding == false, ExcNotImplemented());
810
811 AssertDimension(v.size(), this->n());
812 AssertDimension(w.size(), this->m());
813
814 const char diag = 'N';
815 const char trans = 'T';
816 const char uplo =
818
819 w = v;
820
821 const types::blas_int N = mm;
822 const types::blas_int lda = N;
823 const types::blas_int incx = 1;
824
825 trmv(
826 &uplo, &trans, &diag, &N, this->values.data(), &lda, w.data(), &incx);
827
828 return;
829 }
830
831
832 switch (state)
833 {
834 case matrix:
835 case inverse_matrix:
836 {
837 AssertDimension(w.size(), this->n());
838 AssertDimension(v.size(), this->m());
839
840 gemv("T",
841 &mm,
842 &nn,
843 &alpha,
844 this->values.data(),
845 &mm,
846 v.data(),
847 &one,
848 &beta,
849 w.data(),
850 &one);
851 break;
852 }
853 case svd:
854 {
855 std::lock_guard<std::mutex> lock(mutex);
856 AssertDimension(w.size(), this->n());
857 AssertDimension(v.size(), this->m());
858
859 // Compute U^T v
860 work.resize(std::max(mm, nn));
861 gemv("T",
862 &mm,
863 &mm,
864 &alpha,
865 svd_u->values.data(),
866 &mm,
867 v.data(),
868 &one,
869 &null,
870 work.data(),
871 &one);
872 // Multiply by singular values
873 for (size_type i = 0; i < wr.size(); ++i)
874 work[i] *= wr[i];
875 // Multiply with V
876 gemv("T",
877 &nn,
878 &nn,
879 &alpha,
880 svd_vt->values.data(),
881 &nn,
882 work.data(),
883 &one,
884 &beta,
885 w.data(),
886 &one);
887 break;
888 }
889 case inverse_svd:
890 {
891 std::lock_guard<std::mutex> lock(mutex);
892 AssertDimension(v.size(), this->n());
893 AssertDimension(w.size(), this->m());
894
895 // Compute V^T v
896 work.resize(std::max(mm, nn));
897 gemv("N",
898 &nn,
899 &nn,
900 &alpha,
901 svd_vt->values.data(),
902 &nn,
903 v.data(),
904 &one,
905 &null,
906 work.data(),
907 &one);
908 // Multiply by singular values
909 for (size_type i = 0; i < wr.size(); ++i)
910 work[i] *= wr[i];
911 // Multiply with U
912 gemv("N",
913 &mm,
914 &mm,
915 &alpha,
916 svd_u->values.data(),
917 &mm,
918 work.data(),
919 &one,
920 &beta,
921 w.data(),
922 &one);
923 break;
925 default:
926 Assert(false, ExcState(state));
927 }
928}
929
930
932template <typename number>
933void
935 const Vector<number> &v) const
936{
937 vmult(w, v, true);
938}
939
940
941
942template <typename number>
943void
945 const Vector<number> &v) const
946{
947 Tvmult(w, v, true);
948}
949
950
951
952template <typename number>
953void
956 const bool adding) const
957{
960 Assert(C.state == matrix || C.state == inverse_matrix, ExcState(C.state));
961 Assert(this->n() == B.m(), ExcDimensionMismatch(this->n(), B.m()));
962 Assert(C.n() == B.n(), ExcDimensionMismatch(C.n(), B.n()));
963 Assert(C.m() == this->m(), ExcDimensionMismatch(this->m(), C.m()));
964 const types::blas_int mm = this->m();
965 const types::blas_int nn = B.n();
966 const types::blas_int kk = this->n();
967 const number alpha = 1.;
968 const number beta = (adding ? 1. : 0.);
969
970 gemm("N",
971 "N",
972 &mm,
973 &nn,
974 &kk,
975 &alpha,
976 this->values.data(),
977 &mm,
978 B.values.data(),
979 &kk,
980 &beta,
981 C.values.data(),
982 &mm);
983}
984
985
986
987template <typename number>
988void
991 const bool adding) const
992{
995 Assert(this->n() == B.m(), ExcDimensionMismatch(this->n(), B.m()));
996 Assert(C.n() == B.n(), ExcDimensionMismatch(C.n(), B.n()));
997 Assert(C.m() == this->m(), ExcDimensionMismatch(this->m(), C.m()));
998 const types::blas_int mm = this->m();
999 const types::blas_int nn = B.n();
1000 const types::blas_int kk = this->n();
1001 const number alpha = 1.;
1002 const number beta = (adding ? 1. : 0.);
1003
1004 // since FullMatrix stores the matrix in transposed order compared to this
1005 // matrix, compute B^T * A^T = (A * B)^T
1006 gemm("T",
1007 "T",
1008 &nn,
1009 &mm,
1010 &kk,
1011 &alpha,
1012 B.values.data(),
1013 &kk,
1014 this->values.data(),
1015 &mm,
1016 &beta,
1017 &C(0, 0),
1018 &nn);
1019}
1020
1021
1022
1023template <typename number>
1024void
1026 const LAPACKFullMatrix<number> &B,
1027 const Vector<number> &V,
1028 const bool adding) const
1029{
1032 Assert(C.state == matrix || C.state == inverse_matrix, ExcState(C.state));
1033
1034 const LAPACKFullMatrix<number> &A = *this;
1035
1036 Assert(A.m() == B.m(), ExcDimensionMismatch(A.m(), B.m()));
1037 Assert(C.n() == B.n(), ExcDimensionMismatch(C.n(), B.n()));
1038 Assert(C.m() == A.n(), ExcDimensionMismatch(A.n(), C.m()));
1039 Assert(V.size() == A.m(), ExcDimensionMismatch(A.m(), V.size()));
1040
1041 const types::blas_int mm = A.n();
1042 const types::blas_int nn = B.n();
1043 const types::blas_int kk = B.m();
1044
1045 // lapack does not have any triple product routines, including the case of
1046 // diagonal matrix in the middle, see
1047 // https://stackoverflow.com/questions/3548069/multiplying-three-matrices-in-blas-with-the-middle-one-being-diagonal
1048 // http://mathforum.org/kb/message.jspa?messageID=3546564
1049
1050 std::lock_guard<std::mutex> lock(mutex);
1051 // First, get V*B into "work" array
1052 work.resize(kk * nn);
1053 // following http://icl.cs.utk.edu/lapack-forum/viewtopic.php?f=2&t=768#p2577
1054 // do left-multiplication manually. Note that Xscal would require to first
1055 // copy the input vector as multiplication is done inplace.
1056 for (types::blas_int j = 0; j < nn; ++j)
1057 for (types::blas_int i = 0; i < kk; ++i)
1058 {
1059 Assert(j * kk + i < static_cast<types::blas_int>(work.size()),
1061 work[j * kk + i] = V(i) * B(i, j);
1062 }
1063
1064 // Now do the standard Tmmult:
1065 const number alpha = 1.;
1066 const number beta = (adding ? 1. : 0.);
1067
1068 gemm("T",
1069 "N",
1070 &mm,
1071 &nn,
1072 &kk,
1073 &alpha,
1074 this->values.data(),
1075 &kk,
1076 work.data(),
1077 &kk,
1078 &beta,
1079 C.values.data(),
1080 &mm);
1081}
1082
1083
1084
1085template <typename number>
1086void
1088{
1089 const LAPACKFullMatrix<number> &A = *this;
1090 AssertDimension(A.m(), B.n());
1091 AssertDimension(A.n(), B.m());
1092 const types::blas_int m = B.m();
1093 const types::blas_int n = B.n();
1094#ifdef DEAL_II_LAPACK_WITH_MKL
1095 const number one = 1.;
1096 omatcopy('C', 'C', n, m, one, A.values.data(), n, B.values.data(), m);
1097#else
1098 for (types::blas_int i = 0; i < m; ++i)
1099 for (types::blas_int j = 0; j < n; ++j)
1101#endif
1102}
1103
1104
1105
1106template <typename number>
1107void
1109{
1110 LAPACKFullMatrix<number> &A = *this;
1112 Assert(V.size() == A.m(), ExcDimensionMismatch(A.m(), V.size()));
1113
1114 const types::blas_int nn = A.n();
1115 const types::blas_int kk = A.m();
1116 for (types::blas_int j = 0; j < nn; ++j)
1117 for (types::blas_int i = 0; i < kk; ++i)
1118 A(i, j) *= V(i);
1119}
1120
1121
1122
1123template <typename number>
1124void
1126 const LAPACKFullMatrix<number> &B,
1127 const bool adding) const
1128{
1131 Assert(C.state == matrix || C.state == inverse_matrix, ExcState(C.state));
1132 Assert(this->m() == B.m(), ExcDimensionMismatch(this->m(), B.m()));
1133 Assert(C.n() == B.n(), ExcDimensionMismatch(C.n(), B.n()));
1134 Assert(C.m() == this->n(), ExcDimensionMismatch(this->n(), C.m()));
1135 const types::blas_int mm = this->n();
1136 const types::blas_int nn = B.n();
1137 const types::blas_int kk = B.m();
1138 const number alpha = 1.;
1139 const number beta = (adding ? 1. : 0.);
1140
1141 if (PointerComparison::equal(this, &B))
1142 {
1144 "T",
1145 &nn,
1146 &kk,
1147 &alpha,
1148 this->values.data(),
1149 &kk,
1150 &beta,
1151 C.values.data(),
1152 &nn);
1153
1154 // fill-in lower triangular part
1155 for (types::blas_int j = 0; j < nn; ++j)
1156 for (types::blas_int i = 0; i < j; ++i)
1157 C(j, i) = C(i, j);
1158
1159 C.property = symmetric;
1160 }
1161 else
1162 {
1163 gemm("T",
1164 "N",
1165 &mm,
1166 &nn,
1167 &kk,
1168 &alpha,
1169 this->values.data(),
1170 &kk,
1171 B.values.data(),
1172 &kk,
1173 &beta,
1174 C.values.data(),
1175 &mm);
1176 }
1177}
1178
1179
1180
1181template <typename number>
1182void
1184 const LAPACKFullMatrix<number> &B,
1185 const bool adding) const
1186{
1189 Assert(this->m() == B.m(), ExcDimensionMismatch(this->m(), B.m()));
1190 Assert(C.n() == B.n(), ExcDimensionMismatch(C.n(), B.n()));
1191 Assert(C.m() == this->n(), ExcDimensionMismatch(this->n(), C.m()));
1192 const types::blas_int mm = this->n();
1193 const types::blas_int nn = B.n();
1194 const types::blas_int kk = B.m();
1195 const number alpha = 1.;
1196 const number beta = (adding ? 1. : 0.);
1197
1198 // since FullMatrix stores the matrix in transposed order compared to this
1199 // matrix, compute B^T * A = (A^T * B)^T
1200 gemm("T",
1201 "N",
1202 &nn,
1203 &mm,
1204 &kk,
1205 &alpha,
1206 B.values.data(),
1207 &kk,
1208 this->values.data(),
1209 &kk,
1210 &beta,
1211 &C(0, 0),
1212 &nn);
1213}
1214
1215
1216
1217template <typename number>
1218void
1220 const LAPACKFullMatrix<number> &B,
1221 const bool adding) const
1222{
1225 Assert(C.state == matrix || C.state == inverse_matrix, ExcState(C.state));
1226 Assert(this->n() == B.n(), ExcDimensionMismatch(this->n(), B.n()));
1227 Assert(C.n() == B.m(), ExcDimensionMismatch(C.n(), B.m()));
1228 Assert(C.m() == this->m(), ExcDimensionMismatch(this->m(), C.m()));
1229 const types::blas_int mm = this->m();
1230 const types::blas_int nn = B.m();
1231 const types::blas_int kk = B.n();
1232 const number alpha = 1.;
1233 const number beta = (adding ? 1. : 0.);
1234
1235 if (PointerComparison::equal(this, &B))
1236 {
1238 "N",
1239 &nn,
1240 &kk,
1241 &alpha,
1242 this->values.data(),
1243 &nn,
1244 &beta,
1245 C.values.data(),
1246 &nn);
1247
1248 // fill-in lower triangular part
1249 for (types::blas_int j = 0; j < nn; ++j)
1250 for (types::blas_int i = 0; i < j; ++i)
1251 C(j, i) = C(i, j);
1252
1253 C.property = symmetric;
1254 }
1255 else
1256 {
1257 gemm("N",
1258 "T",
1259 &mm,
1260 &nn,
1261 &kk,
1262 &alpha,
1263 this->values.data(),
1264 &mm,
1265 B.values.data(),
1266 &nn,
1267 &beta,
1268 C.values.data(),
1269 &mm);
1270 }
1271}
1272
1273
1274
1275template <typename number>
1276void
1278 const LAPACKFullMatrix<number> &B,
1279 const bool adding) const
1280{
1283 Assert(this->n() == B.n(), ExcDimensionMismatch(this->n(), B.n()));
1284 Assert(C.n() == B.m(), ExcDimensionMismatch(C.n(), B.m()));
1285 Assert(C.m() == this->m(), ExcDimensionMismatch(this->m(), C.m()));
1286 const types::blas_int mm = this->m();
1287 const types::blas_int nn = B.m();
1288 const types::blas_int kk = B.n();
1289 const number alpha = 1.;
1290 const number beta = (adding ? 1. : 0.);
1291
1292 // since FullMatrix stores the matrix in transposed order compared to this
1293 // matrix, compute B * A^T = (A * B^T)^T
1294 gemm("N",
1295 "T",
1296 &nn,
1297 &mm,
1298 &kk,
1299 &alpha,
1300 B.values.data(),
1301 &nn,
1302 this->values.data(),
1303 &mm,
1304 &beta,
1305 &C(0, 0),
1306 &nn);
1307}
1308
1309
1310
1311template <typename number>
1312void
1314 const LAPACKFullMatrix<number> &B,
1315 const bool adding) const
1316{
1319 Assert(C.state == matrix || C.state == inverse_matrix, ExcState(C.state));
1320 Assert(this->m() == B.n(), ExcDimensionMismatch(this->m(), B.n()));
1321 Assert(C.n() == B.m(), ExcDimensionMismatch(C.n(), B.m()));
1322 Assert(C.m() == this->n(), ExcDimensionMismatch(this->n(), C.m()));
1323 const types::blas_int mm = this->n();
1324 const types::blas_int nn = B.m();
1325 const types::blas_int kk = B.n();
1326 const number alpha = 1.;
1327 const number beta = (adding ? 1. : 0.);
1328
1329 gemm("T",
1330 "T",
1331 &mm,
1332 &nn,
1333 &kk,
1334 &alpha,
1335 this->values.data(),
1336 &kk,
1337 B.values.data(),
1338 &nn,
1339 &beta,
1340 C.values.data(),
1341 &mm);
1342}
1343
1344
1345
1346template <typename number>
1347void
1349 const LAPACKFullMatrix<number> &B,
1350 const bool adding) const
1351{
1354 Assert(this->m() == B.n(), ExcDimensionMismatch(this->m(), B.n()));
1355 Assert(C.n() == B.m(), ExcDimensionMismatch(C.n(), B.m()));
1356 Assert(C.m() == this->n(), ExcDimensionMismatch(this->n(), C.m()));
1357 const types::blas_int mm = this->n();
1358 const types::blas_int nn = B.m();
1359 const types::blas_int kk = B.n();
1360 const number alpha = 1.;
1361 const number beta = (adding ? 1. : 0.);
1362
1363 // since FullMatrix stores the matrix in transposed order compared to this
1364 // matrix, compute B * A = (A^T * B^T)^T
1365 gemm("N",
1366 "N",
1367 &nn,
1368 &mm,
1369 &kk,
1370 &alpha,
1371 B.values.data(),
1372 &nn,
1373 this->values.data(),
1374 &kk,
1375 &beta,
1376 &C(0, 0),
1377 &nn);
1378}
1379
1380
1381
1382template <typename number>
1383void
1385{
1388
1389 const types::blas_int mm = this->m();
1390 const types::blas_int nn = this->n();
1391 number *const values = this->values.data();
1392 ipiv.resize(mm);
1393 types::blas_int info = 0;
1394 getrf(&mm, &nn, values, &mm, ipiv.data(), &info);
1395
1396 Assert(info >= 0, ExcInternalError());
1397
1398 // if info >= 0, the factorization has been completed
1399 state = lu;
1400
1402}
1403
1404
1405
1406template <typename number>
1407void
1409{
1410 property = p;
1411}
1412
1413
1414
1415template <typename number>
1416number
1418{
1419 const char type('O');
1420 return norm(type);
1421}
1422
1423
1424
1425template <typename number>
1426number
1428{
1429 const char type('I');
1430 return norm(type);
1431}
1432
1433
1434
1435template <typename number>
1436number
1438{
1439 const char type('F');
1440 return norm(type);
1441}
1442
1443
1444
1445template <typename number>
1446number
1448{
1449 std::lock_guard<std::mutex> lock(mutex);
1450
1453 ExcMessage("norms can be called in matrix state only."));
1454
1455 const types::blas_int N = this->n();
1456 const types::blas_int M = this->m();
1457 const number *const values = this->values.data();
1458 if (property == symmetric)
1459 {
1461 const types::blas_int lwork =
1462 (type == 'I' || type == 'O') ? std::max<types::blas_int>(1, N) : 0;
1463 work.resize(lwork);
1464 return lansy(&type, &LAPACKSupport::L, &N, values, &lda, work.data());
1465 }
1466 else
1467 {
1469 const types::blas_int lwork =
1470 (type == 'I') ? std::max<types::blas_int>(1, M) : 0;
1471 work.resize(lwork);
1472 return lange(&type, &M, &N, values, &lda, work.data());
1473 }
1474}
1475
1476
1477
1478template <typename number>
1479number
1481{
1484 ExcMessage("Trace can be called in matrix state only."));
1485 Assert(this->n() == this->m(), ExcDimensionMismatch(this->n(), this->m()));
1486
1487 number tr = 0;
1488 for (size_type i = 0; i < this->m(); ++i)
1489 tr += (*this)(i, i);
1490
1491 return tr;
1492}
1493
1494
1495
1496template <typename number>
1497void
1499{
1503
1504 const types::blas_int mm = this->m();
1505 const types::blas_int nn = this->n();
1506 Assert(mm == nn, ExcDimensionMismatch(mm, nn));
1507
1508 number *const values = this->values.data();
1509 types::blas_int info = 0;
1510 const types::blas_int lda = std::max<types::blas_int>(1, nn);
1511 potrf(&LAPACKSupport::L, &nn, values, &lda, &info);
1512
1513 // info < 0 : the info-th argument had an illegal value
1514 Assert(info >= 0, ExcInternalError());
1515
1516 state = cholesky;
1518}
1519
1520
1521
1522template <typename number>
1523number
1525{
1526 std::lock_guard<std::mutex> lock(mutex);
1528 number rcond = 0.;
1529
1530 const types::blas_int N = this->m();
1531 const number *values = this->values.data();
1532 types::blas_int info = 0;
1534 work.resize(3 * N);
1535 iwork.resize(N);
1536
1537 // use the same uplo as in Cholesky
1539 &N,
1540 values,
1541 &lda,
1542 &a_norm,
1543 &rcond,
1544 work.data(),
1545 iwork.data(),
1546 &info);
1547
1548 Assert(info >= 0, ExcInternalError());
1549
1550 return rcond;
1551}
1552
1553
1554
1555template <typename number>
1556number
1558{
1559 std::lock_guard<std::mutex> lock(mutex);
1562 number rcond = 0.;
1563
1564 const types::blas_int N = this->m();
1565 const number *const values = this->values.data();
1566 types::blas_int info = 0;
1568 work.resize(3 * N);
1569 iwork.resize(N);
1570
1571 const char norm = '1';
1572 const char diag = 'N';
1573 const char uplo =
1575 trcon(&norm,
1576 &uplo,
1577 &diag,
1578 &N,
1579 values,
1580 &lda,
1581 &rcond,
1582 work.data(),
1583 iwork.data(),
1584 &info);
1585
1586 Assert(info >= 0, ExcInternalError());
1587
1588 return rcond;
1589}
1590
1591
1592
1593template <typename number>
1594void
1596{
1599
1600 const types::blas_int mm = this->m();
1601 const types::blas_int nn = this->n();
1602 wr.resize(std::max(mm, nn));
1603 std::fill(wr.begin(), wr.end(), 0.);
1604 ipiv.resize(8 * mm);
1605
1606 svd_u = std::make_unique<LAPACKFullMatrix<number>>(mm, mm);
1607 svd_vt = std::make_unique<LAPACKFullMatrix<number>>(nn, nn);
1608 types::blas_int info = 0;
1609
1610 // First determine optimal workspace size
1611 work.resize(1);
1612 types::blas_int lwork = -1;
1613
1614 // TODO double check size
1615 std::vector<typename numbers::NumberTraits<number>::real_type> real_work;
1617 {
1618 // This array is only used by the complex versions.
1619 std::size_t min = std::min(this->m(), this->n());
1620 std::size_t max = std::max(this->m(), this->n());
1621 real_work.resize(
1622 std::max(5 * min * min + 5 * min, 2 * max * min + 2 * min * min + min));
1623 }
1624
1625 // make sure that the first entry in the work array is clear, in case the
1626 // routine does not completely overwrite the memory:
1627 work[0] = number();
1629 mm,
1630 nn,
1631 this->values,
1632 wr,
1633 svd_u->values,
1634 svd_vt->values,
1635 work,
1636 real_work,
1637 ipiv,
1638 lwork,
1639 info);
1640
1641 AssertThrow(info == 0, LAPACKSupport::ExcErrorCode("gesdd", info));
1642 // Resize the work array. Add one to the size computed by LAPACK to be on
1643 // the safe side.
1644 lwork = static_cast<types::blas_int>(std::abs(work[0]) + 1);
1645
1646 work.resize(lwork);
1647 // Do the actual SVD.
1649 mm,
1650 nn,
1651 this->values,
1652 wr,
1653 svd_u->values,
1654 svd_vt->values,
1655 work,
1656 real_work,
1657 ipiv,
1658 lwork,
1659 info);
1660 AssertThrow(info == 0, LAPACKSupport::ExcErrorCode("gesdd", info));
1661
1662 work.resize(0);
1663 ipiv.resize(0);
1664
1666}
1667
1668
1669
1670template <typename number>
1671void
1673{
1675 compute_svd();
1676
1678
1680 const double lim = std::abs(wr[0]) * threshold;
1681 for (size_type i = 0; i < wr.size(); ++i)
1682 {
1683 if (std::abs(wr[i]) > lim)
1684 wr[i] = one / wr[i];
1685 else
1686 wr[i] = 0.;
1687 }
1689}
1690
1691
1692
1693template <typename number>
1694void
1696 const unsigned int kernel_size)
1697{
1699 compute_svd();
1700
1702
1704 const unsigned int n_wr = wr.size();
1705 for (size_type i = 0; i < n_wr - kernel_size; ++i)
1706 wr[i] = one / wr[i];
1707 for (size_type i = n_wr - kernel_size; i < n_wr; ++i)
1708 wr[i] = 0.;
1710}
1711
1712
1713
1714template <typename number>
1715void
1717{
1719 const types::blas_int mm = this->m();
1720 const types::blas_int nn = this->n();
1721 Assert(nn == mm, ExcNotQuadratic());
1722
1723 number *const values = this->values.data();
1724 types::blas_int info = 0;
1725
1726 if (property != symmetric)
1727 {
1728 if (state == matrix)
1730
1731 ipiv.resize(mm);
1732 inv_work.resize(mm);
1733 getri(&mm, values, &mm, ipiv.data(), inv_work.data(), &mm, &info);
1734 }
1735 else
1736 {
1737 if (state == matrix)
1739
1740 const types::blas_int lda = std::max<types::blas_int>(1, nn);
1741 potri(&LAPACKSupport::L, &nn, values, &lda, &info);
1742 // inverse is stored in lower diagonal, set the upper diagonal
1743 // appropriately:
1744 for (types::blas_int i = 0; i < nn; ++i)
1745 for (types::blas_int j = i + 1; j < nn; ++j)
1746 this->el(i, j) = this->el(j, i);
1747 }
1748
1749 Assert(info >= 0, ExcInternalError());
1751
1753}
1754
1755
1756
1757template <typename number>
1758void
1759LAPACKFullMatrix<number>::solve(Vector<number> &v, const bool transposed) const
1760{
1761 Assert(this->m() == this->n(), LACExceptions::ExcNotQuadratic());
1762 AssertDimension(this->m(), v.size());
1763 const char *trans = transposed ? &T : &N;
1764 const types::blas_int nn = this->n();
1765 const number *const values = this->values.data();
1766 const types::blas_int n_rhs = 1;
1767 types::blas_int info = 0;
1768
1769 if (state == lu)
1770 {
1771 getrs(
1772 trans, &nn, &n_rhs, values, &nn, ipiv.data(), v.begin(), &nn, &info);
1773 }
1774 else if (state == cholesky)
1775 {
1776 potrs(&LAPACKSupport::L, &nn, &n_rhs, values, &nn, v.begin(), &nn, &info);
1777 }
1779 {
1780 const char uplo =
1782
1783 const types::blas_int lda = nn;
1784 const types::blas_int ldb = nn;
1785 trtrs(
1786 &uplo, trans, "N", &nn, &n_rhs, values, &lda, v.begin(), &ldb, &info);
1787 }
1788 else
1789 {
1790 Assert(false,
1791 ExcMessage(
1792 "The matrix has to be either factorized or triangular."));
1793 }
1794
1795 Assert(info == 0, ExcInternalError());
1796}
1797
1798
1799
1800template <typename number>
1801void
1803 const bool transposed) const
1804{
1805 Assert(B.state == matrix, ExcState(B.state));
1806
1807 Assert(this->m() == this->n(), LACExceptions::ExcNotQuadratic());
1808 AssertDimension(this->m(), B.m());
1809 const char *trans = transposed ? &T : &N;
1810 const types::blas_int nn = this->n();
1811 const number *const values = this->values.data();
1812 const types::blas_int n_rhs = B.n();
1813 types::blas_int info = 0;
1814
1815 if (state == lu)
1816 {
1817 getrs(trans,
1818 &nn,
1819 &n_rhs,
1820 values,
1821 &nn,
1822 ipiv.data(),
1823 B.values.data(),
1824 &nn,
1825 &info);
1826 }
1827 else if (state == cholesky)
1828 {
1830 &nn,
1831 &n_rhs,
1832 values,
1833 &nn,
1834 B.values.data(),
1835 &nn,
1836 &info);
1837 }
1839 {
1840 const char uplo =
1842
1843 const types::blas_int lda = nn;
1844 const types::blas_int ldb = nn;
1845 trtrs(&uplo,
1846 trans,
1847 "N",
1848 &nn,
1849 &n_rhs,
1850 values,
1851 &lda,
1852 B.values.data(),
1853 &ldb,
1854 &info);
1855 }
1856 else
1857 {
1858 Assert(false,
1859 ExcMessage(
1860 "The matrix has to be either factorized or triangular."));
1861 }
1862
1863 Assert(info == 0, ExcInternalError());
1864}
1865
1866
1867
1868template <typename number>
1869number
1871{
1872 Assert(this->m() == this->n(), LACExceptions::ExcNotQuadratic());
1873
1874 // LAPACK doesn't offer a function to compute a matrix determinant.
1875 // This is due to the difficulty in maintaining numerical accuracy, as the
1876 // calculations are likely to overflow or underflow. See
1877 // http://www.netlib.org/lapack/faq.html#_are_there_routines_in_lapack_to_compute_determinants
1878 //
1879 // However, after a PLU decomposition one can compute this by multiplication
1880 // of the diagonal entries with one another. One must take into consideration
1881 // the number of permutations (row swaps) stored in the P matrix.
1882 //
1883 // See the implementations in the blaze library (detNxN)
1884 // https://bitbucket.org/blaze-lib/blaze
1885 // and also
1886 // https://dualm.wordpress.com/2012/01/06/computing-determinant-in-fortran/
1887 // http://icl.cs.utk.edu/lapack-forum/viewtopic.php?p=341&#p336
1888 // for further information.
1890 Assert(ipiv.size() == this->m(), ExcInternalError());
1891 number det = 1.0;
1892 for (size_type i = 0; i < this->m(); ++i)
1893 det *=
1894 (ipiv[i] == types::blas_int(i + 1) ? this->el(i, i) : -this->el(i, i));
1895 return det;
1896}
1897
1898
1899
1900template <typename number>
1901void
1902LAPACKFullMatrix<number>::compute_eigenvalues(const bool right, const bool left)
1903{
1905 const types::blas_int nn = this->n();
1906 wr.resize(nn);
1907 wi.resize(nn);
1908 if (right)
1909 vr.resize(nn * nn);
1910 if (left)
1911 vl.resize(nn * nn);
1912
1913 types::blas_int info = 0;
1914 types::blas_int lwork = 1;
1915 const char jobvr = (right) ? V : N;
1916 const char jobvl = (left) ? V : N;
1917
1918 /*
1919 * The LAPACK routine xGEEV requires a sufficiently large work array; the
1920 * minimum requirement is
1921 *
1922 * work.size >= 4*nn.
1923 *
1924 * However, for better performance, a larger work array may be needed. The
1925 * first call determines the optimal work size and the second does the work.
1926 */
1927 lwork = -1;
1928 work.resize(1);
1929
1930 std::vector<typename numbers::NumberTraits<number>::real_type> real_work;
1932 // This array is only used by the complex versions.
1933 real_work.resize(2 * this->m());
1935 jobvr,
1936 this->values,
1937 this->m(),
1938 wr,
1939 wi,
1940 vl,
1941 vr,
1942 work,
1943 real_work,
1944 lwork,
1945 info);
1946
1947 // geev returns info=0 on success. Since we only queried the optimal size
1948 // for work, everything else would not be acceptable.
1949 Assert(info == 0, ExcInternalError());
1950 // Allocate working array according to suggestion (same strategy as was
1951 // noted in compute_svd).
1952 lwork = static_cast<types::blas_int>(std::abs(work[0]) + 1);
1953
1954 // resize workspace array
1955 work.resize(lwork);
1956 real_work.resize(lwork);
1957
1958 // Finally compute the eigenvalues.
1960 jobvr,
1961 this->values,
1962 this->m(),
1963 wr,
1964 wi,
1965 vl,
1966 vr,
1967 work,
1968 real_work,
1969 lwork,
1970 info);
1971
1972 Assert(info >= 0, ExcInternalError());
1973 if (info < 0)
1974 {
1975 AssertThrow(info == 0,
1976 ExcMessage("Lapack error in geev: the " +
1977 std::to_string(-info) +
1978 "-th"
1979 " parameter had an illegal value."));
1980 }
1981 else
1982 {
1984 info == 0,
1985 ExcMessage(
1986 "Lapack error in geev: the QR algorithm failed to compute "
1987 "all the eigenvalues, and no eigenvectors have been computed."));
1988 }
1989
1991}
1992
1993
1994
1995namespace
1996{
1997 // This function extracts complex eigenvectors from the underlying 'number'
1998 // array 'vr' of the LAPACK eigenvalue routine. For real-valued matrices
1999 // addressed by this function specialization, we might get complex
2000 // eigenvalues, which come in complex-conjugate pairs. In LAPACK, a compact
2001 // storage scheme is applied that stores the real and imaginary part of
2002 // eigenvectors only once, putting the real parts in one column and the
2003 // imaginary part in the next of a real-valued array. Here, we do the
2004 // unpacking into the usual complex values.
2005 template <typename RealNumber>
2006 void
2007 unpack_lapack_eigenvector_and_increment_index(
2008 const std::vector<RealNumber> &vr,
2009 const std::complex<RealNumber> &eigenvalue,
2010 FullMatrix<std::complex<RealNumber>> &result,
2011 unsigned int &index)
2012 {
2013 const std::size_t n = result.n();
2014 if (eigenvalue.imag() != 0.)
2015 {
2016 for (std::size_t j = 0; j < n; ++j)
2017 {
2018 result(j, index).real(vr[index * n + j]);
2019 result(j, index + 1).real(vr[index * n + j]);
2020 result(j, index).imag(vr[(index + 1) * n + j]);
2021 result(j, index + 1).imag(-vr[(index + 1) * n + j]);
2022 }
2023
2024 // we filled two columns with the complex-conjugate pair, so increment
2025 // returned index by 2
2026 index += 2;
2027 }
2028 else
2029 {
2030 for (unsigned int j = 0; j < n; ++j)
2031 result(j, index).real(vr[index * n + j]);
2032
2033 // real-valued case, we only filled one column
2034 ++index;
2035 }
2036 }
2037
2038 // This specialization fills the eigenvectors for complex-valued matrices,
2039 // in which case we simply read off the entry in the 'vr' array.
2040 template <typename ComplexNumber>
2041 void
2042 unpack_lapack_eigenvector_and_increment_index(
2043 const std::vector<ComplexNumber> &vr,
2044 const ComplexNumber &,
2046 unsigned int &index)
2047 {
2048 const std::size_t n = result.n();
2049 for (unsigned int j = 0; j < n; ++j)
2050 result(j, index) = vr[index * n + j];
2051
2052 // complex-valued case always only fills a single column
2053 ++index;
2054 }
2055} // namespace
2056
2057
2058
2059template <typename number>
2062{
2064 Assert(vr.size() == this->n_rows() * this->n_cols(),
2065 ExcMessage("Right eigenvectors are not available! Did you "
2066 "set the associated flag in compute_eigenvalues()?"));
2067
2069 result(n(), n());
2070
2071 for (unsigned int i = 0; i < n();)
2072 unpack_lapack_eigenvector_and_increment_index(vr, eigenvalue(i), result, i);
2073
2074 return result;
2075}
2076
2077
2078
2079template <typename number>
2082{
2084 Assert(vl.size() == this->n_rows() * this->n_cols(),
2085 ExcMessage("Left eigenvectors are not available! Did you "
2086 "set the associated flag in compute_eigenvalues()?"));
2087
2089 result(n(), n());
2090
2091 for (unsigned int i = 0; i < n();)
2092 unpack_lapack_eigenvector_and_increment_index(vl, eigenvalue(i), result, i);
2093
2094 return result;
2095}
2096
2097
2098
2099template <typename number>
2100void
2102 const number lower_bound,
2103 const number upper_bound,
2104 const number abs_accuracy,
2107{
2109 const types::blas_int nn = (this->n() > 0 ? this->n() : 1);
2110 Assert(static_cast<size_type>(nn) == this->m(), ExcNotQuadratic());
2111
2112 wr.resize(nn);
2113 LAPACKFullMatrix<number> matrix_eigenvectors(nn, nn);
2114
2115 number *const values_A = this->values.data();
2116 number *const values_eigenvectors = matrix_eigenvectors.values.data();
2117
2118 types::blas_int info(0), lwork(-1), n_eigenpairs(0);
2119 const char *const jobz(&V);
2120 const char *const uplo(&U);
2121 const char *const range(&V);
2122 const types::blas_int *const dummy(&one);
2123 std::vector<types::blas_int> iwork(static_cast<size_type>(5 * nn));
2124 std::vector<types::blas_int> ifail(static_cast<size_type>(nn));
2125
2126
2127 /*
2128 * The LAPACK routine xSYEVX requires a sufficiently large work array; the
2129 * minimum requirement is
2130 *
2131 * work.size >= 8*nn.
2132 *
2133 * However, for better performance, a larger work array may be needed. The
2134 * first call determines the optimal work size and the second does the work.
2135 */
2136 work.resize(1);
2137
2138 syevx(jobz,
2139 range,
2140 uplo,
2141 &nn,
2142 values_A,
2143 &nn,
2144 &lower_bound,
2145 &upper_bound,
2146 dummy,
2147 dummy,
2148 &abs_accuracy,
2149 &n_eigenpairs,
2150 wr.data(),
2151 values_eigenvectors,
2152 &nn,
2153 work.data(),
2154 &lwork,
2155 iwork.data(),
2156 ifail.data(),
2157 &info);
2158 // syevx returns info=0 on success. Since we only queried the optimal size
2159 // for work, everything else would not be acceptable.
2160 Assert(info == 0, ExcInternalError());
2161 // Allocate working array according to suggestion (same strategy as was noted
2162 // in compute_svd).
2163 lwork = static_cast<types::blas_int>(std::abs(work[0]) + 1);
2164 work.resize(static_cast<size_type>(lwork));
2165
2166 // Finally compute the eigenvalues.
2167 syevx(jobz,
2168 range,
2169 uplo,
2170 &nn,
2171 values_A,
2172 &nn,
2173 &lower_bound,
2174 &upper_bound,
2175 dummy,
2176 dummy,
2177 &abs_accuracy,
2178 &n_eigenpairs,
2179 wr.data(),
2180 values_eigenvectors,
2181 &nn,
2182 work.data(),
2183 &lwork,
2184 iwork.data(),
2185 ifail.data(),
2186 &info);
2187
2188 // Negative return value implies a wrong argument. This should be internal.
2189 Assert(info >= 0, ExcInternalError());
2190 if (info < 0)
2191 {
2192 AssertThrow(info == 0,
2193 ExcMessage("Lapack error in syevx: the " +
2194 std::to_string(-info) +
2195 "-th"
2196 " parameter had an illegal value."));
2197 }
2198 else if ((info > 0) && (info <= nn))
2199 {
2200 AssertThrow(info == 0,
2201 ExcMessage(
2202 "Lapack error in syevx: " + std::to_string(info) +
2203 " eigenvectors failed to converge."
2204 " (You may need to scale the abs_accuracy according"
2205 " to your matrix norm.)"));
2206 }
2207 else
2208 {
2209 AssertThrow(info == 0,
2210 ExcMessage("Lapack error in syevx: unknown error."));
2211 }
2212
2213 eigenvalues.reinit(n_eigenpairs);
2214 eigenvectors.reinit(nn, n_eigenpairs, true);
2215
2216 for (size_type i = 0; i < static_cast<size_type>(n_eigenpairs); ++i)
2217 {
2218 eigenvalues(i) = wr[i];
2219 size_type col_begin(i * nn);
2220 for (size_type j = 0; j < static_cast<size_type>(nn); ++j)
2221 {
2222 eigenvectors(j, i) = values_eigenvectors[col_begin + j];
2223 }
2224 }
2225
2227}
2228
2229
2230
2231template <typename number>
2232void
2235 const number lower_bound,
2236 const number upper_bound,
2237 const number abs_accuracy,
2239 std::vector<Vector<number>> &eigenvectors,
2240 const types::blas_int itype)
2241{
2243 const types::blas_int nn = (this->n() > 0 ? this->n() : 1);
2244 Assert(static_cast<size_type>(nn) == this->m(), ExcNotQuadratic());
2245 Assert(B.m() == B.n(), ExcNotQuadratic());
2246 Assert(static_cast<size_type>(nn) == B.n(), ExcDimensionMismatch(nn, B.n()));
2247
2248 wr.resize(nn);
2249 LAPACKFullMatrix<number> matrix_eigenvectors(nn, nn);
2250
2251 number *const values_A = this->values.data();
2252 number *const values_B = B.values.data();
2253 number *const values_eigenvectors = matrix_eigenvectors.values.data();
2254
2255 types::blas_int info(0), lwork(-1), n_eigenpairs(0);
2256 const char *const jobz(&V);
2257 const char *const uplo(&U);
2258 const char *const range(&V);
2259 const types::blas_int *const dummy(&one);
2260 iwork.resize(static_cast<size_type>(5 * nn));
2261 std::vector<types::blas_int> ifail(static_cast<size_type>(nn));
2262
2263
2264 /*
2265 * The LAPACK routine xSYGVX requires a sufficiently large work array; the
2266 * minimum requirement is
2267 *
2268 * work.size >= 8*nn.
2269 *
2270 * However, for better performance, a larger work array may be needed. The
2271 * first call determines the optimal work size and the second does the work.
2272 */
2273 work.resize(1);
2274
2275 sygvx(&itype,
2276 jobz,
2277 range,
2278 uplo,
2279 &nn,
2280 values_A,
2281 &nn,
2282 values_B,
2283 &nn,
2284 &lower_bound,
2285 &upper_bound,
2286 dummy,
2287 dummy,
2288 &abs_accuracy,
2289 &n_eigenpairs,
2290 wr.data(),
2291 values_eigenvectors,
2292 &nn,
2293 work.data(),
2294 &lwork,
2295 iwork.data(),
2296 ifail.data(),
2297 &info);
2298 // sygvx returns info=0 on success. Since we only queried the optimal size
2299 // for work, everything else would not be acceptable.
2300 Assert(info == 0, ExcInternalError());
2301 // Allocate working array according to suggestion (same strategy as was
2302 // noted in compute_svd).
2303 lwork = static_cast<types::blas_int>(std::abs(work[0]) + 1);
2304
2305 // resize workspace arrays
2306 work.resize(static_cast<size_type>(lwork));
2307
2308 // Finally compute the generalized eigenvalues.
2309 sygvx(&itype,
2310 jobz,
2311 range,
2312 uplo,
2313 &nn,
2314 values_A,
2315 &nn,
2316 values_B,
2317 &nn,
2318 &lower_bound,
2319 &upper_bound,
2320 dummy,
2321 dummy,
2322 &abs_accuracy,
2323 &n_eigenpairs,
2324 wr.data(),
2325 values_eigenvectors,
2326 &nn,
2327 work.data(),
2328 &lwork,
2329 iwork.data(),
2330 ifail.data(),
2331 &info);
2332
2333 // Negative return value implies a wrong argument. This should be internal.
2334 Assert(info >= 0, ExcInternalError());
2335 if (info < 0)
2336 {
2337 AssertThrow(info == 0,
2338 ExcMessage("Lapack error in sygvx: the " +
2339 std::to_string(-info) +
2340 "-th"
2341 " parameter had an illegal value."));
2342 }
2343 else if ((info > 0) && (info <= nn))
2344 {
2346 info == 0,
2347 ExcMessage(
2348 "Lapack error in sygvx: ssyevx/dsyevx failed to converge, and " +
2349 std::to_string(info) +
2350 " eigenvectors failed to converge."
2351 " (You may need to scale the abs_accuracy"
2352 " according to the norms of matrices A and B.)"));
2353 }
2354 else if ((info > nn) && (info <= 2 * nn))
2355 {
2356 AssertThrow(info == 0,
2357 ExcMessage(
2358 "Lapack error in sygvx: the leading minor of order " +
2359 std::to_string(info - nn) +
2360 " of matrix B is not positive-definite."
2361 " The factorization of B could not be completed and"
2362 " no eigenvalues or eigenvectors were computed."));
2363 }
2364 else
2365 {
2366 AssertThrow(info == 0,
2367 ExcMessage("Lapack error in sygvx: unknown error."));
2368 }
2369
2370 eigenvalues.reinit(n_eigenpairs);
2371 eigenvectors.resize(n_eigenpairs);
2372
2373 for (size_type i = 0; i < static_cast<size_type>(n_eigenpairs); ++i)
2374 {
2375 eigenvalues(i) = wr[i];
2376 size_type col_begin(i * nn);
2377 eigenvectors[i].reinit(nn, true);
2378 for (size_type j = 0; j < static_cast<size_type>(nn); ++j)
2379 {
2380 eigenvectors[i](j) = values_eigenvectors[col_begin + j];
2381 }
2382 }
2383
2385}
2386
2387
2388
2389template <typename number>
2390void
2393 std::vector<Vector<number>> &eigenvectors,
2394 const types::blas_int itype)
2395{
2397 const types::blas_int nn = this->n();
2398 Assert(static_cast<size_type>(nn) == this->m(), ExcNotQuadratic());
2399 Assert(B.m() == B.n(), ExcNotQuadratic());
2400 Assert(static_cast<size_type>(nn) == B.n(), ExcDimensionMismatch(nn, B.n()));
2401 Assert(eigenvectors.size() <= static_cast<size_type>(nn),
2402 ExcMessage("eigenvectors.size() > matrix.n()"));
2403
2404 wr.resize(nn);
2405 wi.resize(nn); // This is set purely for consistency reasons with the
2406 // eigenvalues() function.
2407
2408 number *const values_A = this->values.data();
2409 number *const values_B = B.values.data();
2410
2411 types::blas_int info = 0;
2412 types::blas_int lwork = -1;
2413 const char *const jobz = (eigenvectors.size() > 0) ? (&V) : (&N);
2414 const char *const uplo = (&U);
2415
2416 /*
2417 * The LAPACK routine xSYGV requires a sufficiently large work array; the
2418 * minimum requirement is
2419 *
2420 * work.size >= 3*nn - 1.
2421 *
2422 * However, for better performance, a larger work array may be needed. The
2423 * first call determines the optimal work size and the second does the work.
2424 */
2425 work.resize(1);
2426
2427 sygv(&itype,
2428 jobz,
2429 uplo,
2430 &nn,
2431 values_A,
2432 &nn,
2433 values_B,
2434 &nn,
2435 wr.data(),
2436 work.data(),
2437 &lwork,
2438 &info);
2439 // sygv returns info=0 on success. Since we only queried the optimal size
2440 // for work, everything else would not be acceptable.
2441 Assert(info == 0, ExcInternalError());
2442 // Allocate working array according to suggestion (same strategy as was
2443 // noted in compute_svd).
2444 lwork = static_cast<types::blas_int>(std::abs(work[0]) + 1);
2445
2446 // resize workspace array
2447 work.resize(static_cast<size_type>(lwork));
2448
2449 // Finally compute the generalized eigenvalues.
2450 sygv(&itype,
2451 jobz,
2452 uplo,
2453 &nn,
2454 values_A,
2455 &nn,
2456 values_B,
2457 &nn,
2458 wr.data(),
2459 work.data(),
2460 &lwork,
2461 &info);
2462 // Negative return value implies a wrong argument. This should be internal.
2463
2464 Assert(info >= 0, ExcInternalError());
2465 if (info < 0)
2466 {
2467 AssertThrow(info == 0,
2468 ExcMessage("Lapack error in sygv: the " +
2469 std::to_string(-info) +
2470 "-th"
2471 " parameter had an illegal value."));
2472 }
2473 else if ((info > 0) && (info <= nn))
2474 {
2476 info == 0,
2477 ExcMessage(
2478 "Lapack error in sygv: ssyev/dsyev failed to converge, and " +
2479 std::to_string(info) +
2480 " off-diagonal elements of an intermediate "
2481 " tridiagonal did not converge to zero."
2482 " (You may need to scale the abs_accuracy"
2483 " according to the norms of matrices A and B.)"));
2484 }
2485 else if ((info > nn) && (info <= 2 * nn))
2486 {
2487 AssertThrow(info == 0,
2488 ExcMessage(
2489 "Lapack error in sygv: the leading minor of order " +
2490 std::to_string(info - nn) +
2491 " of matrix B is not positive-definite."
2492 " The factorization of B could not be completed and"
2493 " no eigenvalues or eigenvectors were computed."));
2494 }
2495 else
2496 {
2497 AssertThrow(info == 0,
2498 ExcMessage("Lapack error in sygv: unknown error."));
2499 }
2500
2501 for (size_type i = 0; i < eigenvectors.size(); ++i)
2502 {
2503 size_type col_begin(i * nn);
2504 eigenvectors[i].reinit(nn, true);
2505 for (size_type j = 0; j < static_cast<size_type>(nn); ++j)
2506 {
2507 eigenvectors[i](j) = values_A[col_begin + j];
2508 }
2509 }
2511}
2512
2513
2514
2515template <typename number>
2516void
2518 const unsigned int precision,
2519 const bool scientific,
2520 const unsigned int width_,
2521 const char *zero_string,
2522 const double denominator,
2523 const double threshold,
2524 const char *separator) const
2525{
2526 unsigned int width = width_;
2527
2528 Assert((!this->empty()) || (this->n() + this->m() == 0), ExcInternalError());
2532 ExcState(state));
2533
2534 // set output format, but store old
2535 // state
2536 std::ios::fmtflags old_flags = out.flags();
2537 std::streamsize old_precision = out.precision(precision);
2538
2539 if (scientific)
2540 {
2541 out.setf(std::ios::scientific, std::ios::floatfield);
2542 if (width == 0u)
2543 width = precision + 7;
2544 }
2545 else
2546 {
2547 out.setf(std::ios::fixed, std::ios::floatfield);
2548 if (width == 0u)
2549 width = precision + 2;
2550 }
2551
2552 for (size_type i = 0; i < this->m(); ++i)
2553 {
2554 // Cholesky is stored in lower triangular, so just output this part:
2555 const size_type nc = state == LAPACKSupport::cholesky ? i + 1 : this->n();
2556 for (size_type j = 0; j < nc; ++j)
2557 // we might have complex numbers, so use abs also to check for nan
2558 // since there is no isnan on complex numbers
2559 if (numbers::is_nan(std::abs((*this)(i, j))))
2560 out << std::setw(width) << (*this)(i, j) << separator;
2561 else if (std::abs(this->el(i, j)) > threshold)
2562 out << std::setw(width) << this->el(i, j) * denominator << separator;
2563 else
2564 out << std::setw(width) << zero_string << separator;
2565 out << std::endl;
2566 }
2567
2568 AssertThrow(out.fail() == false, ExcIO());
2569 // reset output format
2570 out.flags(old_flags);
2571 out.precision(old_precision);
2572}
2573
2574
2575
2576template <typename number>
2579{
2580 return this->state;
2581}
2582
2583
2584//----------------------------------------------------------------------//
2585
2586template <typename number>
2587void
2589{
2590 matrix = &M;
2591 mem = nullptr;
2592}
2593
2594
2595template <typename number>
2596void
2603
2604
2605template <typename number>
2606void
2608 const Vector<number> &src) const
2609{
2610 dst = src;
2611 matrix->solve(dst, false);
2612}
2613
2614
2615template <typename number>
2616void
2618 const Vector<number> &src) const
2619{
2620 dst = src;
2621 matrix->solve(dst, true);
2622}
2623
2624
2625template <typename number>
2626void
2628 const BlockVector<number> &src) const
2629{
2630 Assert(mem != nullptr, ExcNotInitialized());
2631 Vector<number> *aux = mem->alloc();
2632 *aux = src;
2633 matrix->solve(*aux, false);
2634 dst = *aux;
2635}
2636
2637
2638template <typename number>
2639void
2641 const BlockVector<number> &src) const
2642{
2643 Assert(mem != nullptr, ExcNotInitialized());
2644 Vector<number> *aux = mem->alloc();
2645 *aux = src;
2646 matrix->solve(*aux, true);
2647 dst = *aux;
2648}
2649
2650
2651
2652#include "lac/lapack_full_matrix.inst"
2653
2654
void omatcopy(char, char, ::types::blas_int, ::types::blas_int, const number1, const number2 *, ::types::blas_int, number3 *, ::types::blas_int)
pointer data()
EnableObserverPointer & operator=(const EnableObserverPointer &)
size_type n() const
size_type m() const
LAPACKFullMatrix< number > & operator*=(const number factor)
number reciprocal_condition_number() const
void Tmmult(LAPACKFullMatrix< number > &C, const LAPACKFullMatrix< number > &B, const bool adding=false) const
void scale_rows(const Vector< number > &V)
FullMatrix< std::complex< typename numbers::NumberTraits< number >::real_type > > get_right_eigenvectors() const
void add(const number a, const LAPACKFullMatrix< number > &B)
void Tvmult(Vector< number2 > &w, const Vector< number2 > &v, const bool adding=false) const
void transpose(LAPACKFullMatrix< number > &B) const
void mTmult(LAPACKFullMatrix< number > &C, const LAPACKFullMatrix< number > &B, const bool adding=false) const
std::vector< typename numbers::NumberTraits< typename MatrixType::value_type >::real_type > wr
LAPACKSupport::State get_state() const
void compute_eigenvalues_symmetric(const number lower_bound, const number upper_bound, const number abs_accuracy, Vector< number > &eigenvalues, FullMatrix< number > &eigenvectors)
void vmult(Vector< number2 > &w, const Vector< number2 > &v, const bool adding=false) const
void reinit(const size_type size)
LAPACKFullMatrix< number > & operator=(const LAPACKFullMatrix< number > &)
FullMatrix< std::complex< typename numbers::NumberTraits< number >::real_type > > get_left_eigenvectors() const
std::unique_ptr< LAPACKFullMatrix< typename MatrixType::value_type > > svd_vt
std::vector< typename MatrixType::value_type > work
void grow_or_shrink(const size_type size)
void apply_givens_rotation(const std::array< number, 3 > &csr, const size_type i, const size_type k, const bool left=true)
void set_property(const LAPACKSupport::Property property)
std::complex< typename numbers::NumberTraits< number >::real_type > eigenvalue(const size_type i) const
number norm(const char type) const
void solve(Vector< number > &v, const bool transposed=false) const
void compute_eigenvalues(const bool right_eigenvectors=false, const bool left_eigenvectors=false)
LAPACKSupport::State state
std::make_unsigned_t< types::blas_int > size_type
std::vector< number > inv_work
number frobenius_norm() const
LAPACKFullMatrix(const size_type size=0)
LAPACKSupport::Property property
std::vector< number > wi
size_type m() const
void compute_inverse_svd(const double threshold=0.)
void compute_generalized_eigenvalues_symmetric(LAPACKFullMatrix< number > &B, const number lower_bound, const number upper_bound, const number abs_accuracy, Vector< number > &eigenvalues, std::vector< Vector< number > > &eigenvectors, const types::blas_int itype=1)
void vmult_add(Vector< number2 > &w, const Vector< number2 > &v) const
size_type n() const
number linfty_norm() const
void TmTmult(LAPACKFullMatrix< number > &C, const LAPACKFullMatrix< number > &B, const bool adding=false) const
std::unique_ptr< LAPACKFullMatrix< typename MatrixType::value_type > > svd_u
std::vector< typename MatrixType::value_type > vr
void compute_inverse_svd_with_kernel(const unsigned int kernel_size)
void Tvmult_add(Vector< number2 > &w, const Vector< number2 > &v) const
std::vector< types::blas_int > iwork
void rank1_update(const number a, const Vector< number > &v)
std::vector< types::blas_int > ipiv
void remove_row_and_column(const size_type row, const size_type col)
LAPACKFullMatrix< number > & operator/=(const number factor)
number determinant() const
std::vector< typename MatrixType::value_type > vl
void mmult(LAPACKFullMatrix< number > &C, const LAPACKFullMatrix< number > &B, const bool adding=false) const
void print_formatted(std::ostream &out, const unsigned int precision=3, const bool scientific=true, const unsigned int width=0, const char *zero_string=" ", const double denominator=1., const double threshold=0., const char *separator=" ") const
void vmult(Vector< number > &, const Vector< number > &) const
void initialize(const LAPACKFullMatrix< number > &)
void Tvmult(Vector< number > &, const Vector< number > &) const
ObserverPointer< VectorMemory< Vector< number > >, PreconditionLU< number > > mem
ObserverPointer< const LAPACKFullMatrix< number >, PreconditionLU< number > > matrix
number el(const size_type i, const size_type j) const
size_type n() const
size_type m() const
AlignedVector< number > values
Definition table.h:795
size_type n_elements() const
void reinit(const size_type size1, const size_type size2, const bool omit_default_initialization=false)
reference el(const size_type i, const size_type j)
pointer data()
virtual size_type size() const override
iterator begin()
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:40
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:41
static ::ExceptionBase & ExcIO()
static ::ExceptionBase & ExcZero()
static ::ExceptionBase & ExcNotImplemented()
static ::ExceptionBase & ExcErrorCode(std::string arg1, types::blas_int arg2)
static ::ExceptionBase & ExcScalarAssignmentOnlyForZeroValue()
static ::ExceptionBase & ExcProperty(Property arg1)
#define Assert(cond, exc)
static ::ExceptionBase & ExcSingular()
#define AssertIsFinite(number)
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcNotInitialized()
static ::ExceptionBase & ExcNotQuadratic()
static ::ExceptionBase & ExcInvalidState()
static ::ExceptionBase & ExcMessage(std::string arg1)
static ::ExceptionBase & ExcState(State arg1)
#define AssertThrow(cond, exc)
void getrs(const char *, const ::types::blas_int *, const ::types::blas_int *, const number1 *, const ::types::blas_int *, const ::types::blas_int *, number2 *, const ::types::blas_int *, ::types::blas_int *)
void syrk(const char *, const char *, const ::types::blas_int *, const ::types::blas_int *, const number1 *, const number2 *, const ::types::blas_int *, const number3 *, number4 *, const ::types::blas_int *)
void geev(const char *, const char *, const ::types::blas_int *, number1 *, const ::types::blas_int *, number2 *, number3 *, number4 *, const ::types::blas_int *, number5 *, const ::types::blas_int *, number6 *, const ::types::blas_int *, ::types::blas_int *)
void gemm(const char *, const char *, const ::types::blas_int *, const ::types::blas_int *, const ::types::blas_int *, const number1 *, const number2 *, const ::types::blas_int *, const number3 *, const ::types::blas_int *, const number4 *, number5 *, const ::types::blas_int *)
void pocon(const char *, const ::types::blas_int *, const number1 *, const ::types::blas_int *, const number2 *, number3 *, number4 *, ::types::blas_int *, ::types::blas_int *)
void lascl(const char *, const ::types::blas_int *, const ::types::blas_int *, const number1 *, const number2 *, const ::types::blas_int *, const ::types::blas_int *, number3 *, const ::types::blas_int *, ::types::blas_int *)
void syr(const char *, const ::types::blas_int *, const number1 *, const number2 *, const ::types::blas_int *, number3 *, const ::types::blas_int *)
void trtrs(const char *, const char *, const char *, const ::types::blas_int *, const ::types::blas_int *, const number1 *, const ::types::blas_int *, number2 *, const ::types::blas_int *, ::types::blas_int *)
void syevx(const char *, const char *, const char *, const ::types::blas_int *, number1 *, const ::types::blas_int *, const number2 *, const number3 *, const ::types::blas_int *, const ::types::blas_int *, const number4 *, ::types::blas_int *, number5 *, number6 *, const ::types::blas_int *, number7 *, const ::types::blas_int *, ::types::blas_int *, ::types::blas_int *, ::types::blas_int *)
void trcon(const char *, const char *, const char *, const ::types::blas_int *, const number1 *, const ::types::blas_int *, number2 *, number3 *, ::types::blas_int *, ::types::blas_int *)
void gesdd(const char *, const ::types::blas_int *, const ::types::blas_int *, number1 *, const ::types::blas_int *, number2 *, number3 *, const ::types::blas_int *, number4 *, const ::types::blas_int *, number5 *, const ::types::blas_int *, ::types::blas_int *, ::types::blas_int *)
void axpy(const ::types::blas_int *, const number1 *, const number2 *, const ::types::blas_int *, number3 *, const ::types::blas_int *)
void trmv(const char *, const char *, const char *, const ::types::blas_int *, const number1 *, const ::types::blas_int *, number2 *, const ::types::blas_int *)
void sygvx(const ::types::blas_int *, const char *, const char *, const char *, const ::types::blas_int *, number1 *, const ::types::blas_int *, number2 *, const ::types::blas_int *, const number3 *, const number4 *, const ::types::blas_int *, const ::types::blas_int *, const number5 *, ::types::blas_int *, number6 *, number7 *, const ::types::blas_int *, number8 *, const ::types::blas_int *, ::types::blas_int *, ::types::blas_int *, ::types::blas_int *)
void gemv(const char *, const ::types::blas_int *, const ::types::blas_int *, const number1 *, const number2 *, const ::types::blas_int *, const number3 *, const ::types::blas_int *, const number4 *, number5 *, const ::types::blas_int *)
number1 lange(const char *, const ::types::blas_int *, const ::types::blas_int *, const number1 *, const ::types::blas_int *, number2 *)
void potrs(const char *, const ::types::blas_int *, const ::types::blas_int *, const number1 *, const ::types::blas_int *, number2 *, const ::types::blas_int *, ::types::blas_int *)
void sygv(const ::types::blas_int *, const char *, const char *, const ::types::blas_int *, number1 *, const ::types::blas_int *, number2 *, const ::types::blas_int *, number3 *, number4 *, const ::types::blas_int *, ::types::blas_int *)
void potri(const char *, const ::types::blas_int *, number1 *, const ::types::blas_int *, ::types::blas_int *)
number1 lansy(const char *, const char *, const ::types::blas_int *, const number1 *, const ::types::blas_int *, number2 *)
void potrf(const char *, const ::types::blas_int *, number1 *, const ::types::blas_int *, ::types::blas_int *)
void getrf(const ::types::blas_int *, const ::types::blas_int *, number1 *, const ::types::blas_int *, ::types::blas_int *, ::types::blas_int *)
void getri(const ::types::blas_int *, number1 *, const ::types::blas_int *, const ::types::blas_int *, number2 *, const ::types::blas_int *, ::types::blas_int *)
@ cholesky
Contents is a Cholesky decomposition.
@ lu
Contents is an LU decomposition.
@ matrix
Contents is actually a matrix.
@ unusable
Contents is something useless.
@ inverse_matrix
Contents is the inverse of a matrix.
@ svd
Matrix contains singular value decomposition,.
@ inverse_svd
Matrix is the inverse of a singular value decomposition.
@ eigenvalues
Eigenvalue vector is filled.
@ symmetric
Matrix is symmetric.
@ upper_triangular
Matrix is upper triangular.
@ lower_triangular
Matrix is lower triangular.
@ general
No special properties.
constexpr char L
constexpr char N
constexpr char U
constexpr char T
constexpr char V
constexpr char A
constexpr types::blas_int one
std::array< NumberType, 3 > givens_rotation(const NumberType &x, const NumberType &y)
std::array< NumberType, 3 > hyperbolic_rotation(const NumberType &x, const NumberType &y)
void gesdd_helper(const char job, const types::blas_int n_rows, const types::blas_int n_cols, AlignedVector< T > &matrix, std::vector< T > &singular_values, AlignedVector< T > &left_vectors, AlignedVector< T > &right_vectors, std::vector< T > &real_work, std::vector< T > &, std::vector< types::blas_int > &integer_work, const types::blas_int work_flag, types::blas_int &info)
void geev_helper(const char vl, const char vr, AlignedVector< T > &matrix, const types::blas_int n_rows, std::vector< T > &real_part_eigenvalues, std::vector< T > &imag_part_eigenvalues, std::vector< T > &left_eigenvectors, std::vector< T > &right_eigenvectors, std::vector< T > &real_work, std::vector< T > &, const types::blas_int work_flag, types::blas_int &info)
bool is_nan(const double x)
Definition numbers.h:500
::VectorizedArray< Number, width > min(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)
int blas_int
static bool equal(const T *p1, const T *p2)
static constexpr const number & conjugate(const number &x)
Definition numbers.h:540
static constexpr bool is_complex
Definition numbers.h:404
std::array< std::pair< Number, Tensor< 1, dim, Number > >, dim > eigenvectors(const SymmetricTensor< 2, dim, Number > &T, const SymmetricTensorEigenvectorMethod method=SymmetricTensorEigenvectorMethod::ql_implicit_shifts)