dune-composites (2.5.1)

arpackpp_fork.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_FORK
4#define DUNE_ISTL_EIGENVALUE_ARPACKPP_HH_FORK
5
6#if HAVE_ARPACKPP
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#include "argsym.h" // provides ARSymGenEig
29#include "argnsym.h" // provides ARSymGenEig
30
31
32namespace ArpackFork
33{
34
50 template <class BCRSMatrix>
51 class ArPackPlusPlus_BCRSMatrixWrapper
52 {
53 public:
55 typedef typename BCRSMatrix::field_type Real;
56
57 public:
59 ArPackPlusPlus_BCRSMatrixWrapper (const BCRSMatrix& A)
60 : A_(A),
61 a_(A_.M() * mBlock), n_(A_.N() * nBlock)
62 {
63 // assert that BCRSMatrix type has blocklevel 2
64 static_assert
65 (BCRSMatrix::blocklevel == 2,
66 "Only BCRSMatrices with blocklevel 2 are supported.");
67
68 // allocate memory for auxiliary block vector objects
69 // which are compatible to matrix rows / columns
70 domainBlockVector.resize(A_.N(),false);
71 rangeBlockVector.resize(A_.M(),false);
72 }
73
75 inline void multMv (Real* v, Real* w)
76 {
77 // get vector v as an object of appropriate type
78 arrayToDomainBlockVector(v,domainBlockVector);
79
80 // perform matrix-vector product
81 A_.mv(domainBlockVector,rangeBlockVector);
82
83 // get vector w from object of appropriate type
84 rangeBlockVectorToArray(rangeBlockVector,w);
85 };
86
88 inline void multMtMv (Real* v, Real* w)
89 {
90 // get vector v as an object of appropriate type
91 arrayToDomainBlockVector(v,domainBlockVector);
92
93 // perform matrix-vector product
94 A_.mv(domainBlockVector,rangeBlockVector);
95 A_.mtv(rangeBlockVector,domainBlockVector);
96
97 // get vector w from object of appropriate type
98 domainBlockVectorToArray(domainBlockVector,w);
99 };
100
102 inline void multMMtv (Real* v, Real* w)
103 {
104 // get vector v as an object of appropriate type
105 arrayToRangeBlockVector(v,rangeBlockVector);
106
107 // perform matrix-vector product
108 A_.mtv(rangeBlockVector,domainBlockVector);
109 A_.mv(domainBlockVector,rangeBlockVector);
110
111 // get vector w from object of appropriate type
112 rangeBlockVectorToArray(rangeBlockVector,w);
113 };
114
116 inline int nrows () const { return a_; }
117
119 inline int ncols () const { return n_; }
120
121 protected:
122 // Number of rows and columns in each block of the matrix
123 constexpr static int mBlock = BCRSMatrix::block_type::rows;
124 constexpr static int nBlock = BCRSMatrix::block_type::cols;
125
126 // Type of vectors in the domain of the linear map associated with
127 // the matrix, i.e. block vectors compatible to matrix rows
128 constexpr static int dbvBlockSize = nBlock;
129 typedef Dune::FieldVector<Real,dbvBlockSize> DomainBlockVectorBlock;
130 typedef Dune::BlockVector<DomainBlockVectorBlock> DomainBlockVector;
131
132 // Type of vectors in the range of the linear map associated with
133 // the matrix, i.e. block vectors compatible to matrix columns
134 constexpr static int rbvBlockSize = mBlock;
135 typedef Dune::FieldVector<Real,rbvBlockSize> RangeBlockVectorBlock;
136 typedef Dune::BlockVector<RangeBlockVectorBlock> RangeBlockVector;
137
138 // Types for vector index access
139 typedef typename DomainBlockVector::size_type dbv_size_type;
140 typedef typename RangeBlockVector::size_type rbv_size_type;
141 typedef typename DomainBlockVectorBlock::size_type dbvb_size_type;
142 typedef typename RangeBlockVectorBlock::size_type rbvb_size_type;
143
144 // Get vector v from a block vector object which is compatible to
145 // matrix rows
146 static inline void
147 domainBlockVectorToArray (const DomainBlockVector& dbv, Real* v)
148 {
149 for (dbv_size_type block = 0; block < dbv.N(); ++block)
150 for (dbvb_size_type iBlock = 0; iBlock < dbvBlockSize; ++iBlock)
151 v[block*dbvBlockSize + iBlock] = dbv[block][iBlock];
152 }
153
154 // Get vector v from a block vector object which is compatible to
155 // matrix columns
156 static inline void
157 rangeBlockVectorToArray (const RangeBlockVector& rbv, Real* v)
158 {
159 for (rbv_size_type block = 0; block < rbv.N(); ++block)
160 for (rbvb_size_type iBlock = 0; iBlock < rbvBlockSize; ++iBlock)
161 v[block*rbvBlockSize + iBlock] = rbv[block][iBlock];
162 }
163
164 public:
167 static inline void arrayToDomainBlockVector (const Real* v,
168 DomainBlockVector& dbv)
169 {
170 for (dbv_size_type block = 0; block < dbv.N(); ++block)
171 for (dbvb_size_type iBlock = 0; iBlock < dbvBlockSize; ++iBlock)
172 dbv[block][iBlock] = v[block*dbvBlockSize + iBlock];
173 }
174
177 static inline void arrayToRangeBlockVector (const Real* v,
178 RangeBlockVector& rbv)
179 {
180 for (rbv_size_type block = 0; block < rbv.N(); ++block)
181 for (rbvb_size_type iBlock = 0; iBlock < rbvBlockSize; ++iBlock)
182 rbv[block][iBlock] = v[block*rbvBlockSize + iBlock];
183 }
184
185 protected:
186 // The DUNE-ISTL BCRSMatrix
187 const BCRSMatrix& A_;
188
189 // Number of rows and columns in the matrix
190 const int a_, n_;
191
192 // Auxiliary block vector objects which are
193 // compatible to matrix rows / columns
194 mutable DomainBlockVector domainBlockVector;
195 mutable RangeBlockVector rangeBlockVector;
196 };
197
198
199
200 template <class BCRSMatrix>
201 class ArPackPlusPlus_BCRSMatrixWrapper2
202 {
203 public:
205 typedef typename BCRSMatrix::field_type Real;
206
207 public:
209 ArPackPlusPlus_BCRSMatrixWrapper2 (const BCRSMatrix& A, const BCRSMatrix& B)
210 : A_(A), B_(B),
211 A_solver(A_),
212 a_(A_.M() * mBlock), n_(A_.N() * nBlock)
213 {
214 // assert that BCRSMatrix type has blocklevel 2
215 static_assert
216 (BCRSMatrix::blocklevel == 2,
217 "Only BCRSMatrices with blocklevel 2 are supported.");
218
219 // allocate memory for auxiliary block vector objects
220 // which are compatible to matrix rows / columns
221 domainBlockVector.resize(A_.N(),false);
222 rangeBlockVector.resize(A_.M(),false);
223 }
224
226 inline void multMv (Real* v, Real* w)
227 {
228 // get vector v as an object of appropriate type
229 arrayToDomainBlockVector(v,domainBlockVector);
230
231 // perform matrix-vector product
232 B_.mv(domainBlockVector,rangeBlockVector);
233
234 Dune::InverseOperatorResult result;
235 auto rangeBlockVector2 = rangeBlockVector;
236 A_solver.apply(rangeBlockVector2, rangeBlockVector, result);
237
238 // get vector w from object of appropriate type
239 rangeBlockVectorToArray(rangeBlockVector2,w);
240 };
241
242 inline void multMvB (Real* v, Real* w)
243 {
244 // get vector v as an object of appropriate type
245 arrayToDomainBlockVector(v,domainBlockVector);
246
247 // perform matrix-vector product
248 B_.mv(domainBlockVector,rangeBlockVector);
249
250 // get vector w from object of appropriate type
251 rangeBlockVectorToArray(rangeBlockVector,w);
252 };
253
254
256 /*inline void multMtMv (Real* v, Real* w)
257 {
258 // get vector v as an object of appropriate type
259 arrayToDomainBlockVector(v,domainBlockVector);
260
261 // perform matrix-vector product
262 A_.mv(domainBlockVector,rangeBlockVector);
263 A_.mtv(rangeBlockVector,domainBlockVector);
264
265 // get vector w from object of appropriate type
266 domainBlockVectorToArray(domainBlockVector,w);
267 };
268
270 inline void multMMtv (Real* v, Real* w)
271 {
272 // get vector v as an object of appropriate type
273 arrayToRangeBlockVector(v,rangeBlockVector);
274
275 // perform matrix-vector product
276 A_.mtv(rangeBlockVector,domainBlockVector);
277 A_.mv(domainBlockVector,rangeBlockVector);
278
279 // get vector w from object of appropriate type
280 rangeBlockVectorToArray(rangeBlockVector,w);
281 };*/
282
284 inline int nrows () const { return a_; }
285
287 inline int ncols () const { return n_; }
288
289 protected:
290 // Number of rows and columns in each block of the matrix
291 constexpr static int mBlock = BCRSMatrix::block_type::rows;
292 constexpr static int nBlock = BCRSMatrix::block_type::cols;
293
294 // Type of vectors in the domain of the linear map associated with
295 // the matrix, i.e. block vectors compatible to matrix rows
296 constexpr static int dbvBlockSize = nBlock;
297 typedef Dune::FieldVector<Real,dbvBlockSize> DomainBlockVectorBlock;
298 typedef Dune::BlockVector<DomainBlockVectorBlock> DomainBlockVector;
299
300 // Type of vectors in the range of the linear map associated with
301 // the matrix, i.e. block vectors compatible to matrix columns
302 constexpr static int rbvBlockSize = mBlock;
303 typedef Dune::FieldVector<Real,rbvBlockSize> RangeBlockVectorBlock;
304 typedef Dune::BlockVector<RangeBlockVectorBlock> RangeBlockVector;
305
306 // Types for vector index access
307 typedef typename DomainBlockVector::size_type dbv_size_type;
308 typedef typename RangeBlockVector::size_type rbv_size_type;
309 typedef typename DomainBlockVectorBlock::size_type dbvb_size_type;
310 typedef typename RangeBlockVectorBlock::size_type rbvb_size_type;
311
312 // Get vector v from a block vector object which is compatible to
313 // matrix rows
314 static inline void
315 domainBlockVectorToArray (const DomainBlockVector& dbv, Real* v)
316 {
317 for (dbv_size_type block = 0; block < dbv.N(); ++block)
318 for (dbvb_size_type iBlock = 0; iBlock < dbvBlockSize; ++iBlock)
319 v[block*dbvBlockSize + iBlock] = dbv[block][iBlock];
320 }
321
322 // Get vector v from a block vector object which is compatible to
323 // matrix columns
324 static inline void
325 rangeBlockVectorToArray (const RangeBlockVector& rbv, Real* v)
326 {
327 for (rbv_size_type block = 0; block < rbv.N(); ++block)
328 for (rbvb_size_type iBlock = 0; iBlock < rbvBlockSize; ++iBlock)
329 v[block*rbvBlockSize + iBlock] = rbv[block][iBlock];
330 }
331
332 public:
335 static inline void arrayToDomainBlockVector (const Real* v,
336 DomainBlockVector& dbv)
337 {
338 for (dbv_size_type block = 0; block < dbv.N(); ++block)
339 for (dbvb_size_type iBlock = 0; iBlock < dbvBlockSize; ++iBlock)
340 dbv[block][iBlock] = v[block*dbvBlockSize + iBlock];
341 }
342
345 static inline void arrayToRangeBlockVector (const Real* v,
346 RangeBlockVector& rbv)
347 {
348 for (rbv_size_type block = 0; block < rbv.N(); ++block)
349 for (rbvb_size_type iBlock = 0; iBlock < rbvBlockSize; ++iBlock)
350 rbv[block][iBlock] = v[block*rbvBlockSize + iBlock];
351 }
352
353 protected:
354 // The DUNE-ISTL BCRSMatrix
355 const BCRSMatrix& A_;
356 const BCRSMatrix& B_;
357
358 Dune::UMFPack<BCRSMatrix> A_solver;
359
360 // Number of rows and columns in the matrix
361 const int a_, n_;
362
363 // Auxiliary block vector objects which are
364 // compatible to matrix rows / columns
365 mutable DomainBlockVector domainBlockVector;
366 mutable RangeBlockVector rangeBlockVector;
367 };
368
369
407 template <typename BCRSMatrix, typename BlockVector>
408 class ArPackPlusPlus_Algorithms
409 {
410 public:
411 typedef typename BlockVector::field_type Real;
412
413 public:
432 ArPackPlusPlus_Algorithms (const BCRSMatrix& m,
433 const unsigned int nIterationsMax = 100000,
434 const unsigned int verbosity_level = 0)
435 : a_(m), nIterationsMax_(nIterationsMax),
436 verbosity_level_(verbosity_level),
437 nIterations_(0),
438 title_(" ArPackPlusPlus_Algorithms: "),
439 blank_(title_.length(),' ')
440 {}
441
453 inline void computeSymMaxMagnitude (const Real& epsilon,
454 BlockVector& x, Real& lambda) const
455 {
456 // print verbosity information
457 if (verbosity_level_ > 0)
458 std::cout << title_ << "Computing an approximation of "
459 << "the dominant eigenvalue of a matrix which "
460 << "is assumed to be symmetric." << std::endl;
461
462 // use type ArPackPlusPlus_BCRSMatrixWrapper to store matrix information
463 // and to perform the product A*v (LU decomposition is not used)
464 typedef ArPackPlusPlus_BCRSMatrixWrapper<BCRSMatrix> WrappedMatrix;
465 WrappedMatrix A(a_);
466
467 // get number of rows and columns in A
468 const int nrows = A.nrows();
469 const int ncols = A.ncols();
470
471 // assert that A is square
472 if (nrows != ncols)
473 DUNE_THROW(Dune::ISTLError,"Matrix is not square ("
474 << nrows << "x" << ncols << ").");
475
476 // allocate memory for variables, set parameters
477 const int nev = 1; // Number of eigenvalues to compute
478 const int ncv = 0; // Number of Arnoldi vectors generated at each iteration (0 == auto)
479 const Real tol = epsilon; // Stopping tolerance (relative accuracy of Ritz values) (0 == machine precision)
480 const int maxit = nIterationsMax_*nev; // Maximum number of Arnoldi update iterations allowed (0 == 100*nev)
481 Real* ev = new Real[nev]; // Computed eigenvalues of A
482 const bool ivec = true; // Flag deciding if eigenvectors shall be determined
483 int nconv; // Number of converged eigenvalues
484
485 // define what we need: eigenvalues with largest magnitude
486 char which[] = "LM";
487 ARSymStdEig<Real,WrappedMatrix>
488 dprob(nrows, nev, &A, &WrappedMatrix::multMv, which, ncv, tol, maxit);
489
490 // set ARPACK verbosity mode if requested
491 if (verbosity_level_ > 3) dprob.Trace();
492
493 // find eigenvalues and eigenvectors of A, obtain the eigenvalues
494 nconv = dprob.Eigenvalues(ev,ivec);
495
496 // obtain approximated dominant eigenvalue of A
497 lambda = ev[nev-1];
498
499 // obtain associated approximated eigenvector of A
500 Real* x_raw = dprob.RawEigenvector(nev-1);
501 WrappedMatrix::arrayToDomainBlockVector(x_raw,x);
502
503 // obtain number of Arnoldi update iterations actually taken
504 nIterations_ = dprob.GetIter();
505
506 // compute residual norm
507 BlockVector r(x);
508 Real* Ax_raw = new Real[nrows];
509 A.multMv(x_raw,Ax_raw);
510 WrappedMatrix::arrayToDomainBlockVector(Ax_raw,r);
511 r.axpy(-lambda,x);
512 const Real r_norm = r.two_norm();
513
514 // print verbosity information
515 if (verbosity_level_ > 0)
516 {
517 if (verbosity_level_ > 1)
518 {
519 // print some information about the problem
520 std::cout << blank_ << "Obtained eigenvalues of A by solving "
521 << "A*x = λ*x using the ARPACK++ class ARSym"
522 << "StdEig:" << std::endl;
523 std::cout << blank_ << " converged eigenvalues of A: "
524 << nconv << " / " << nev << std::endl;
525 std::cout << blank_ << " dominant eigenvalue of A: "
526 << lambda << std::endl;
527 }
528 std::cout << blank_ << "Result ("
529 << "#iterations = " << nIterations_ << ", "
530 << "║residual║_2 = " << r_norm << "): "
531 << "λ = " << lambda << std::endl;
532 if (verbosity_level_ > 2)
533 {
534 // print approximated eigenvector via DUNE-ISTL I/O methods
535 Dune::printvector(std::cout,x,blank_+"x",blank_+"row");
536 }
537 }
538
539 // free dynamically allocated memory
540 delete[] Ax_raw;
541 delete[] ev;
542 }
543
544
545 inline void computeGenSymMaxMagnitude (const BCRSMatrix& b_, const Real& epsilon,
546 BlockVector& x, Real& lambda) const
547 {
548 // print verbosity information
549 if (verbosity_level_ > 0)
550 std::cout << title_ << "Computing an approximation of "
551 << "the dominant eigenvalue of a matrix which "
552 << "is assumed to be symmetric." << std::endl;
553
554 // use type ArPackPlusPlus_BCRSMatrixWrapper to store matrix information
555 // and to perform the product A*v (LU decomposition is not used)
556 typedef ArPackPlusPlus_BCRSMatrixWrapper<BCRSMatrix> WrappedMatrix;
557 WrappedMatrix A(a_);
558 WrappedMatrix B(b_);
559
560 // get number of rows and columns in A
561 const int nrows = A.nrows();
562 const int ncols = A.ncols();
563
564 // assert that A is square
565 if (nrows != ncols)
566 DUNE_THROW(Dune::ISTLError,"Matrix is not square ("
567 << nrows << "x" << ncols << ").");
568
569 // allocate memory for variables, set parameters
570 const int nev = 1; // Number of eigenvalues to compute
571 const int ncv = 0; // Number of Arnoldi vectors generated at each iteration (0 == auto)
572 const Real tol = epsilon; // Stopping tolerance (relative accuracy of Ritz values) (0 == machine precision)
573 const int maxit = nIterationsMax_*nev; // Maximum number of Arnoldi update iterations allowed (0 == 100*nev)
574 Real* ev = new Real[nev]; // Computed eigenvalues of A
575 const bool ivec = true; // Flag deciding if eigenvectors shall be determined
576 int nconv; // Number of converged eigenvalues
577
578 // define what we need: eigenvalues with largest magnitude
579 char which[] = "LM";
580 ARSymGenEig<Real,WrappedMatrix,WrappedMatrix>
581 dprob(nrows, nev, &A, &WrappedMatrix::multMv, &B, &WrappedMatrix::multMv, which, ncv, tol, maxit);
582
583 // set ARPACK verbosity mode if requested
584 if (verbosity_level_ > 3) dprob.Trace();
585
586 // find eigenvalues and eigenvectors of A, obtain the eigenvalues
587 nconv = dprob.Eigenvalues(ev,ivec);
588
589 // obtain approximated dominant eigenvalue of A
590 lambda = ev[nev-1];
591
592 // obtain associated approximated eigenvector of A
593 Real* x_raw = dprob.RawEigenvector(nev-1);
594 WrappedMatrix::arrayToDomainBlockVector(x_raw,x);
595
596 // obtain number of Arnoldi update iterations actually taken
597 nIterations_ = dprob.GetIter();
598
599 // compute residual norm
600 BlockVector r(x);
601 Real* Ax_raw = new Real[nrows];
602 A.multMv(x_raw,Ax_raw);
603 WrappedMatrix::arrayToDomainBlockVector(Ax_raw,r);
604 r.axpy(-lambda,x);
605 const Real r_norm = r.two_norm();
606
607 // print verbosity information
608 if (verbosity_level_ > 0)
609 {
610 if (verbosity_level_ > 1)
611 {
612 // print some information about the problem
613 std::cout << blank_ << "Obtained eigenvalues of A by solving "
614 << "A*x = λ*x using the ARPACK++ class ARSym"
615 << "StdEig:" << std::endl;
616 std::cout << blank_ << " converged eigenvalues of A: "
617 << nconv << " / " << nev << std::endl;
618 std::cout << blank_ << " dominant eigenvalue of A: "
619 << lambda << std::endl;
620 }
621 std::cout << blank_ << "Result ("
622 << "#iterations = " << nIterations_ << ", "
623 << "║residual║_2 = " << r_norm << "): "
624 << "λ = " << lambda << std::endl;
625 if (verbosity_level_ > 2)
626 {
627 // print approximated eigenvector via DUNE-ISTL I/O methods
628 Dune::printvector(std::cout,x,blank_+"x",blank_+"row");
629 }
630 }
631
632 // free dynamically allocated memory
633 delete[] Ax_raw;
634 delete[] ev;
635 }
636
637
649 inline void computeSymMinMagnitude (const Real& epsilon,
650 BlockVector& x, Real& lambda) const
651 {
652 // print verbosity information
653 if (verbosity_level_ > 0)
654 std::cout << title_ << "Computing an approximation of the "
655 << "least dominant eigenvalue of a matrix which "
656 << "is assumed to be symmetric." << std::endl;
657
658 // use type ArPackPlusPlus_BCRSMatrixWrapper to store matrix information
659 // and to perform the product A*v (LU decomposition is not used)
660 typedef ArPackPlusPlus_BCRSMatrixWrapper<BCRSMatrix> WrappedMatrix;
661 WrappedMatrix A(a_);
662
663 // get number of rows and columns in A
664 const int nrows = A.nrows();
665 const int ncols = A.ncols();
666
667 // assert that A is square
668 if (nrows != ncols)
669 DUNE_THROW(Dune::ISTLError,"Matrix is not square ("
670 << nrows << "x" << ncols << ").");
671
672 // allocate memory for variables, set parameters
673 const int nev = 1; // Number of eigenvalues to compute
674 const int ncv = 20; // Number of Arnoldi vectors generated at each iteration (0 == auto)
675 const Real tol = epsilon; // Stopping tolerance (relative accuracy of Ritz values) (0 == machine precision)
676 const int maxit = nIterationsMax_*nev; // Maximum number of Arnoldi update iterations allowed (0 == 100*nev)
677 Real* ev = new Real[nev]; // Computed eigenvalues of A
678 const bool ivec = true; // Flag deciding if eigenvectors shall be determined
679 int nconv; // Number of converged eigenvalues
680 //const Real sigma = 0.0;
681
682 // define what we need: eigenvalues with smallest magnitude
683 char which[] = "SM";
684 ARSymStdEig<Real,WrappedMatrix>
685 //dprob(nrows, nev, &A, &WrappedMatrix::multMv, sigma, which, ncv, tol, maxit);
686 dprob(nrows, nev, &A, &WrappedMatrix::multMv, which, ncv, tol, maxit);
687
688 // set ARPACK verbosity mode if requested
689 if (verbosity_level_ > 3) dprob.Trace();
690
691 // find eigenvalues and eigenvectors of A, obtain the eigenvalues
692 nconv = dprob.Eigenvalues(ev,ivec);
693
694 // obtain approximated least dominant eigenvalue of A
695 lambda = ev[nev-1];
696
697 // obtain associated approximated eigenvector of A
698 Real* x_raw = dprob.RawEigenvector(nev-1);
699 WrappedMatrix::arrayToDomainBlockVector(x_raw,x);
700
701 // obtain number of Arnoldi update iterations actually taken
702 nIterations_ = dprob.GetIter();
703
704 // compute residual norm
705 BlockVector r(x);
706 Real* Ax_raw = new Real[nrows];
707 A.multMv(x_raw,Ax_raw);
708 WrappedMatrix::arrayToDomainBlockVector(Ax_raw,r);
709 r.axpy(-lambda,x);
710 const Real r_norm = r.two_norm();
711
712 // print verbosity information
713 if (verbosity_level_ > 0)
714 {
715 if (verbosity_level_ > 1)
716 {
717 // print some information about the problem
718 std::cout << blank_ << "Obtained eigenvalues of A by solving "
719 << "A*x = λ*x using the ARPACK++ class ARSym"
720 << "StdEig:" << std::endl;
721 std::cout << blank_ << " converged eigenvalues of A: "
722 << nconv << " / " << nev << std::endl;
723 std::cout << blank_ << " least dominant eigenvalue of A: "
724 << lambda << std::endl;
725 }
726 std::cout << blank_ << "Result ("
727 << "#iterations = " << nIterations_ << ", "
728 << "║residual║_2 = " << r_norm << "): "
729 << "λ = " << lambda << std::endl;
730 if (verbosity_level_ > 2)
731 {
732 // print approximated eigenvector via DUNE-ISTL I/O methods
733 Dune::printvector(std::cout,x,blank_+"x",blank_+"row");
734 }
735 }
736
737 // free dynamically allocated memory
738 delete[] Ax_raw;
739 delete[] ev;
740 }
741
742 inline void computeGenSymMinMagnitude (const BCRSMatrix& b_, const Real& epsilon, int& nev,
743 std::vector<BlockVector>& x, std::vector<Real>& lambda, Real sigma) const
744 {
745 // print verbosity information
746 if (verbosity_level_ > 0)
747 std::cout << title_ << "Computing an approximation of the "
748 << "least dominant eigenvalue of a matrix which "
749 << "is assumed to be symmetric." << std::endl;
750
751
752
753 // allocate memory for variables, set parameters
754 //const int nev = 1; // Number of eigenvalues to compute
755 const int ncv = 0; // Number of Arnoldi vectors generated at each iteration (0 == auto)
756 const Real tol = epsilon; // Stopping tolerance (relative accuracy of Ritz values) (0 == machine precision)
757 const int maxit = nIterationsMax_*nev; // Maximum number of Arnoldi update iterations allowed (0 == 100*nev)
758 Real* ev = new Real[nev]; // Computed eigenvalues of A
759 Real* ev_imag = new Real[nev]; // Computed eigenvalues of A
760 const bool ivec = true; // Flag deciding if eigenvectors shall be determined
761 int nconv; // Number of converged eigenvalues
762 char invertmodep = 'S';
763
764 BCRSMatrix ashiftb(a_);
765 ashiftb.axpy(-sigma,b_);
766
767 // use type ArPackPlusPlus_BCRSMatrixWrapper to store matrix information
768 // and to perform the product (A-sigma B)^-1 v (LU decomposition is not used)
769 typedef ArPackPlusPlus_BCRSMatrixWrapper2<BCRSMatrix> WrappedMatrix;
770 WrappedMatrix A(ashiftb,b_);
771
772 // get number of rows and columns in A
773 const int nrows = A.nrows();
774 const int ncols = A.ncols();
775
776 // assert that A is square
777 if (nrows != ncols)
778 DUNE_THROW(Dune::ISTLError,"Matrix is not square ("
779 << nrows << "x" << ncols << ").");
780
781 // define what we need: eigenvalues with smallest magnitude
782 char which[] = "LM";
783 //ARNonSymGenEig<Real,WrappedMatrix,WrappedMatrix>
784 // dprob(nrows, nev, &A, &WrappedMatrix::multMv, &A, &WrappedMatrix::multMvB, sigma, which, ncv, tol, maxit);
785
786 //Non generalised problem
787 ARSymStdEig<Real,WrappedMatrix>
788 dprob(nrows, nev, &A, &WrappedMatrix::multMv, which, ncv, tol, maxit);
789
790
791 // set ARPACK verbosity mode if requested
792 if (verbosity_level_ > 3) dprob.Trace();
793
794 // find eigenvalues and eigenvectors of A
795 nconv = dprob.Eigenvalues(ev,ev_imag,ivec);
796
797 // Get sorting permutation for un-shifted eigenvalues
798 std::vector<int> index(nev, 0);
799 for (int i = 0 ; i != index.size() ; i++) {
800 index[i] = i;
801 }
802 std::sort(index.begin(), index.end(),
803 [&](const int& a, const int& b) {
804 return (sigma+1./ev[a] < sigma+1./ev[b]);
805 }
806 );
807
808 // obtain approximated least dominant eigenvalue of A
809 lambda.resize(nev);
810 // obtain associated approximated eigenvector of A
811 x.resize(nev);
812 for (int i = 0; i < nev; i++) {
813 lambda[i] = sigma+1./ev[index[i]];
814 Real* x_raw = dprob.RawEigenvector(index[i]);
815 WrappedMatrix::arrayToDomainBlockVector(x_raw,x[i]);
816 }
817
818 // obtain number of Arnoldi update iterations actually taken
819 nIterations_ = dprob.GetIter();
820
821 // compute residual norm
822 /*BlockVector r(x);
823 Real* Ax_raw = new Real[nrows];
824 A.multMv(x_raw,Ax_raw);
825 WrappedMatrix::arrayToDomainBlockVector(Ax_raw,r);
826 r.axpy(-lambda,x);
827 const Real r_norm = r.two_norm();
828
829 // print verbosity information
830 if (verbosity_level_ > 0)
831 {
832 if (verbosity_level_ > 1)
833 {
834 // print some information about the problem
835 std::cout << blank_ << "Obtained eigenvalues of A by solving "
836 << "A*x = λ*x using the ARPACK++ class ARSym"
837 << "StdEig:" << std::endl;
838 std::cout << blank_ << " converged eigenvalues of A: "
839 << nconv << " / " << nev << std::endl;
840 std::cout << blank_ << " least dominant eigenvalue of A: "
841 << lambda << std::endl;
842 }
843 std::cout << blank_ << "Result ("
844 << "#iterations = " << nIterations_ << ", "
845 << "║residual║_2 = " << r_norm << "): "
846 << "λ = " << lambda << std::endl;
847 if (verbosity_level_ > 2)
848 {
849 // print approximated eigenvector via DUNE-ISTL I/O methods
850 Dune::printvector(std::cout,x,blank_+"x",blank_+"row");
851 }
852 }*/
853
854 // free dynamically allocated memory
855 //delete[] Ax_raw;
856 delete[] ev;
857 }
858
859
860 inline void computeGenNonSymMinMagnitude (const BCRSMatrix& b_, const Real& epsilon, int& nev,
861 std::vector<BlockVector>& x, std::vector<Real>& lambda, Real sigma) const
862 {
863 // print verbosity information
864 if (verbosity_level_ > 0)
865 std::cout << title_ << "Computing an approximation of the "
866 << "least dominant eigenvalue of a matrix which "
867 << "is assumed to be symmetric." << std::endl;
868
869
870
871 // allocate memory for variables, set parameters
872 //const int nev = 1; // Number of eigenvalues to compute
873 const int ncv = 0; // Number of Arnoldi vectors generated at each iteration (0 == auto)
874 const Real tol = epsilon; // Stopping tolerance (relative accuracy of Ritz values) (0 == machine precision)
875 const int maxit = nIterationsMax_*nev; // Maximum number of Arnoldi update iterations allowed (0 == 100*nev)
876 Real* ev = new Real[nev]; // Computed eigenvalues of A
877 Real* ev_imag = new Real[nev]; // Computed eigenvalues of A
878 const bool ivec = true; // Flag deciding if eigenvectors shall be determined
879 int nconv; // Number of converged eigenvalues
880 //char invertmodep = 'S';
881
882 BCRSMatrix ashiftb(a_);
883 ashiftb.axpy(-sigma,b_);
884
885 // use type ArPackPlusPlus_BCRSMatrixWrapper to store matrix information
886 // and to perform the product (A-sigma B)^-1 v (LU decomposition is not used)
887 typedef ArPackPlusPlus_BCRSMatrixWrapper2<BCRSMatrix> WrappedMatrix;
888 WrappedMatrix A(ashiftb,b_);
889
890 // get number of rows and columns in A
891 const int nrows = A.nrows();
892 const int ncols = A.ncols();
893
894 // assert that A is square
895 if (nrows != ncols)
896 DUNE_THROW(Dune::ISTLError,"Matrix is not square ("
897 << nrows << "x" << ncols << ").");
898
899 // define what we need: eigenvalues with smallest magnitude
900 char which[] = "LM";
901 //ARNonSymGenEig<Real,WrappedMatrix,WrappedMatrix>
902 // dprob(nrows, nev, &A, &WrappedMatrix::multMv, &A, &WrappedMatrix::multMvB, sigma, which, ncv, tol, maxit);
903
904 //Non generalised problem
905 ARNonSymStdEig<Real,WrappedMatrix>
906 dprob(nrows, nev, &A, &WrappedMatrix::multMv, which, ncv, tol, maxit);
907
908
909 // set ARPACK verbosity mode if requested
910 if (verbosity_level_ > 3) dprob.Trace();
911
912 // find eigenvalues and eigenvectors of A
913 nconv = dprob.Eigenvalues(ev,ev_imag,ivec);
914 if (verbosity_level_ > 3) std::cout << nconv << std::endl;
915
916 // Get sorting permutation for un-shifted eigenvalues
917 std::vector<int> index(nev, 0);
918 for (unsigned int i = 0 ; i != index.size() ; i++) {
919 index[i] = i;
920 }
921 std::sort(index.begin(), index.end(),
922 [&](const int& a, const int& b) {
923 return (sigma+1./ev[a] < sigma+1./ev[b]);
924 }
925 );
926
927 // obtain approximated least dominant eigenvalue of A
928 lambda.resize(nev);
929 // obtain associated approximated eigenvector of A
930 x.resize(nev);
931 for (int i = 0; i < nev; i++) {
932 lambda[i] = sigma+1./ev[index[i]];
933 Real* x_raw = dprob.RawEigenvector(index[i]);
934 WrappedMatrix::arrayToDomainBlockVector(x_raw,x[i]);
935 }
936
937 // obtain number of Arnoldi update iterations actually taken
938 nIterations_ = dprob.GetIter();
939
940 // compute residual norm
941 /*BlockVector r(x);
942 Real* Ax_raw = new Real[nrows];
943 A.multMv(x_raw,Ax_raw);
944 WrappedMatrix::arrayToDomainBlockVector(Ax_raw,r);
945 r.axpy(-lambda,x);
946 const Real r_norm = r.two_norm();
947
948 // print verbosity information
949 if (verbosity_level_ > 0)
950 {
951 if (verbosity_level_ > 1)
952 {
953 // print some information about the problem
954 std::cout << blank_ << "Obtained eigenvalues of A by solving "
955 << "A*x = λ*x using the ARPACK++ class ARSym"
956 << "StdEig:" << std::endl;
957 std::cout << blank_ << " converged eigenvalues of A: "
958 << nconv << " / " << nev << std::endl;
959 std::cout << blank_ << " least dominant eigenvalue of A: "
960 << lambda << std::endl;
961 }
962 std::cout << blank_ << "Result ("
963 << "#iterations = " << nIterations_ << ", "
964 << "║residual║_2 = " << r_norm << "): "
965 << "λ = " << lambda << std::endl;
966 if (verbosity_level_ > 2)
967 {
968 // print approximated eigenvector via DUNE-ISTL I/O methods
969 Dune::printvector(std::cout,x,blank_+"x",blank_+"row");
970 }
971 }*/
972
973 // free dynamically allocated memory
974 //delete[] Ax_raw;
975 delete[] ev;
976 }
977
989 inline void computeSymCond2 (const Real& epsilon, Real& cond_2) const
990 {
991 // print verbosity information
992 if (verbosity_level_ > 0)
993 std::cout << title_ << "Computing an approximation of the "
994 << "spectral condition number of a matrix which "
995 << "is assumed to be symmetric." << std::endl;
996
997 // use type ArPackPlusPlus_BCRSMatrixWrapper to store matrix information
998 // and to perform the product A*v (LU decomposition is not used)
999 typedef ArPackPlusPlus_BCRSMatrixWrapper<BCRSMatrix> WrappedMatrix;
1000 WrappedMatrix A(a_);
1001
1002 // get number of rows and columns in A
1003 const int nrows = A.nrows();
1004 const int ncols = A.ncols();
1005
1006 // assert that A is square
1007 if (nrows != ncols)
1008 DUNE_THROW(Dune::ISTLError,"Matrix is not square ("
1009 << nrows << "x" << ncols << ").");
1010
1011 // allocate memory for variables, set parameters
1012 const int nev = 2; // Number of eigenvalues to compute
1013 const int ncv = 20; // Number of Arnoldi vectors generated at each iteration (0 == auto)
1014 const Real tol = epsilon; // Stopping tolerance (relative accuracy of Ritz values) (0 == machine precision)
1015 const int maxit = nIterationsMax_*nev; // Maximum number of Arnoldi update iterations allowed (0 == 100*nev)
1016 Real* ev = new Real[nev]; // Computed eigenvalues of A
1017 const bool ivec = true; // Flag deciding if eigenvectors shall be determined
1018 int nconv; // Number of converged eigenvalues
1019
1020 // define what we need: eigenvalues from both ends of the spectrum
1021 char which[] = "BE";
1022 ARSymStdEig<Real,WrappedMatrix>
1023 dprob(nrows, nev, &A, &WrappedMatrix::multMv, which, ncv, tol, maxit);
1024
1025 // set ARPACK verbosity mode if requested
1026 if (verbosity_level_ > 3) dprob.Trace();
1027
1028 // find eigenvalues and eigenvectors of A, obtain the eigenvalues
1029 nconv = dprob.Eigenvalues(ev,ivec);
1030
1031 // obtain approximated dominant and least dominant eigenvalue of A
1032 const Real& lambda_max = ev[nev-1];
1033 const Real& lambda_min = ev[0];
1034
1035 // obtain associated approximated eigenvectors of A
1036 Real* x_max_raw = dprob.RawEigenvector(nev-1);
1037 Real* x_min_raw = dprob.RawEigenvector(0);
1038
1039 // obtain approximated spectral condition number of A
1040 cond_2 = std::abs(lambda_max / lambda_min);
1041
1042 // obtain number of Arnoldi update iterations actually taken
1043 nIterations_ = dprob.GetIter();
1044
1045 // compute each residual norm
1046 Real* Ax_max_raw = new Real[nrows];
1047 Real* Ax_min_raw = new Real[nrows];
1048 A.multMv(x_max_raw,Ax_max_raw);
1049 A.multMv(x_min_raw,Ax_min_raw);
1050 Real r_max_norm = 0.0;
1051 Real r_min_norm = 0.0;
1052 for (int i = 0; i < nrows; ++i)
1053 {
1054 r_max_norm += std::pow(Ax_max_raw[i] - lambda_max * x_max_raw[i],2);
1055 r_min_norm += std::pow(Ax_min_raw[i] - lambda_min * x_min_raw[i],2);
1056 }
1057 r_max_norm = std::sqrt(r_max_norm);
1058 r_min_norm = std::sqrt(r_min_norm);
1059
1060 // print verbosity information
1061 if (verbosity_level_ > 0)
1062 {
1063 if (verbosity_level_ > 1)
1064 {
1065 // print some information about the problem
1066 std::cout << blank_ << "Obtained eigenvalues of A by solving "
1067 << "A*x = λ*x using the ARPACK++ class ARSym"
1068 << "StdEig:" << std::endl;
1069 std::cout << blank_ << " converged eigenvalues of A: "
1070 << nconv << " / " << nev << std::endl;
1071 std::cout << blank_ << " dominant eigenvalue of A: "
1072 << lambda_max << std::endl;
1073 std::cout << blank_ << " least dominant eigenvalue of A: "
1074 << lambda_min << std::endl;
1075 std::cout << blank_ << " spectral condition number of A: "
1076 << cond_2 << std::endl;
1077 }
1078 std::cout << blank_ << "Result ("
1079 << "#iterations = " << nIterations_ << ", "
1080 << "║residual║_2 = {" << r_max_norm << ","
1081 << r_min_norm << "}, " << "λ = {"
1082 << lambda_max << "," << lambda_min
1083 << "}): cond_2 = " << cond_2 << std::endl;
1084 }
1085
1086 // free dynamically allocated memory
1087 delete[] Ax_min_raw;
1088 delete[] Ax_max_raw;
1089 delete[] ev;
1090 }
1091
1105 inline void computeNonSymMax (const Real& epsilon,
1106 BlockVector& x, Real& sigma) const
1107 {
1108 // print verbosity information
1109 if (verbosity_level_ > 0)
1110 std::cout << title_ << "Computing an approximation of the "
1111 << "largest singular value of a matrix which "
1112 << "is assumed to be nonsymmetric." << std::endl;
1113
1114 // use type ArPackPlusPlus_BCRSMatrixWrapper to store matrix information
1115 // and to perform the product A^T*A*v (LU decomposition is not used)
1116 typedef ArPackPlusPlus_BCRSMatrixWrapper<BCRSMatrix> WrappedMatrix;
1117 WrappedMatrix A(a_);
1118
1119 // get number of rows and columns in A
1120 const int nrows = A.nrows();
1121 const int ncols = A.ncols();
1122
1123 // assert that A has more rows than columns (extend code later to the opposite case!)
1124 if (nrows < ncols)
1125 DUNE_THROW(Dune::ISTLError,"Matrix has less rows than "
1126 << "columns (" << nrows << "x" << ncols << ")."
1127 << " This case is not implemented, yet.");
1128
1129 // allocate memory for variables, set parameters
1130 const int nev = 1; // Number of eigenvalues to compute
1131 const int ncv = 20; // Number of Arnoldi vectors generated at each iteration (0 == auto)
1132 const Real tol = epsilon; // Stopping tolerance (relative accuracy of Ritz values) (0 == machine precision)
1133 const int maxit = nIterationsMax_*nev; // Maximum number of Arnoldi update iterations allowed (0 == 100*nev)
1134 Real* ev = new Real[nev]; // Computed eigenvalues of A^T*A
1135 const bool ivec = true; // Flag deciding if eigenvectors shall be determined
1136 int nconv; // Number of converged eigenvalues
1137
1138 // define what we need: eigenvalues with largest algebraic value
1139 char which[] = "LA";
1140 ARSymStdEig<Real,WrappedMatrix>
1141 dprob(ncols, nev, &A, &WrappedMatrix::multMtMv, which, ncv, tol, maxit);
1142
1143 // set ARPACK verbosity mode if requested
1144 if (verbosity_level_ > 3) dprob.Trace();
1145
1146 // find eigenvalues and eigenvectors of A^T*A, obtain the eigenvalues
1147 nconv = dprob.Eigenvalues(ev,ivec);
1148
1149 // obtain approximated largest eigenvalue of A^T*A
1150 const Real& lambda = ev[nev-1];
1151
1152 // obtain associated approximated eigenvector of A^T*A
1153 Real* x_raw = dprob.RawEigenvector(nev-1);
1154 WrappedMatrix::arrayToDomainBlockVector(x_raw,x);
1155
1156 // obtain number of Arnoldi update iterations actually taken
1157 nIterations_ = dprob.GetIter();
1158
1159 // compute residual norm
1160 BlockVector r(x);
1161 Real* AtAx_raw = new Real[ncols];
1162 A.multMtMv(x_raw,AtAx_raw);
1163 WrappedMatrix::arrayToDomainBlockVector(AtAx_raw,r);
1164 r.axpy(-lambda,x);
1165 const Real r_norm = r.two_norm();
1166
1167 // calculate largest singular value of A (note that
1168 // x is right-singular / left-singular vector of A)
1169 sigma = std::sqrt(lambda);
1170
1171 // print verbosity information
1172 if (verbosity_level_ > 0)
1173 {
1174 if (verbosity_level_ > 1)
1175 {
1176 // print some information about the problem
1177 std::cout << blank_ << "Obtained singular values of A by sol"
1178 << "ving (A^T*A)*x = σ²*x using the ARPACK++ "
1179 << "class ARSymStdEig:" << std::endl;
1180 std::cout << blank_ << " converged eigenvalues of A^T*A: "
1181 << nconv << " / " << nev << std::endl;
1182 std::cout << blank_ << " largest eigenvalue of A^T*A: "
1183 << lambda << std::endl;
1184 std::cout << blank_ << " => largest singular value of A: "
1185 << sigma << std::endl;
1186 }
1187 std::cout << blank_ << "Result ("
1188 << "#iterations = " << nIterations_ << ", "
1189 << "║residual║_2 = " << r_norm << "): "
1190 << "σ = " << sigma << std::endl;
1191 if (verbosity_level_ > 2)
1192 {
1193 // print approximated right-singular / left-singular vector
1194 // via DUNE-ISTL I/O methods
1195 Dune::printvector(std::cout,x,blank_+"x",blank_+"row");
1196 }
1197 }
1198
1199 // free dynamically allocated memory
1200 delete[] AtAx_raw;
1201 delete[] ev;
1202 }
1203
1217 inline void computeNonSymMin (const BCRSMatrix& b_, const Real& epsilon,
1218 BlockVector& x, Real& sigma) const
1219 {
1220 // print verbosity information
1221 if (verbosity_level_ > 0)
1222 std::cout << title_ << "Computing an approximation of the "
1223 << "smallest singular value of a matrix which "
1224 << "is assumed to be nonsymmetric." << std::endl;
1225
1226 // use type ArPackPlusPlus_BCRSMatrixWrapper to store matrix information
1227 // and to perform the product A^T*A*v (LU decomposition is not used)
1228 typedef ArPackPlusPlus_BCRSMatrixWrapper<BCRSMatrix> WrappedMatrix;
1229 WrappedMatrix A(a_);
1230 WrappedMatrix B(b_);
1231
1232 // get number of rows and columns in A
1233 const int nrows = A.nrows();
1234 const int ncols = A.ncols();
1235
1236 // assert that A has more rows than columns (extend code later to the opposite case!)
1237 if (nrows < ncols)
1238 DUNE_THROW(Dune::ISTLError,"Matrix has less rows than "
1239 << "columns (" << nrows << "x" << ncols << ")."
1240 << " This case is not implemented, yet.");
1241
1242 // allocate memory for variables, set parameters
1243 const int nev = 2; // Number of eigenvalues to compute
1244 const int ncv = 20; // Number of Arnoldi vectors generated at each iteration (0 == auto)
1245 const Real tol = epsilon; // Stopping tolerance (relative accuracy of Ritz values) (0 == machine precision)
1246 const int maxit = nIterationsMax_*nev; // Maximum number of Arnoldi update iterations allowed (0 == 100*nev)
1247 Real* ev = new Real[nev]; // Computed eigenvalues of A^T*A
1248 Real* ev_imag = new Real[nev]; // Computed eigenvalues of A^T*A
1249 const bool ivec = true; // Flag deciding if eigenvectors shall be determined
1250 int nconv; // Number of converged eigenvalues
1251
1252 // define what we need: eigenvalues with smallest algebraic value
1253 char which[] = "SM";
1254 ARNonSymGenEig<Real,WrappedMatrix,WrappedMatrix>
1255 dprob(ncols, nev, &A, &WrappedMatrix::multMv, &B, &WrappedMatrix::multMv, which, ncv, tol, maxit);
1256
1257 // set ARPACK verbosity mode if requested
1258 if (verbosity_level_ > 3) dprob.Trace();
1259
1260 // find eigenvalues and eigenvectors of A^T*A, obtain the eigenvalues
1261 nconv = dprob.Eigenvalues(ev,ev_imag,ivec);
1262
1263 // obtain approximated smallest eigenvalue of A^T*A
1264 const Real& lambda = ev[nev-1];
1265
1266 // obtain associated approximated eigenvector of A^T*A
1267 Real* x_raw = dprob.RawEigenvector(nev-1);
1268 WrappedMatrix::arrayToDomainBlockVector(x_raw,x);
1269
1270 // obtain number of Arnoldi update iterations actually taken
1271 nIterations_ = dprob.GetIter();
1272
1273 // compute residual norm
1274 BlockVector r(x);
1275 Real* AtAx_raw = new Real[ncols];
1276 A.multMtMv(x_raw,AtAx_raw);
1277 WrappedMatrix::arrayToDomainBlockVector(AtAx_raw,r);
1278 r.axpy(-lambda,x);
1279 const Real r_norm = r.two_norm();
1280
1281 // calculate smallest singular value of A (note that
1282 // x is right-singular / left-singular vector of A)
1283 sigma = std::sqrt(lambda);
1284
1285 // print verbosity information
1286 if (verbosity_level_ > 0)
1287 {
1288 if (verbosity_level_ > 1)
1289 {
1290 // print some information about the problem
1291 std::cout << blank_ << "Obtained singular values of A by sol"
1292 << "ving (A^T*A)*x = σ²*x using the ARPACK++ "
1293 << "class ARSymStdEig:" << std::endl;
1294 std::cout << blank_ << " converged eigenvalues of A^T*A: "
1295 << nconv << " / " << nev << std::endl;
1296 std::cout << blank_ << " smallest eigenvalue of A^T*A: "
1297 << lambda << std::endl;
1298 std::cout << blank_ << " => smallest singular value of A: "
1299 << sigma << std::endl;
1300 }
1301 std::cout << blank_ << "Result ("
1302 << "#iterations = " << nIterations_ << ", "
1303 << "║residual║_2 = " << r_norm << "): "
1304 << "σ = " << sigma << std::endl;
1305 if (verbosity_level_ > 2)
1306 {
1307 // print approximated right-singular / left-singular vector
1308 // via DUNE-ISTL I/O methods
1309 Dune::printvector(std::cout,x,blank_+"x",blank_+"row");
1310 }
1311 }
1312
1313 // free dynamically allocated memory
1314 delete[] AtAx_raw;
1315 delete[] ev;
1316 }
1317
1328 inline void computeNonSymCond2 (const Real& epsilon, Real& cond_2) const
1329 {
1330 // print verbosity information
1331 if (verbosity_level_ > 0)
1332 std::cout << title_ << "Computing an approximation of the "
1333 << "spectral condition number of a matrix which "
1334 << "is assumed to be nonsymmetric." << std::endl;
1335
1336 // use type ArPackPlusPlus_BCRSMatrixWrapper to store matrix information
1337 // and to perform the product A^T*A*v (LU decomposition is not used)
1338 typedef ArPackPlusPlus_BCRSMatrixWrapper<BCRSMatrix> WrappedMatrix;
1339 WrappedMatrix A(a_);
1340
1341 // get number of rows and columns in A
1342 const int nrows = A.nrows();
1343 const int ncols = A.ncols();
1344
1345 // assert that A has more rows than columns (extend code later to the opposite case!)
1346 if (nrows < ncols)
1347 DUNE_THROW(Dune::ISTLError,"Matrix has less rows than "
1348 << "columns (" << nrows << "x" << ncols << ")."
1349 << " This case is not implemented, yet.");
1350
1351 // allocate memory for variables, set parameters
1352 const int nev = 2; // Number of eigenvalues to compute
1353 const int ncv = 20; // Number of Arnoldi vectors generated at each iteration (0 == auto)
1354 const Real tol = epsilon; // Stopping tolerance (relative accuracy of Ritz values) (0 == machine precision)
1355 const int maxit = nIterationsMax_*nev; // Maximum number of Arnoldi update iterations allowed (0 == 100*nev)
1356 Real* ev = new Real[nev]; // Computed eigenvalues of A^T*A
1357 const bool ivec = true; // Flag deciding if eigenvectors shall be determined
1358 int nconv; // Number of converged eigenvalues
1359
1360 // define what we need: eigenvalues from both ends of the spectrum
1361 char which[] = "BE";
1362 ARSymStdEig<Real,WrappedMatrix>
1363 dprob(ncols, nev, &A, &WrappedMatrix::multMtMv, which, ncv, tol, maxit);
1364
1365 // set ARPACK verbosity mode if requested
1366 if (verbosity_level_ > 3) dprob.Trace();
1367
1368 // find eigenvalues and eigenvectors of A^T*A, obtain the eigenvalues
1369 nconv = dprob.Eigenvalues(ev,ivec);
1370
1371 // obtain approximated largest and smallest eigenvalue of A^T*A
1372 const Real& lambda_max = ev[nev-1];
1373 const Real& lambda_min = ev[0];
1374
1375 // obtain associated approximated eigenvectors of A^T*A
1376 Real* x_max_raw = dprob.RawEigenvector(nev-1);
1377 Real* x_min_raw = dprob.RawEigenvector(0);
1378
1379 // obtain number of Arnoldi update iterations actually taken
1380 nIterations_ = dprob.GetIter();
1381
1382 // compute each residual norm
1383 Real* AtAx_max_raw = new Real[ncols];
1384 Real* AtAx_min_raw = new Real[ncols];
1385 A.multMtMv(x_max_raw,AtAx_max_raw);
1386 A.multMtMv(x_min_raw,AtAx_min_raw);
1387 Real r_max_norm = 0.0;
1388 Real r_min_norm = 0.0;
1389 for (int i = 0; i < ncols; ++i)
1390 {
1391 r_max_norm += std::pow(AtAx_max_raw[i] - lambda_max * x_max_raw[i],2);
1392 r_min_norm += std::pow(AtAx_min_raw[i] - lambda_min * x_min_raw[i],2);
1393 }
1394 r_max_norm = std::sqrt(r_max_norm);
1395 r_min_norm = std::sqrt(r_min_norm);
1396
1397 // calculate largest and smallest singular value of A
1398 const Real sigma_max = std::sqrt(lambda_max);
1399 const Real sigma_min = std::sqrt(lambda_min);
1400
1401 // obtain approximated spectral condition number of A
1402 cond_2 = sigma_max / sigma_min;
1403
1404 // print verbosity information
1405 if (verbosity_level_ > 0)
1406 {
1407 if (verbosity_level_ > 1)
1408 {
1409 // print some information about the problem
1410 std::cout << blank_ << "Obtained singular values of A by sol"
1411 << "ving (A^T*A)*x = σ²*x using the ARPACK++ "
1412 << "class ARSymStdEig:" << std::endl;
1413 std::cout << blank_ << " converged eigenvalues of A^T*A: "
1414 << nconv << " / " << nev << std::endl;
1415 std::cout << blank_ << " largest eigenvalue of A^T*A: "
1416 << lambda_max << std::endl;
1417 std::cout << blank_ << " smallest eigenvalue of A^T*A: "
1418 << lambda_min << std::endl;
1419 std::cout << blank_ << " => largest singular value of A: "
1420 << sigma_max << std::endl;
1421 std::cout << blank_ << " => smallest singular value of A: "
1422 << sigma_min << std::endl;
1423 }
1424 std::cout << blank_ << "Result ("
1425 << "#iterations = " << nIterations_ << ", "
1426 << "║residual║_2 = {" << r_max_norm << ","
1427 << r_min_norm << "}, " << "σ = {"
1428 << sigma_max << "," << sigma_min
1429 << "}): cond_2 = " << cond_2 << std::endl;
1430 }
1431
1432 // free dynamically allocated memory
1433 delete[] AtAx_min_raw;
1434 delete[] AtAx_max_raw;
1435 delete[] ev;
1436 }
1437
1442 inline unsigned int getIterationCount () const
1443 {
1444 if (nIterations_ == 0)
1445 DUNE_THROW(Dune::ISTLError,"No algorithm applied, yet.");
1446
1447 return nIterations_;
1448 }
1449
1450 protected:
1451 // parameters related to iterative eigenvalue algorithms
1452 const BCRSMatrix& a_;
1453 const unsigned int nIterationsMax_;
1454
1455 // verbosity setting
1456 const unsigned int verbosity_level_;
1457
1458 // memory for storing temporary variables (mutable as they shall
1459 // just be effectless auxiliary variables of the const apply*(...)
1460 // methods)
1461 mutable unsigned int nIterations_;
1462
1463 // constants for printing verbosity information
1464 const std::string title_;
1465 const std::string blank_;
1466 };
1467
1468} // namespace Dune
1469
1470
1471#endif // HAVE_ARPACKPP
1472
1473#endif // DUNE_ISTL_EIGENVALUE_ARPACKPP_HH_FORK
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden & Uni Heidelberg  |  generated with Hugo v0.111.3 (Sep 5, 22:35, 2025)