Reference documentation for deal.II version 8.4.2
arpack_solver.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2010 - 2016 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii__arpack_solver_h
17 #define dealii__arpack_solver_h
18 
19 #include <deal.II/base/config.h>
20 #include <deal.II/base/smartpointer.h>
21 #include <deal.II/lac/solver_control.h>
22 
23 #include <cstring>
24 
25 
26 #ifdef DEAL_II_WITH_ARPACK
27 
28 DEAL_II_NAMESPACE_OPEN
29 
30 
31 extern "C" void dnaupd_(int *ido, char *bmat, const unsigned int *n, char *which,
32  const unsigned int *nev, const double *tol, double *resid, int *ncv,
33  double *v, int *ldv, int *iparam, int *ipntr,
34  double *workd, double *workl, int *lworkl,
35  int *info);
36 
37 extern "C" void dneupd_(int *rvec, char *howmany, int *select, double *d,
38  double *di, double *z, int *ldz, double *sigmar,
39  double *sigmai, double *workev, char *bmat,const unsigned int *n, char *which,
40  const unsigned int *nev, const double *tol, double *resid, int *ncv,
41  double *v, int *ldv, int *iparam, int *ipntr,
42  double *workd, double *workl, int *lworkl, int *info);
43 
90 class ArpackSolver : public Subscriptor
91 {
92 public:
97 
98 
104  {
105  algebraically_largest,
106  algebraically_smallest,
107  largest_magnitude,
108  smallest_magnitude,
109  largest_real_part,
110  smallest_real_part,
111  largest_imaginary_part,
112  smallest_imaginary_part,
113  both_ends
114  };
115 
121  {
122  const unsigned int number_of_arnoldi_vectors;
123  const WhichEigenvalues eigenvalue_of_interest;
125  const unsigned int number_of_arnoldi_vectors = 15,
126  const WhichEigenvalues eigenvalue_of_interest = largest_magnitude);
127  };
128 
132  SolverControl &control () const;
133 
137  ArpackSolver(SolverControl &control,
138  const AdditionalData &data = AdditionalData());
139 
175  template <typename VectorType, typename MatrixType1,
176  typename MatrixType2, typename INVERSE>
177  void solve (const MatrixType1 &A,
178  const MatrixType2 &B,
179  const INVERSE &inverse,
180  std::vector<std::complex<double> > &eigenvalues,
181  std::vector<VectorType> &eigenvectors,
182  const unsigned int n_eigenvalues = 0);
183 
184 protected:
185 
191 
196 
197 private:
198 
202  DeclException2 (ExcInvalidNumberofEigenvalues, int, int,
203  << "Number of wanted eigenvalues " << arg1
204  << " is larger that the size of the matrix " << arg2);
205 
206  DeclException2 (ExcInvalidNumberofArnoldiVectors, int, int,
207  << "Number of Arnoldi vectors " << arg1
208  << " is larger that the size of the matrix " << arg2);
209 
210  DeclException2 (ExcSmallNumberofArnoldiVectors, int, int,
211  << "Number of Arnoldi vectors " << arg1
212  << " is too small to obtain " << arg2
213  << " eigenvalues");
214 
215  DeclException1 (ExcArpackIdo, int, << "This ido " << arg1
216  << " is not supported. Check documentation of ARPACK");
217 
218  DeclException1 (ExcArpackMode, int, << "This mode " << arg1
219  << " is not supported. Check documentation of ARPACK");
220 
221  DeclException1 (ExcArpackInfodsaupd, int,
222  << "Error with dsaupd, info " << arg1
223  << ". Check documentation of ARPACK");
224 
225  DeclException1 (ExcArpackInfodneupd, int,
226  << "Error with dneupd, info " << arg1
227  << ". Check documentation of ARPACK");
228 
229  DeclException1 (ExcArpackInfoMaxIt, int,
230  << "Maximum number " << arg1
231  << " of iterations reached.");
232 
233  DeclExceptionMsg (ExcArpackNoShifts,
234  "No shifts could be applied during implicit"
235  " Arnoldi update, try increasing the number of"
236  " Arnoldi vectors.");
237 };
238 
239 
240 inline
241 ArpackSolver::AdditionalData::
242 AdditionalData (const unsigned int number_of_arnoldi_vectors,
243  const WhichEigenvalues eigenvalue_of_interest)
244  :
245  number_of_arnoldi_vectors(number_of_arnoldi_vectors),
246  eigenvalue_of_interest(eigenvalue_of_interest)
247 {}
248 
249 
250 inline
252  const AdditionalData &data)
253  :
254  solver_control (control),
255  additional_data (data)
256 
257 {}
258 
259 
260 template <typename VectorType, typename MatrixType1,
261  typename MatrixType2, typename INVERSE>
262 inline
263 void ArpackSolver::solve (const MatrixType1 &/*system_matrix*/,
264  const MatrixType2 &mass_matrix,
265  const INVERSE &inverse,
266  std::vector<std::complex<double> > &eigenvalues,
267  std::vector<VectorType> &eigenvectors,
268  const unsigned int n_eigenvalues)
269 {
270  //inside the routines of ARPACK the
271  //values change magically, so store
272  //them here
273 
274  const unsigned int n = eigenvectors[0].size();
275  const unsigned int n_inside_arpack = eigenvectors[0].size();
276 
277  // Number of eigenvalues for arpack
278  const unsigned int nev = (n_eigenvalues == 0) ? eigenvalues.size() : n_eigenvalues;
279  AssertIndexRange(eigenvalues.size()-1, nev);
280  /*
281  if(n < 0 || nev <0 || p < 0 || maxit < 0 )
282  std:cout << "All input parameters have to be positive.\n";
283  */
284  Assert (n_eigenvalues < n,
285  ExcInvalidNumberofEigenvalues(nev, n));
286 
287  Assert (additional_data.number_of_arnoldi_vectors < n,
288  ExcInvalidNumberofArnoldiVectors(
289  additional_data.number_of_arnoldi_vectors, n));
290 
291  Assert (additional_data.number_of_arnoldi_vectors > 2*nev+1,
292  ExcSmallNumberofArnoldiVectors(
293  additional_data.number_of_arnoldi_vectors, nev));
294  // ARPACK mode for dnaupd, here only mode 3
295  int mode = 3;
296 
297  // reverse communication parameter
298  int ido = 0;
299 
303  char bmat[2] = "G";
304 
312  char which[3];
313  switch (additional_data.eigenvalue_of_interest)
314  {
315  case algebraically_largest:
316  std::strcpy (which, "LA");
317  break;
318  case algebraically_smallest:
319  std::strcpy (which, "SA");
320  break;
321  case largest_magnitude:
322  std::strcpy (which, "LM");
323  break;
324  case smallest_magnitude:
325  std::strcpy (which, "SM");
326  break;
327  case largest_real_part:
328  std::strcpy (which, "LR");
329  break;
330  case smallest_real_part:
331  std::strcpy (which, "SR");
332  break;
333  case largest_imaginary_part:
334  std::strcpy (which, "LI");
335  break;
336  case smallest_imaginary_part:
337  std::strcpy (which, "SI");
338  break;
339  case both_ends:
340  std::strcpy (which, "BE");
341  break;
342  }
343 
344  // tolerance for ARPACK
345  const double tol = control().tolerance();
346 
347  // if the starting vector is used it has to be in resid
348  std::vector<double> resid(n, 1.);
349 
350  // number of Arnoldi basis vectors specified
351  // in additional_data
352  int ncv = additional_data.number_of_arnoldi_vectors;
353 
354  int ldv = n;
355  std::vector<double> v (ldv*ncv, 0.0);
356 
357  //information to the routines
358  std::vector<int> iparam (11, 0);
359 
360  iparam[0] = 1; // shift strategy
361 
362  // maximum number of iterations
363  iparam[2] = control().max_steps();
364 
370  iparam[6] = mode;
371  std::vector<int> ipntr (14, 0);
372 
373  // work arrays for ARPACK
374  double *workd;
375  workd = new double[3*n];
376 
377  for (unsigned int i=0; i<3*n; ++i)
378  workd[i] = 0.0;
379 
380  int lworkl = 3*ncv*(ncv + 6);
381  std::vector<double> workl (lworkl, 0.);
382  //information out of the iteration
383  int info = 1;
384 
385  while (ido != 99)
386  {
387  // call of ARPACK dnaupd routine
388  dnaupd_(&ido, bmat, &n_inside_arpack, which, &nev, &tol,
389  &resid[0], &ncv, &v[0], &ldv, &iparam[0], &ipntr[0],
390  workd, &workl[0], &lworkl, &info);
391 
392  if (ido == 99)
393  break;
394 
395  switch (mode)
396  {
397  case 3:
398  {
399  switch (ido)
400  {
401  case -1:
402  {
403 
404  VectorType src,dst,tmp;
405  src.reinit(eigenvectors[0]);
406  dst.reinit(src);
407  tmp.reinit(src);
408 
409 
410  for (size_type i=0; i<src.size(); ++i)
411  src(i) = *(workd+ipntr[0]-1+i);
412 
413  // multiplication with mass matrix M
414  mass_matrix.vmult(tmp, src);
415  // solving linear system
416  inverse.vmult(dst,tmp);
417 
418  for (size_type i=0; i<dst.size(); ++i)
419  *(workd+ipntr[1]-1+i) = dst(i);
420  }
421  break;
422 
423  case 1:
424  {
425 
426  VectorType src,dst,tmp, tmp2;
427  src.reinit(eigenvectors[0]);
428  dst.reinit(src);
429  tmp.reinit(src);
430  tmp2.reinit(src);
431 
432  for (size_type i=0; i<src.size(); ++i)
433  {
434  src(i) = *(workd+ipntr[2]-1+i);
435  tmp(i) = *(workd+ipntr[0]-1+i);
436  }
437  // solving linear system
438  inverse.vmult(dst,src);
439 
440  for (size_type i=0; i<dst.size(); ++i)
441  *(workd+ipntr[1]-1+i) = dst(i);
442  }
443  break;
444 
445  case 2:
446  {
447 
448  VectorType src,dst;
449  src.reinit(eigenvectors[0]);
450  dst.reinit(src);
451 
452  for (size_type i=0; i<src.size(); ++i)
453  src(i) = *(workd+ipntr[0]-1+i);
454 
455  // Multiplication with mass matrix M
456  mass_matrix.vmult(dst, src);
457 
458  for (size_type i=0; i<dst.size(); ++i)
459  *(workd+ipntr[1]-1+i) = dst(i);
460 
461  }
462  break;
463 
464  default:
465  Assert (false, ExcArpackIdo(ido));
466  break;
467  }
468  }
469  break;
470  default:
471  Assert (false, ExcArpackMode(mode));
472  break;
473  }
474  }
475 
476  if (info<0)
477  {
478  Assert (false, ExcArpackInfodsaupd(info));
479  }
480  else
481  {
485  int rvec = 1;
486 
487  // which eigenvectors
488  char howmany = 'A';
489 
490  std::vector<int> select (ncv, 1);
491 
492  int ldz = n;
493 
494  std::vector<double> z (ldz*ncv, 0.);
495 
496  double sigmar = 0.0; // real part of the shift
497  double sigmai = 0.0; // imaginary part of the shift
498 
499  int lworkev = 3*ncv;
500  std::vector<double> workev (lworkev, 0.);
501 
502  std::vector<double> eigenvalues_real (nev, 0.);
503  std::vector<double> eigenvalues_im (nev, 0.);
504 
505  // call of ARPACK dneupd routine
506  dneupd_(&rvec, &howmany, &select[0], &eigenvalues_real[0],
507  &eigenvalues_im[0], &z[0], &ldz, &sigmar, &sigmai,
508  &workev[0], bmat, &n_inside_arpack, which, &nev, &tol,
509  &resid[0], &ncv, &v[0], &ldv,
510  &iparam[0], &ipntr[0], workd, &workl[0], &lworkl, &info);
511 
512  if (info == 1)
513  {
514  Assert (false, ExcArpackInfoMaxIt(control().max_steps()));
515  }
516  else if (info == 3)
517  {
518  Assert (false, ExcArpackNoShifts());
519  }
520  else if (info!=0)
521  {
522  Assert (false, ExcArpackInfodneupd(info));
523  }
524 
525 
526  const unsigned int n_eigenvecs = eigenvectors.size();
527  for (size_type i=0; i<n_eigenvecs; ++i)
528  for (unsigned int j=0; j<n; ++j)
529  eigenvectors[i](j) = v[i*n+j];
530 
531  delete[] workd;
532 
533  AssertDimension (eigenvalues.size(), eigenvalues_real.size());
534  AssertDimension (eigenvalues.size(), eigenvalues_im.size());
535 
536  for (size_type i=0; i<eigenvalues.size(); ++i)
537  eigenvalues[i] = std::complex<double> (eigenvalues_real[i],
538  eigenvalues_im[i]);
539  }
540 }
541 
542 
543 inline
545 {
546  return solver_control;
547 }
548 
549 DEAL_II_NAMESPACE_CLOSE
550 
551 
552 #endif
553 #endif
SolverControl & solver_control
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1052
#define AssertIndexRange(index, range)
Definition: exceptions.h:1081
types::global_dof_index size_type
Definition: arpack_solver.h:96
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:542
DeclException2(ExcInvalidNumberofEigenvalues, int, int,<< "Number of wanted eigenvalues "<< arg1<< " is larger that the size of the matrix "<< arg2)
unsigned int global_dof_index
Definition: types.h:88
#define Assert(cond, exc)
Definition: exceptions.h:294
#define DeclExceptionMsg(Exception, defaulttext)
Definition: exceptions.h:533
ArpackSolver(SolverControl &control, const AdditionalData &data=AdditionalData())
double tolerance() const
const AdditionalData additional_data
void solve(const MatrixType1 &A, const MatrixType2 &B, const INVERSE &inverse, std::vector< std::complex< double > > &eigenvalues, std::vector< VectorType > &eigenvectors, const unsigned int n_eigenvalues=0)
unsigned int max_steps() const
SolverControl & control() const