Dune Core Modules (2.6.0)

arpackpp.hh
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3#ifndef DUNE_ISTL_EIGENVALUE_ARPACKPP_HH
4#define DUNE_ISTL_EIGENVALUE_ARPACKPP_HH
5
6#if HAVE_ARPACKPP || defined DOXYGEN
7
8#include <cmath> // provides std::abs, std::pow, std::sqrt
9
10#include <iostream> // provides std::cout, std::endl
11#include <string> // provides std::string
12
13#include <dune/common/fvector.hh> // provides Dune::FieldVector
14#include <dune/common/exceptions.hh> // provides DUNE_THROW(...)
15
16#include <dune/istl/bvector.hh> // provides Dune::BlockVector
17#include <dune/istl/istlexception.hh> // provides Dune::ISTLError
18#include <dune/istl/io.hh> // provides Dune::printvector(...)
19
20#ifdef Status
21#undef Status // prevent preprocessor from damaging the ARPACK++
22 // code when "X11/Xlib.h" is included (the latter
23 // defines Status as "#define Status int" and
24 // ARPACK++ provides a class with a method called
25 // Status)
26#endif
27#include "arssym.h" // provides ARSymStdEig
28
29namespace Dune
30{
31
36 namespace Impl {
52 template <class BCRSMatrix>
53 class ArPackPlusPlus_BCRSMatrixWrapper
54 {
55 public:
57 typedef typename BCRSMatrix::field_type Real;
58
59 public:
61 ArPackPlusPlus_BCRSMatrixWrapper (const BCRSMatrix& A)
62 : A_(A),
63 m_(A_.M() * mBlock), n_(A_.N() * nBlock)
64 {
65 // assert that BCRSMatrix type has blocklevel 2
66 static_assert
68 "Only BCRSMatrices with blocklevel 2 are supported.");
69
70 // allocate memory for auxiliary block vector objects
71 // which are compatible to matrix rows / columns
72 domainBlockVector.resize(A_.N(),false);
73 rangeBlockVector.resize(A_.M(),false);
74 }
75
77 inline void multMv (Real* v, Real* w)
78 {
79 // get vector v as an object of appropriate type
80 arrayToDomainBlockVector(v,domainBlockVector);
81
82 // perform matrix-vector product
83 A_.mv(domainBlockVector,rangeBlockVector);
84
85 // get vector w from object of appropriate type
86 rangeBlockVectorToArray(rangeBlockVector,w);
87 };
88
90 inline void multMtMv (Real* v, Real* w)
91 {
92 // get vector v as an object of appropriate type
93 arrayToDomainBlockVector(v,domainBlockVector);
94
95 // perform matrix-vector product
96 A_.mv(domainBlockVector,rangeBlockVector);
97 A_.mtv(rangeBlockVector,domainBlockVector);
98
99 // get vector w from object of appropriate type
100 domainBlockVectorToArray(domainBlockVector,w);
101 };
102
104 inline void multMMtv (Real* v, Real* w)
105 {
106 // get vector v as an object of appropriate type
107 arrayToRangeBlockVector(v,rangeBlockVector);
108
109 // perform matrix-vector product
110 A_.mtv(rangeBlockVector,domainBlockVector);
111 A_.mv(domainBlockVector,rangeBlockVector);
112
113 // get vector w from object of appropriate type
114 rangeBlockVectorToArray(rangeBlockVector,w);
115 };
116
118 inline int nrows () const { return m_; }
119
121 inline int ncols () const { return n_; }
122
123 protected:
124 // Number of rows and columns in each block of the matrix
125 constexpr static int mBlock = BCRSMatrix::block_type::rows;
126 constexpr static int nBlock = BCRSMatrix::block_type::cols;
127
128 // Type of vectors in the domain of the linear map associated with
129 // the matrix, i.e. block vectors compatible to matrix rows
130 constexpr static int dbvBlockSize = nBlock;
131 typedef Dune::FieldVector<Real,dbvBlockSize> DomainBlockVectorBlock;
132 typedef Dune::BlockVector<DomainBlockVectorBlock> DomainBlockVector;
133
134 // Type of vectors in the range of the linear map associated with
135 // the matrix, i.e. block vectors compatible to matrix columns
136 constexpr static int rbvBlockSize = mBlock;
137 typedef Dune::FieldVector<Real,rbvBlockSize> RangeBlockVectorBlock;
138 typedef Dune::BlockVector<RangeBlockVectorBlock> RangeBlockVector;
139
140 // Types for vector index access
141 typedef typename DomainBlockVector::size_type dbv_size_type;
142 typedef typename RangeBlockVector::size_type rbv_size_type;
143 typedef typename DomainBlockVectorBlock::size_type dbvb_size_type;
144 typedef typename RangeBlockVectorBlock::size_type rbvb_size_type;
145
146 // Get vector v from a block vector object which is compatible to
147 // matrix rows
148 static inline void
149 domainBlockVectorToArray (const DomainBlockVector& dbv, Real* v)
150 {
151 for (dbv_size_type block = 0; block < dbv.N(); ++block)
152 for (dbvb_size_type iBlock = 0; iBlock < dbvBlockSize; ++iBlock)
153 v[block*dbvBlockSize + iBlock] = dbv[block][iBlock];
154 }
155
156 // Get vector v from a block vector object which is compatible to
157 // matrix columns
158 static inline void
159 rangeBlockVectorToArray (const RangeBlockVector& rbv, Real* v)
160 {
161 for (rbv_size_type block = 0; block < rbv.N(); ++block)
162 for (rbvb_size_type iBlock = 0; iBlock < rbvBlockSize; ++iBlock)
163 v[block*rbvBlockSize + iBlock] = rbv[block][iBlock];
164 }
165
166 public:
169 static inline void arrayToDomainBlockVector (const Real* v,
170 DomainBlockVector& dbv)
171 {
172 for (dbv_size_type block = 0; block < dbv.N(); ++block)
173 for (dbvb_size_type iBlock = 0; iBlock < dbvBlockSize; ++iBlock)
174 dbv[block][iBlock] = v[block*dbvBlockSize + iBlock];
175 }
176
179 static inline void arrayToRangeBlockVector (const Real* v,
180 RangeBlockVector& rbv)
181 {
182 for (rbv_size_type block = 0; block < rbv.N(); ++block)
183 for (rbvb_size_type iBlock = 0; iBlock < rbvBlockSize; ++iBlock)
184 rbv[block][iBlock] = v[block*rbvBlockSize + iBlock];
185 }
186
187 protected:
188 // The DUNE-ISTL BCRSMatrix
189 const BCRSMatrix& A_;
190
191 // Number of rows and columns in the matrix
192 const int m_, n_;
193
194 // Auxiliary block vector objects which are
195 // compatible to matrix rows / columns
196 mutable DomainBlockVector domainBlockVector;
197 mutable RangeBlockVector rangeBlockVector;
198 };
199 } // end namespace Impl
200
240 template <typename BCRSMatrix, typename BlockVector>
242 {
243 public:
244 typedef typename BlockVector::field_type Real;
245
246 public:
266 const unsigned int nIterationsMax = 100000,
267 const unsigned int verbosity_level = 0)
268 : m_(m), nIterationsMax_(nIterationsMax),
269 verbosity_level_(verbosity_level),
270 nIterations_(0),
271 title_(" ArPackPlusPlus_Algorithms: "),
272 blank_(title_.length(),' ')
273 {}
274
286 inline void computeSymMaxMagnitude (const Real& epsilon,
287 BlockVector& x, Real& lambda) const
288 {
289 // print verbosity information
290 if (verbosity_level_ > 0)
291 std::cout << title_ << "Computing an approximation of "
292 << "the dominant eigenvalue of a matrix which "
293 << "is assumed to be symmetric." << std::endl;
294
295 // use type ArPackPlusPlus_BCRSMatrixWrapper to store matrix information
296 // and to perform the product A*v (LU decomposition is not used)
297 typedef Impl::ArPackPlusPlus_BCRSMatrixWrapper<BCRSMatrix> WrappedMatrix;
298 WrappedMatrix A(m_);
299
300 // get number of rows and columns in A
301 const int nrows = A.nrows();
302 const int ncols = A.ncols();
303
304 // assert that A is square
305 if (nrows != ncols)
306 DUNE_THROW(Dune::ISTLError,"Matrix is not square ("
307 << nrows << "x" << ncols << ").");
308
309 // allocate memory for variables, set parameters
310 const int nev = 1; // Number of eigenvalues to compute
311 int ncv = std::min(20, nrows); // Number of Arnoldi vectors generated at each iteration (0 == auto)
312 const Real tol = epsilon; // Stopping tolerance (relative accuracy of Ritz values) (0 == machine precision)
313 const int maxit = nIterationsMax_*nev; // Maximum number of Arnoldi update iterations allowed (0 == 100*nev)
314 Real* ev = new Real[nev]; // Computed eigenvalues of A
315 const bool ivec = true; // Flag deciding if eigenvectors shall be determined
316 int nconv; // Number of converged eigenvalues
317
318 // define what we need: eigenvalues with largest magnitude
319 char which[] = "LM";
320 ARSymStdEig<Real,WrappedMatrix>
321 dprob(nrows, nev, &A, &WrappedMatrix::multMv, which, ncv, tol, maxit);
322
323 // set ARPACK verbosity mode if requested
324 if (verbosity_level_ > 3) dprob.Trace();
325
326 // find eigenvalues and eigenvectors of A, obtain the eigenvalues
327 nconv = dprob.Eigenvalues(ev,ivec);
328
329 // obtain approximated dominant eigenvalue of A
330 lambda = ev[nev-1];
331
332 // obtain associated approximated eigenvector of A
333 Real* x_raw = dprob.RawEigenvector(nev-1);
334 WrappedMatrix::arrayToDomainBlockVector(x_raw,x);
335
336 // obtain number of Arnoldi update iterations actually taken
337 nIterations_ = dprob.GetIter();
338
339 // compute residual norm
340 BlockVector r(x);
341 Real* Ax_raw = new Real[nrows];
342 A.multMv(x_raw,Ax_raw);
343 WrappedMatrix::arrayToDomainBlockVector(Ax_raw,r);
344 r.axpy(-lambda,x);
345 const Real r_norm = r.two_norm();
346
347 // print verbosity information
348 if (verbosity_level_ > 0)
349 {
350 if (verbosity_level_ > 1)
351 {
352 // print some information about the problem
353 std::cout << blank_ << "Obtained eigenvalues of A by solving "
354 << "A*x = λ*x using the ARPACK++ class ARSym"
355 << "StdEig:" << std::endl;
356 std::cout << blank_ << " converged eigenvalues of A: "
357 << nconv << " / " << nev << std::endl;
358 std::cout << blank_ << " dominant eigenvalue of A: "
359 << lambda << std::endl;
360 }
361 std::cout << blank_ << "Result ("
362 << "#iterations = " << nIterations_ << ", "
363 << "║residual║_2 = " << r_norm << "): "
364 << "λ = " << lambda << std::endl;
365 if (verbosity_level_ > 2)
366 {
367 // print approximated eigenvector via DUNE-ISTL I/O methods
368 Dune::printvector(std::cout,x,blank_+"x",blank_+"row");
369 }
370 }
371
372 // free dynamically allocated memory
373 delete[] Ax_raw;
374 delete[] ev;
375 }
376
388 inline void computeSymMinMagnitude (const Real& epsilon,
389 BlockVector& x, Real& lambda) const
390 {
391 // print verbosity information
392 if (verbosity_level_ > 0)
393 std::cout << title_ << "Computing an approximation of the "
394 << "least dominant eigenvalue of a matrix which "
395 << "is assumed to be symmetric." << std::endl;
396
397 // use type ArPackPlusPlus_BCRSMatrixWrapper to store matrix information
398 // and to perform the product A*v (LU decomposition is not used)
399 typedef Impl::ArPackPlusPlus_BCRSMatrixWrapper<BCRSMatrix> WrappedMatrix;
400 WrappedMatrix A(m_);
401
402 // get number of rows and columns in A
403 const int nrows = A.nrows();
404 const int ncols = A.ncols();
405
406 // assert that A is square
407 if (nrows != ncols)
408 DUNE_THROW(Dune::ISTLError,"Matrix is not square ("
409 << nrows << "x" << ncols << ").");
410
411 // allocate memory for variables, set parameters
412 const int nev = 1; // Number of eigenvalues to compute
413 int ncv = std::min(20, nrows); // Number of Arnoldi vectors generated at each iteration (0 == auto)
414 const Real tol = epsilon; // Stopping tolerance (relative accuracy of Ritz values) (0 == machine precision)
415 const int maxit = nIterationsMax_*nev; // Maximum number of Arnoldi update iterations allowed (0 == 100*nev)
416 Real* ev = new Real[nev]; // Computed eigenvalues of A
417 const bool ivec = true; // Flag deciding if eigenvectors shall be determined
418 int nconv; // Number of converged eigenvalues
419
420 // define what we need: eigenvalues with smallest magnitude
421 char which[] = "SM";
422 ARSymStdEig<Real,WrappedMatrix>
423 dprob(nrows, nev, &A, &WrappedMatrix::multMv, which, ncv, tol, maxit);
424
425 // set ARPACK verbosity mode if requested
426 if (verbosity_level_ > 3) dprob.Trace();
427
428 // find eigenvalues and eigenvectors of A, obtain the eigenvalues
429 nconv = dprob.Eigenvalues(ev,ivec);
430
431 // obtain approximated least dominant eigenvalue of A
432 lambda = ev[nev-1];
433
434 // obtain associated approximated eigenvector of A
435 Real* x_raw = dprob.RawEigenvector(nev-1);
436 WrappedMatrix::arrayToDomainBlockVector(x_raw,x);
437
438 // obtain number of Arnoldi update iterations actually taken
439 nIterations_ = dprob.GetIter();
440
441 // compute residual norm
442 BlockVector r(x);
443 Real* Ax_raw = new Real[nrows];
444 A.multMv(x_raw,Ax_raw);
445 WrappedMatrix::arrayToDomainBlockVector(Ax_raw,r);
446 r.axpy(-lambda,x);
447 const Real r_norm = r.two_norm();
448
449 // print verbosity information
450 if (verbosity_level_ > 0)
451 {
452 if (verbosity_level_ > 1)
453 {
454 // print some information about the problem
455 std::cout << blank_ << "Obtained eigenvalues of A by solving "
456 << "A*x = λ*x using the ARPACK++ class ARSym"
457 << "StdEig:" << std::endl;
458 std::cout << blank_ << " converged eigenvalues of A: "
459 << nconv << " / " << nev << std::endl;
460 std::cout << blank_ << " least dominant eigenvalue of A: "
461 << lambda << std::endl;
462 }
463 std::cout << blank_ << "Result ("
464 << "#iterations = " << nIterations_ << ", "
465 << "║residual║_2 = " << r_norm << "): "
466 << "λ = " << lambda << std::endl;
467 if (verbosity_level_ > 2)
468 {
469 // print approximated eigenvector via DUNE-ISTL I/O methods
470 Dune::printvector(std::cout,x,blank_+"x",blank_+"row");
471 }
472 }
473
474 // free dynamically allocated memory
475 delete[] Ax_raw;
476 delete[] ev;
477 }
478
490 inline void computeSymCond2 (const Real& epsilon, Real& cond_2) const
491 {
492 // print verbosity information
493 if (verbosity_level_ > 0)
494 std::cout << title_ << "Computing an approximation of the "
495 << "spectral condition number of a matrix which "
496 << "is assumed to be symmetric." << std::endl;
497
498 // use type ArPackPlusPlus_BCRSMatrixWrapper to store matrix information
499 // and to perform the product A*v (LU decomposition is not used)
500 typedef Impl::ArPackPlusPlus_BCRSMatrixWrapper<BCRSMatrix> WrappedMatrix;
501 WrappedMatrix A(m_);
502
503 // get number of rows and columns in A
504 const int nrows = A.nrows();
505 const int ncols = A.ncols();
506
507 // assert that A is square
508 if (nrows != ncols)
509 DUNE_THROW(Dune::ISTLError,"Matrix is not square ("
510 << nrows << "x" << ncols << ").");
511
512 // allocate memory for variables, set parameters
513 const int nev = 2; // Number of eigenvalues to compute
514 int ncv = std::min(20, nrows); // Number of Arnoldi vectors generated at each iteration (0 == auto)
515 const Real tol = epsilon; // Stopping tolerance (relative accuracy of Ritz values) (0 == machine precision)
516 const int maxit = nIterationsMax_*nev; // Maximum number of Arnoldi update iterations allowed (0 == 100*nev)
517 Real* ev = new Real[nev]; // Computed eigenvalues of A
518 const bool ivec = true; // Flag deciding if eigenvectors shall be determined
519 int nconv; // Number of converged eigenvalues
520
521 // define what we need: eigenvalues from both ends of the spectrum
522 char which[] = "BE";
523 ARSymStdEig<Real,WrappedMatrix>
524 dprob(nrows, nev, &A, &WrappedMatrix::multMv, which, ncv, tol, maxit);
525
526 // set ARPACK verbosity mode if requested
527 if (verbosity_level_ > 3) dprob.Trace();
528
529 // find eigenvalues and eigenvectors of A, obtain the eigenvalues
530 nconv = dprob.Eigenvalues(ev,ivec);
531
532 // obtain approximated dominant and least dominant eigenvalue of A
533 const Real& lambda_max = ev[nev-1];
534 const Real& lambda_min = ev[0];
535
536 // obtain associated approximated eigenvectors of A
537 Real* x_max_raw = dprob.RawEigenvector(nev-1);
538 Real* x_min_raw = dprob.RawEigenvector(0);
539
540 // obtain approximated spectral condition number of A
541 cond_2 = std::abs(lambda_max / lambda_min);
542
543 // obtain number of Arnoldi update iterations actually taken
544 nIterations_ = dprob.GetIter();
545
546 // compute each residual norm
547 Real* Ax_max_raw = new Real[nrows];
548 Real* Ax_min_raw = new Real[nrows];
549 A.multMv(x_max_raw,Ax_max_raw);
550 A.multMv(x_min_raw,Ax_min_raw);
551 Real r_max_norm = 0.0;
552 Real r_min_norm = 0.0;
553 for (int i = 0; i < nrows; ++i)
554 {
555 r_max_norm += std::pow(Ax_max_raw[i] - lambda_max * x_max_raw[i],2);
556 r_min_norm += std::pow(Ax_min_raw[i] - lambda_min * x_min_raw[i],2);
557 }
558 r_max_norm = std::sqrt(r_max_norm);
559 r_min_norm = std::sqrt(r_min_norm);
560
561 // print verbosity information
562 if (verbosity_level_ > 0)
563 {
564 if (verbosity_level_ > 1)
565 {
566 // print some information about the problem
567 std::cout << blank_ << "Obtained eigenvalues of A by solving "
568 << "A*x = λ*x using the ARPACK++ class ARSym"
569 << "StdEig:" << std::endl;
570 std::cout << blank_ << " converged eigenvalues of A: "
571 << nconv << " / " << nev << std::endl;
572 std::cout << blank_ << " dominant eigenvalue of A: "
573 << lambda_max << std::endl;
574 std::cout << blank_ << " least dominant eigenvalue of A: "
575 << lambda_min << std::endl;
576 std::cout << blank_ << " spectral condition number of A: "
577 << cond_2 << std::endl;
578 }
579 std::cout << blank_ << "Result ("
580 << "#iterations = " << nIterations_ << ", "
581 << "║residual║_2 = {" << r_max_norm << ","
582 << r_min_norm << "}, " << "λ = {"
583 << lambda_max << "," << lambda_min
584 << "}): cond_2 = " << cond_2 << std::endl;
585 }
586
587 // free dynamically allocated memory
588 delete[] Ax_min_raw;
589 delete[] Ax_max_raw;
590 delete[] ev;
591 }
592
606 inline void computeNonSymMax (const Real& epsilon,
607 BlockVector& x, Real& sigma) const
608 {
609 // print verbosity information
610 if (verbosity_level_ > 0)
611 std::cout << title_ << "Computing an approximation of the "
612 << "largest singular value of a matrix which "
613 << "is assumed to be nonsymmetric." << std::endl;
614
615 // use type ArPackPlusPlus_BCRSMatrixWrapper to store matrix information
616 // and to perform the product A^T*A*v (LU decomposition is not used)
617 typedef Impl::ArPackPlusPlus_BCRSMatrixWrapper<BCRSMatrix> WrappedMatrix;
618 WrappedMatrix A(m_);
619
620 // get number of rows and columns in A
621 const int nrows = A.nrows();
622 const int ncols = A.ncols();
623
624 // assert that A has more rows than columns (extend code later to the opposite case!)
625 if (nrows < ncols)
626 DUNE_THROW(Dune::ISTLError,"Matrix has less rows than "
627 << "columns (" << nrows << "x" << ncols << ")."
628 << " This case is not implemented, yet.");
629
630 // allocate memory for variables, set parameters
631 const int nev = 1; // Number of eigenvalues to compute
632 int ncv = std::min(20, nrows); // Number of Arnoldi vectors generated at each iteration (0 == auto)
633 const Real tol = epsilon; // Stopping tolerance (relative accuracy of Ritz values) (0 == machine precision)
634 const int maxit = nIterationsMax_*nev; // Maximum number of Arnoldi update iterations allowed (0 == 100*nev)
635 Real* ev = new Real[nev]; // Computed eigenvalues of A^T*A
636 const bool ivec = true; // Flag deciding if eigenvectors shall be determined
637 int nconv; // Number of converged eigenvalues
638
639 // define what we need: eigenvalues with largest algebraic value
640 char which[] = "LA";
641 ARSymStdEig<Real,WrappedMatrix>
642 dprob(ncols, nev, &A, &WrappedMatrix::multMtMv, which, ncv, tol, maxit);
643
644 // set ARPACK verbosity mode if requested
645 if (verbosity_level_ > 3) dprob.Trace();
646
647 // find eigenvalues and eigenvectors of A^T*A, obtain the eigenvalues
648 nconv = dprob.Eigenvalues(ev,ivec);
649
650 // obtain approximated largest eigenvalue of A^T*A
651 const Real& lambda = ev[nev-1];
652
653 // obtain associated approximated eigenvector of A^T*A
654 Real* x_raw = dprob.RawEigenvector(nev-1);
655 WrappedMatrix::arrayToDomainBlockVector(x_raw,x);
656
657 // obtain number of Arnoldi update iterations actually taken
658 nIterations_ = dprob.GetIter();
659
660 // compute residual norm
661 BlockVector r(x);
662 Real* AtAx_raw = new Real[ncols];
663 A.multMtMv(x_raw,AtAx_raw);
664 WrappedMatrix::arrayToDomainBlockVector(AtAx_raw,r);
665 r.axpy(-lambda,x);
666 const Real r_norm = r.two_norm();
667
668 // calculate largest singular value of A (note that
669 // x is right-singular / left-singular vector of A)
670 sigma = std::sqrt(lambda);
671
672 // print verbosity information
673 if (verbosity_level_ > 0)
674 {
675 if (verbosity_level_ > 1)
676 {
677 // print some information about the problem
678 std::cout << blank_ << "Obtained singular values of A by sol"
679 << "ving (A^T*A)*x = σ²*x using the ARPACK++ "
680 << "class ARSymStdEig:" << std::endl;
681 std::cout << blank_ << " converged eigenvalues of A^T*A: "
682 << nconv << " / " << nev << std::endl;
683 std::cout << blank_ << " largest eigenvalue of A^T*A: "
684 << lambda << std::endl;
685 std::cout << blank_ << " => largest singular value of A: "
686 << sigma << std::endl;
687 }
688 std::cout << blank_ << "Result ("
689 << "#iterations = " << nIterations_ << ", "
690 << "║residual║_2 = " << r_norm << "): "
691 << "σ = " << sigma << std::endl;
692 if (verbosity_level_ > 2)
693 {
694 // print approximated right-singular / left-singular vector
695 // via DUNE-ISTL I/O methods
696 Dune::printvector(std::cout,x,blank_+"x",blank_+"row");
697 }
698 }
699
700 // free dynamically allocated memory
701 delete[] AtAx_raw;
702 delete[] ev;
703 }
704
718 inline void computeNonSymMin (const Real& epsilon,
719 BlockVector& x, Real& sigma) const
720 {
721 // print verbosity information
722 if (verbosity_level_ > 0)
723 std::cout << title_ << "Computing an approximation of the "
724 << "smallest singular value of a matrix which "
725 << "is assumed to be nonsymmetric." << std::endl;
726
727 // use type ArPackPlusPlus_BCRSMatrixWrapper to store matrix information
728 // and to perform the product A^T*A*v (LU decomposition is not used)
729 typedef Impl::ArPackPlusPlus_BCRSMatrixWrapper<BCRSMatrix> WrappedMatrix;
730 WrappedMatrix A(m_);
731
732 // get number of rows and columns in A
733 const int nrows = A.nrows();
734 const int ncols = A.ncols();
735
736 // assert that A has more rows than columns (extend code later to the opposite case!)
737 if (nrows < ncols)
738 DUNE_THROW(Dune::ISTLError,"Matrix has less rows than "
739 << "columns (" << nrows << "x" << ncols << ")."
740 << " This case is not implemented, yet.");
741
742 // allocate memory for variables, set parameters
743 const int nev = 1; // Number of eigenvalues to compute
744 int ncv = std::min(20, nrows); // Number of Arnoldi vectors generated at each iteration (0 == auto)
745 const Real tol = epsilon; // Stopping tolerance (relative accuracy of Ritz values) (0 == machine precision)
746 const int maxit = nIterationsMax_*nev; // Maximum number of Arnoldi update iterations allowed (0 == 100*nev)
747 Real* ev = new Real[nev]; // Computed eigenvalues of A^T*A
748 const bool ivec = true; // Flag deciding if eigenvectors shall be determined
749 int nconv; // Number of converged eigenvalues
750
751 // define what we need: eigenvalues with smallest algebraic value
752 char which[] = "SA";
753 ARSymStdEig<Real,WrappedMatrix>
754 dprob(ncols, nev, &A, &WrappedMatrix::multMtMv, which, ncv, tol, maxit);
755
756 // set ARPACK verbosity mode if requested
757 if (verbosity_level_ > 3) dprob.Trace();
758
759 // find eigenvalues and eigenvectors of A^T*A, obtain the eigenvalues
760 nconv = dprob.Eigenvalues(ev,ivec);
761
762 // obtain approximated smallest eigenvalue of A^T*A
763 const Real& lambda = ev[nev-1];
764
765 // obtain associated approximated eigenvector of A^T*A
766 Real* x_raw = dprob.RawEigenvector(nev-1);
767 WrappedMatrix::arrayToDomainBlockVector(x_raw,x);
768
769 // obtain number of Arnoldi update iterations actually taken
770 nIterations_ = dprob.GetIter();
771
772 // compute residual norm
773 BlockVector r(x);
774 Real* AtAx_raw = new Real[ncols];
775 A.multMtMv(x_raw,AtAx_raw);
776 WrappedMatrix::arrayToDomainBlockVector(AtAx_raw,r);
777 r.axpy(-lambda,x);
778 const Real r_norm = r.two_norm();
779
780 // calculate smallest singular value of A (note that
781 // x is right-singular / left-singular vector of A)
782 sigma = std::sqrt(lambda);
783
784 // print verbosity information
785 if (verbosity_level_ > 0)
786 {
787 if (verbosity_level_ > 1)
788 {
789 // print some information about the problem
790 std::cout << blank_ << "Obtained singular values of A by sol"
791 << "ving (A^T*A)*x = σ²*x using the ARPACK++ "
792 << "class ARSymStdEig:" << std::endl;
793 std::cout << blank_ << " converged eigenvalues of A^T*A: "
794 << nconv << " / " << nev << std::endl;
795 std::cout << blank_ << " smallest eigenvalue of A^T*A: "
796 << lambda << std::endl;
797 std::cout << blank_ << " => smallest singular value of A: "
798 << sigma << std::endl;
799 }
800 std::cout << blank_ << "Result ("
801 << "#iterations = " << nIterations_ << ", "
802 << "║residual║_2 = " << r_norm << "): "
803 << "σ = " << sigma << std::endl;
804 if (verbosity_level_ > 2)
805 {
806 // print approximated right-singular / left-singular vector
807 // via DUNE-ISTL I/O methods
808 Dune::printvector(std::cout,x,blank_+"x",blank_+"row");
809 }
810 }
811
812 // free dynamically allocated memory
813 delete[] AtAx_raw;
814 delete[] ev;
815 }
816
827 inline void computeNonSymCond2 (const Real& epsilon, Real& cond_2) const
828 {
829 // print verbosity information
830 if (verbosity_level_ > 0)
831 std::cout << title_ << "Computing an approximation of the "
832 << "spectral condition number of a matrix which "
833 << "is assumed to be nonsymmetric." << std::endl;
834
835 // use type ArPackPlusPlus_BCRSMatrixWrapper to store matrix information
836 // and to perform the product A^T*A*v (LU decomposition is not used)
837 typedef Impl::ArPackPlusPlus_BCRSMatrixWrapper<BCRSMatrix> WrappedMatrix;
838 WrappedMatrix A(m_);
839
840 // get number of rows and columns in A
841 const int nrows = A.nrows();
842 const int ncols = A.ncols();
843
844 // assert that A has more rows than columns (extend code later to the opposite case!)
845 if (nrows < ncols)
846 DUNE_THROW(Dune::ISTLError,"Matrix has less rows than "
847 << "columns (" << nrows << "x" << ncols << ")."
848 << " This case is not implemented, yet.");
849
850 // allocate memory for variables, set parameters
851 const int nev = 2; // Number of eigenvalues to compute
852 int ncv = std::min(20, nrows); // Number of Arnoldi vectors generated at each iteration (0 == auto)
853 const Real tol = epsilon; // Stopping tolerance (relative accuracy of Ritz values) (0 == machine precision)
854 const int maxit = nIterationsMax_*nev; // Maximum number of Arnoldi update iterations allowed (0 == 100*nev)
855 Real* ev = new Real[nev]; // Computed eigenvalues of A^T*A
856 const bool ivec = true; // Flag deciding if eigenvectors shall be determined
857 int nconv; // Number of converged eigenvalues
858
859 // define what we need: eigenvalues from both ends of the spectrum
860 char which[] = "BE";
861 ARSymStdEig<Real,WrappedMatrix>
862 dprob(ncols, nev, &A, &WrappedMatrix::multMtMv, which, ncv, tol, maxit);
863
864 // set ARPACK verbosity mode if requested
865 if (verbosity_level_ > 3) dprob.Trace();
866
867 // find eigenvalues and eigenvectors of A^T*A, obtain the eigenvalues
868 nconv = dprob.Eigenvalues(ev,ivec);
869
870 // obtain approximated largest and smallest eigenvalue of A^T*A
871 const Real& lambda_max = ev[nev-1];
872 const Real& lambda_min = ev[0];
873
874 // obtain associated approximated eigenvectors of A^T*A
875 Real* x_max_raw = dprob.RawEigenvector(nev-1);
876 Real* x_min_raw = dprob.RawEigenvector(0);
877
878 // obtain number of Arnoldi update iterations actually taken
879 nIterations_ = dprob.GetIter();
880
881 // compute each residual norm
882 Real* AtAx_max_raw = new Real[ncols];
883 Real* AtAx_min_raw = new Real[ncols];
884 A.multMtMv(x_max_raw,AtAx_max_raw);
885 A.multMtMv(x_min_raw,AtAx_min_raw);
886 Real r_max_norm = 0.0;
887 Real r_min_norm = 0.0;
888 for (int i = 0; i < ncols; ++i)
889 {
890 r_max_norm += std::pow(AtAx_max_raw[i] - lambda_max * x_max_raw[i],2);
891 r_min_norm += std::pow(AtAx_min_raw[i] - lambda_min * x_min_raw[i],2);
892 }
893 r_max_norm = std::sqrt(r_max_norm);
894 r_min_norm = std::sqrt(r_min_norm);
895
896 // calculate largest and smallest singular value of A
897 const Real sigma_max = std::sqrt(lambda_max);
898 const Real sigma_min = std::sqrt(lambda_min);
899
900 // obtain approximated spectral condition number of A
901 cond_2 = sigma_max / sigma_min;
902
903 // print verbosity information
904 if (verbosity_level_ > 0)
905 {
906 if (verbosity_level_ > 1)
907 {
908 // print some information about the problem
909 std::cout << blank_ << "Obtained singular values of A by sol"
910 << "ving (A^T*A)*x = σ²*x using the ARPACK++ "
911 << "class ARSymStdEig:" << std::endl;
912 std::cout << blank_ << " converged eigenvalues of A^T*A: "
913 << nconv << " / " << nev << std::endl;
914 std::cout << blank_ << " largest eigenvalue of A^T*A: "
915 << lambda_max << std::endl;
916 std::cout << blank_ << " smallest eigenvalue of A^T*A: "
917 << lambda_min << std::endl;
918 std::cout << blank_ << " => largest singular value of A: "
919 << sigma_max << std::endl;
920 std::cout << blank_ << " => smallest singular value of A: "
921 << sigma_min << std::endl;
922 }
923 std::cout << blank_ << "Result ("
924 << "#iterations = " << nIterations_ << ", "
925 << "║residual║_2 = {" << r_max_norm << ","
926 << r_min_norm << "}, " << "σ = {"
927 << sigma_max << "," << sigma_min
928 << "}): cond_2 = " << cond_2 << std::endl;
929 }
930
931 // free dynamically allocated memory
932 delete[] AtAx_min_raw;
933 delete[] AtAx_max_raw;
934 delete[] ev;
935 }
936
941 inline unsigned int getIterationCount () const
942 {
943 if (nIterations_ == 0)
944 DUNE_THROW(Dune::ISTLError,"No algorithm applied, yet.");
945
946 return nIterations_;
947 }
948
949 protected:
950 // parameters related to iterative eigenvalue algorithms
951 const BCRSMatrix& m_;
952 const unsigned int nIterationsMax_;
953
954 // verbosity setting
955 const unsigned int verbosity_level_;
956
957 // memory for storing temporary variables (mutable as they shall
958 // just be effectless auxiliary variables of the const apply*(...)
959 // methods)
960 mutable unsigned int nIterations_;
961
962 // constants for printing verbosity information
963 const std::string title_;
964 const std::string blank_;
965 };
966
969} // namespace Dune
970
971#endif // HAVE_ARPACKPP
972
973#endif // DUNE_ISTL_EIGENVALUE_ARPACKPP_HH
This file implements a vector space as a tensor product of a given vector space. The number of compon...
Wrapper to use a range of ARPACK++ eigenvalue solvers.
Definition: arpackpp.hh:242
unsigned int getIterationCount() const
Return the number of iterations in last application of an algorithm.
Definition: arpackpp.hh:941
void computeSymMaxMagnitude(const Real &epsilon, BlockVector &x, Real &lambda) const
Assume the matrix to be square, symmetric and perform IRLM to compute an approximation lambda of its ...
Definition: arpackpp.hh:286
void computeSymMinMagnitude(const Real &epsilon, BlockVector &x, Real &lambda) const
Assume the matrix to be square, symmetric and perform IRLM to compute an approximation lambda of its ...
Definition: arpackpp.hh:388
ArPackPlusPlus_Algorithms(const BCRSMatrix &m, const unsigned int nIterationsMax=100000, const unsigned int verbosity_level=0)
Construct from required parameters.
Definition: arpackpp.hh:265
void computeSymCond2(const Real &epsilon, Real &cond_2) const
Assume the matrix to be square, symmetric and perform IRLM to compute an approximation of its spectra...
Definition: arpackpp.hh:490
void computeNonSymMin(const Real &epsilon, BlockVector &x, Real &sigma) const
Assume the matrix to be nonsymmetric and perform IRLM to compute an approximation sigma of its smalle...
Definition: arpackpp.hh:718
void computeNonSymCond2(const Real &epsilon, Real &cond_2) const
Assume the matrix to be nonsymmetric and perform IRLM to compute an approximation of its spectral con...
Definition: arpackpp.hh:827
void computeNonSymMax(const Real &epsilon, BlockVector &x, Real &sigma) const
Assume the matrix to be nonsymmetric and perform IRLM to compute an approximation sigma of its larges...
Definition: arpackpp.hh:606
A sparse block matrix with compressed row storage.
Definition: bcrsmatrix.hh:423
B::field_type field_type
export the type representing the field
Definition: bcrsmatrix.hh:445
@ blocklevel
The number of blocklevels the matrix contains.
Definition: bcrsmatrix.hh:465
A vector of blocks with memory management.
Definition: bvector.hh:317
A::size_type size_type
The type for the index access.
Definition: bvector.hh:332
B::field_type field_type
export the type representing the field
Definition: bvector.hh:323
vector space out of a tensor product of fields.
Definition: fvector.hh:93
derive error class from the base class in common
Definition: istlexception.hh:16
A few common exception classes.
Implements a vector constructed from a given type representing a field and a compile-time given size.
#define DUNE_THROW(E, m)
Definition: exceptions.hh:216
void printvector(std::ostream &s, const V &v, std::string title, std::string rowtext, int columns=1, int width=10, int precision=2)
Print an ISTL vector.
Definition: io.hh:102
Some generic functions for pretty printing vectors and matrices.
Dune namespace.
Definition: alignedallocator.hh:10
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 13, 23:29, 2024)