/usr/include/deal.II/lac/slepc_solver.h is in libdeal.ii-dev 8.1.0-4.
This file is owned by root:root, with mode 0o644.
The actual contents of the file can be viewed below.
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736 737 738 739 740 741 742 743 744 745 746 747 748 749 750 751 752 753 754 755 756 757 758 759 760 761 762 763 764 765 766 767 768 769 770 771 772 773 774 775 776 777 778 779 780 781 782 783 784 785 786 787 788 789 790 791 792 793 794 795 796 797 798 799 800 801 802 803 804 805 806 807 808 809 810 811 812 813 814 815 816 817 818 819 820 821 822 823 824 825 826 827 828 829 830 831 832 833 834 835 836 837 838 839 840 841 842 843 844 845 846 847 848 849 850 851 852 853 854 855 856 857 858 859 860 861 862 863 864 865 866 867 868 869 870 871 872 873 874 875 876 | // ---------------------------------------------------------------------
// $Id: slepc_solver.h 31932 2013-12-08 02:15:54Z heister $
//
// Copyright (C) 2009 - 2013 by the deal.II authors
//
// This file is part of the deal.II library.
//
// The deal.II library is free software; you can use it, redistribute
// it, and/or modify it under the terms of the GNU Lesser General
// Public License as published by the Free Software Foundation; either
// version 2.1 of the License, or (at your option) any later version.
// The full text of the license can be found in the file LICENSE at
// the top level of the deal.II distribution.
//
// ---------------------------------------------------------------------
#ifndef __deal2__slepc_solver_h
#define __deal2__slepc_solver_h
#include <deal.II/base/config.h>
#ifdef DEAL_II_WITH_SLEPC
# include <deal.II/base/std_cxx1x/shared_ptr.h>
# include <deal.II/lac/exceptions.h>
# include <deal.II/lac/solver_control.h>
# include <deal.II/lac/petsc_matrix_base.h>
# include <deal.II/lac/slepc_spectral_transformation.h>
# include <petscconf.h>
# include <petscksp.h>
# include <slepceps.h>
DEAL_II_NAMESPACE_OPEN
/**
* Base class for solver classes using the SLEPc solvers which are
* selected based on flags passed to the eigenvalue problem solver
* context. Derived classes set the right flags to set the right
* solver. Note that: the <code>AdditionalData</code> structure is a
* dummy structure that currently exists for backward/forward
* compatibility only.
*
* The SLEPc solvers are intended to be used for solving the
* generalized eigenspectrum problem $(A-\lambda B)x=0$, for $x\neq0$;
* where $A$ is a system matrix, $B$ is a mass matrix, and $\lambda,
* x$ are a set of eigenvalues and eigenvectors respectively. The
* emphasis is on methods and techniques appropriate for problems in
* which the associated matrices are sparse. Most of the methods
* offered by the SLEPc library are projection methods or other
* methods with similar properties; and wrappers are provided to
* interface to SLEPc solvers that handle both of these problem sets.
*
* SLEPcWrappers can be implemented in application codes in the
* following way:
* @code
* SolverControl solver_control (1000, 1e-9);
* SolverArnoldi system (solver_control, mpi_communicator);
* system.solve (A, B, lambda, x, size_of_spectrum);
* @endcode
* for the generalized eigenvalue problem $Ax=B\lambda x$, where the
* variable <code>const unsigned int size_of_spectrum</code> tells
* SLEPc the number of eigenvector/eigenvalue pairs to solve
* for. Additional options and solver parameters can be passed to the
* SLEPc solvers before calling <code>solve()</code>. For example, if
* the matrices of the general eigenspectrum problem are not hermitian
* and the lower eigenvalues are wanted only, the following code can
* be implemented before calling <code>solve()</code>:
* @code
* system.set_problem_type (EPS_NHEP);
* system.set_which_eigenpairs (EPS_SMALLEST_REAL);
* @endcode
* These options can also be set at the commandline.
*
* See also <code>step-36</code> for a hands-on example.
*
* An alternative implementation to the one above is to use the API
* internals directly within the application code. In this way the
* calling sequence requires calling several of SolverBase functions
* rather than just one. This freedom is intended for use of the
* SLEPcWrappers that require a greater handle on the eigenvalue
* problem solver context. See also the API of, for example:
@code
template <typename OutputVector>
void
SolverBase::solve (const PETScWrappers::MatrixBase &A,
const PETScWrappers::MatrixBase &B,
std::vector<PetscScalar> &eigenvalues,
std::vector<OutputVector> &eigenvectors,
const unsigned int n_eigenpairs)
{ ... }
@endcode
* as an example on how to do this.
*
* For further information and explanations on handling the @ref
* SLEPcWrappers "SLEPcWrappers", see also the @ref PETScWrappers
* "PETScWrappers", on which they depend.
*
* @ingroup SLEPcWrappers
*
* @author Toby D. Young 2008, 2009, 2010, 2011, 2013; and Rickard
* Armiento 2008.
*
* @note Various tweaks and enhancments contributed by Eloy Romero and
* Jose E. Roman 2009, 2010.
*/
namespace SLEPcWrappers
{
/**
* Base class for solver classes using the SLEPc solvers. Since
* solvers in SLEPc are selected based on flags passed to a generic
* solver object, basically all the actual solver calls happen in
* this class, and derived classes simply set the right flags to
* select one solver or another, or to set certain parameters for
* individual solvers.
*/
class SolverBase
{
public:
/**
* Constructor. Takes the MPI communicator over which parallel
* computations are to happen.
*/
SolverBase (SolverControl &cn,
const MPI_Comm &mpi_communicator);
/**
* Destructor.
*/
virtual ~SolverBase ();
/**
* Composite method that solves the eigensystem $Ax=\lambda
* x$. The eigenvector sent in has to have at least one element
* that we can use as a template when resizing, since we do not
* know the parameters of the specific vector class used
* (i.e. local_dofs for MPI vectors). However, while copying
* eigenvectors, at least twice the memory size of
* <tt>eigenvectors</tt> is being used (and can be more). To avoid
* doing this, the fairly standard calling sequence executed here
* is used: Initialise; Set up matrices for solving; Actually
* solve the system; Gather the solution(s); and reset.
*
* @note Note that the number of converged eigenvectors can be
* larger than the number of eigenvectors requested; this is due
* to a round off error (success) of the eigenproblem solver
* context. If this is found to be the case we simply do not
* bother with more eigenpairs than requested, but handle that it
* may be more than specified by ignoring any extras. By default
* one eigenvector/eigenvalue pair is computed.
*/
template <typename OutputVector>
void
solve (const PETScWrappers::MatrixBase &A,
std::vector<PetscScalar> &eigenvalues,
std::vector<OutputVector> &eigenvectors,
const unsigned int n_eigenpairs = 1);
/**
* Same as above, but here a composite method for solving the
* system $A x=\lambda B x$, for real matrices, vectors, and
* values $A, B, x, \lambda$.
*/
template <typename OutputVector>
void
solve (const PETScWrappers::MatrixBase &A,
const PETScWrappers::MatrixBase &B,
std::vector<PetscScalar> &eigenvalues,
std::vector<OutputVector> &eigenvectors,
const unsigned int n_eigenpairs = 1);
/**
* Same as above, but here a composite method for solving the
* system $A x=\lambda B x$ with real matrices $A, B$ and
* imaginary eigenpairs $x, \lamda$.
*/
template <typename OutputVector>
void
solve (const PETScWrappers::MatrixBase &A,
const PETScWrappers::MatrixBase &B,
std::vector<double> &real_eigenvalues,
std::vector<double> &imag_eigenvalues,
std::vector<OutputVector> &real_eigenvectors,
std::vector<OutputVector> &imag_eigenvectors,
const unsigned int n_eigenpairs = 1);
/**
* Set the initial vector for the solver.
*/
void
set_initial_vector
(const PETScWrappers::VectorBase &this_initial_vector);
/**
* Set the spectral transformation to be used.
*/
void
set_transformation (SLEPcWrappers::TransformationBase &this_transformation);
/**
* Set target eigenvalues in the spectrum to be computed. By
* default, no target is set.
*/
void
set_target_eigenvalue (const PetscScalar &this_target);
/**
* Indicate which part of the spectrum is to be computed. By
* default largest magnitude eigenvalues are computed.
*
* @note For other allowed values see the SLEPc documentation.
*/
void
set_which_eigenpairs (EPSWhich set_which);
/**
* Specify the type of the eigenspectrum problem. This can be used
* to exploit known symmetries of the matrices that make up the
* standard/generalized eigenspectrum problem. By default a
* non-Hermitian problem is assumed.
*
* @note For other allowed values see the SLEPc documentation.
*/
void
set_problem_type (EPSProblemType set_problem);
/**
* Take the information provided from SLEPc and checks it against
* deal.II's own SolverControl objects to see if convergence has
* been reached.
*/
void
get_solver_state (const SolverControl::State state);
/**
* Exception. Standard exception.
*/
DeclException0 (ExcSLEPcWrappersUsageError);
/**
* Exception. SLEPc error with error number.
*/
DeclException1 (ExcSLEPcError,
int,
<< " An error with error number " << arg1
<< " occurred while calling a SLEPc function");
/**
* Exception. Convergence failure on the number of eigenvectors.
*/
DeclException2 (ExcSLEPcEigenvectorConvergenceMismatchError,
int, int,
<< " The number of converged eigenvectors is " << arg1
<< " but " << arg2 << " were requested. ");
/**
* Access to the object that controls convergence.
*/
SolverControl &control () const;
protected:
/**
* Reference to the object that controls convergence of the
* iterative solver.
*/
SolverControl &solver_control;
/**
* Copy of the MPI communicator object to be used for the solver.
*/
const MPI_Comm mpi_communicator;
/**
* Function that takes an Eigenvalue Problem Solver context
* object, and sets the type of solver that is requested by the
* derived class.
*/
virtual void set_solver_type (EPS &eps) const = 0;
/**
* Reset the solver, and return memory for eigenvectors
*/
void
reset ();
/**
* Retrieve the SLEPc solver object that is internally used.
*/
EPS *get_eps ();
/**
* Solve the linear system for <code>n_eigenpairs</code>
* eigenstates. Parameter <code>n_converged</code> contains the
* actual number of eigenstates that have converged; this can
* be both fewer or more than n_eigenpairs, depending on the
* SLEPc eigensolver used.
*/
void
solve (const unsigned int n_eigenpairs,
unsigned int *n_converged);
/**
* Access the real parts of solutions for a solved eigenvector
* problem, pair index solutions, $\text{index}\,\in\,0\hdots
* \text{n\_converged}-1$.
*/
void
get_eigenpair (const unsigned int index,
PetscScalar &eigenvalues,
PETScWrappers::VectorBase &eigenvectors);
/**
* Access the real and imaginary parts of solutions for a solved
* eigenvector problem, pair index solutions,
* $\text{index}\,\in\,0\hdots \text{n\_converged}-1$.
*/
void
get_eigenpair (const unsigned int index,
double &real_eigenvalues,
double &imag_eigenvalues,
PETScWrappers::VectorBase &real_eigenvectors,
PETScWrappers::VectorBase &imag_eigenvectors);
/**
* Initialize solver for the linear system $Ax=\lambda x$. (Note:
* this is required before calling solve ())
*/
void
set_matrices (const PETScWrappers::MatrixBase &A);
/**
* Same as above, but here initialize solver for the linear system
* $A x=\lambda B x$.
*/
void
set_matrices (const PETScWrappers::MatrixBase &A,
const PETScWrappers::MatrixBase &B);
/**
* Target eigenvalue to solve for.
*/
PetscScalar target_eigenvalue;
/**
* Which portion of the spectrum to solve from.
*/
EPSWhich set_which;
/**
* Set the eigenspectrum problem type.
*/
EPSProblemType set_problem;
private:
/**
* The matrix $A$ of the generalized eigenvalue problem
* $Ax=B\lambda x$, or the standard eigenvalue problem $Ax=\lambda
* x$.
*/
const PETScWrappers::MatrixBase *opA;
/**
* The matrix $B$ of the generalized eigenvalue problem
* $Ax=B\lambda x$.
*/
const PETScWrappers::MatrixBase *opB;
/**
* An initial vector used to "feed" some SLEPc solvers.
*/
const PETScWrappers::VectorBase *initial_vector;
/**
* Pointer to an an object that describes transformations that can
* be applied to the eigenvalue problem.
*/
SLEPcWrappers::TransformationBase *transformation;
/**
* A function that can be used in SLEPc as a callback to check on
* convergence.
*
* @note This function is redundant.
*/
static
int
convergence_test (EPS eps,
PetscScalar real_eigenvalue,
PetscScalar imag_eigenvalue,
PetscReal residual_norm,
PetscReal *estimated_error,
void *solver_control);
/**
* Objects of this type are explicitly created, but are destroyed
* when the surrounding solver object goes out of scope, or when
* we assign a new value to the pointer to this object. The
* respective Destroy functions are therefore written into the
* destructor of this object, even though the object does not have
* a constructor.
*/
struct SolverData
{
/**
* Destructor.
*/
~SolverData ();
/**
* Objects for Eigenvalue Problem Solver.
*/
EPS eps;
/**
* Convergence.
*/
EPSConvergedReason reason;
};
/**
* Pointer to the <code>SolverData</code> object.
*/
std_cxx1x::shared_ptr<SolverData> solver_data;
};
/**
* An implementation of the solver interface using the SLEPc
* Krylov-Schur solver. Usage: All spectrum, all problem types,
* complex.
*
* @ingroup SLEPcWrappers
*
* @author Toby D. Young 2008
*/
class SolverKrylovSchur : public SolverBase
{
public:
/**
* Standardized data struct to pipe additional data to the solver,
* should it be needed.
*/
struct AdditionalData
{};
/**
* SLEPc solvers will want to have an MPI communicator context
* over which computations are parallelized. By default, this
* carries the same behaviour has the PETScWrappers, but you can
* change that.
*/
SolverKrylovSchur (SolverControl &cn,
const MPI_Comm &mpi_communicator = PETSC_COMM_SELF,
const AdditionalData &data = AdditionalData());
protected:
/**
* Store a copy of the flags for this particular solver.
*/
const AdditionalData additional_data;
/**
* Function that takes a Eigenvalue Problem Solver context object,
* and sets the type of solver that is appropriate for this class.
*/
virtual void set_solver_type (EPS &eps) const;
};
/**
* An implementation of the solver interface using the SLEPc Arnoldi
* solver. Usage: All spectrum, all problem types, complex.
*
* @ingroup SLEPcWrappers
*
* @author Toby D. Young 2008, 2011
*/
class SolverArnoldi : public SolverBase
{
public:
/**
* Standardized data struct to pipe additional data to the solver,
* should it be needed.
*/
struct AdditionalData
{
/**
* Constructor. By default, set the option of delayed
* reorthogonalization to false, i.e. don't do it.
*/
AdditionalData (const bool delayed_reorthogonalization = false);
/**
* Flag for delayed reorthogonalization.
*/
bool delayed_reorthogonalization;
};
/**
* SLEPc solvers will want to have an MPI communicator context
* over which computations are parallelized. By default, this
* carries the same behaviour has the PETScWrappers, but you can
* change that.
*/
SolverArnoldi (SolverControl &cn,
const MPI_Comm &mpi_communicator = PETSC_COMM_SELF,
const AdditionalData &data = AdditionalData());
protected:
/**
* Store a copy of the flags for this particular solver.
*/
const AdditionalData additional_data;
/**
* Function that takes a Eigenvalue Problem Solver context object,
* and sets the type of solver that is appropriate for this class.
*/
virtual void set_solver_type (EPS &eps) const;
};
/**
* An implementation of the solver interface using the SLEPc Lanczos
* solver. Usage: All spectrum, all problem types, complex.
*
* @ingroup SLEPcWrappers
*
* @author Toby D. Young 2009
*/
class SolverLanczos : public SolverBase
{
public:
/**
* Standardized data struct to pipe additional data to the solver,
* should it be needed.
*/
struct AdditionalData
{};
/**
* SLEPc solvers will want to have an MPI communicator context
* over which computations are parallelized. By default, this
* carries the same behaviour has the PETScWrappers, but you can
* change that.
*/
SolverLanczos (SolverControl &cn,
const MPI_Comm &mpi_communicator = PETSC_COMM_SELF,
const AdditionalData &data = AdditionalData());
protected:
/**
* Store a copy of the flags for this particular solver.
*/
const AdditionalData additional_data;
/**
* Function that takes a Eigenvalue Problem Solver context object,
* and sets the type of solver that is appropriate for this class.
*/
virtual void set_solver_type (EPS &eps) const;
};
/**
* An implementation of the solver interface using the SLEPc Power
* solver. Usage: Largest values of spectrum only, all problem
* types, complex.
*
* @ingroup SLEPcWrappers
*
* @author Toby D. Young 2010
*/
class SolverPower : public SolverBase
{
public:
/**
* Standardized data struct to pipe additional data to the solver,
* should it be needed.
*/
struct AdditionalData
{};
/**
* SLEPc solvers will want to have an MPI communicator context
* over which computations are parallelized. By default, this
* carries the same behaviour has the PETScWrappers, but you can
* change that.
*/
SolverPower (SolverControl &cn,
const MPI_Comm &mpi_communicator = PETSC_COMM_SELF,
const AdditionalData &data = AdditionalData());
protected:
/**
* Store a copy of the flags for this particular solver.
*/
const AdditionalData additional_data;
/**
* Function that takes a Eigenvalue Problem Solver context object,
* and sets the type of solver that is appropriate for this class.
*/
virtual void set_solver_type (EPS &eps) const;
};
/**
* An implementation of the solver interface using the SLEPc
* Davidson solver. Usage (incomplete/untested): All problem types.
*
* @ingroup SLEPcWrappers
*
* @author Toby D. Young 2010
*/
class SolverGeneralizedDavidson : public SolverBase
{
public:
/**
* Standardized data struct to pipe additional data to the solver,
* should it be needed.
*/
struct AdditionalData
{};
/**
* SLEPc solvers will want to have an MPI communicator context
* over which computations are parallelized. By default, this
* carries the same behaviour has the PETScWrappers, but you can
* change that.
*/
SolverGeneralizedDavidson (SolverControl &cn,
const MPI_Comm &mpi_communicator = PETSC_COMM_SELF,
const AdditionalData &data = AdditionalData());
protected:
/**
* Store a copy of the flags for this particular solver.
*/
const AdditionalData additional_data;
/**
* Function that takes a Eigenvalue Problem Solver context object,
* and sets the type of solver that is appropriate for this class.
*/
virtual void set_solver_type (EPS &eps) const;
};
/**
* An implementation of the solver interface using the SLEPc
* Jacobi-Davidson solver. Usage: All problem types.
*
* @ingroup SLEPcWrappers
*
* @author Toby D. Young 2013
*/
class SolverJacobiDavidson : public SolverBase
{
public:
/**
* Standardized data struct to pipe additional data to the solver,
* should it be needed.
*/
struct AdditionalData
{};
/**
* SLEPc solvers will want to have an MPI communicator context
* over which computations are parallelized. By default, this
* carries the same behaviour has the PETScWrappers, but you can
* change that.
*/
SolverJacobiDavidson (SolverControl &cn,
const MPI_Comm &mpi_communicator = PETSC_COMM_SELF,
const AdditionalData &data = AdditionalData());
protected:
/**
* Store a copy of the flags for this particular solver.
*/
const AdditionalData additional_data;
/**
* Function that takes a Eigenvalue Problem Solver context object,
* and sets the type of solver that is appropriate for this class.
*/
virtual void set_solver_type (EPS &eps) const;
};
/**
* An implementation of the solver interface using the SLEPc LAPACK
* direct solver.
*
* @ingroup SLEPcWrappers
*
* @author Toby D. Young 2013
*/
class SolverLAPACK : public SolverBase
{
public:
/**
* Standardized data struct to pipe additional data to the solver,
* should it be needed.
*/
struct AdditionalData
{};
/**
* SLEPc solvers will want to have an MPI communicator context
* over which computations are parallelized. By default, this
* carries the same behaviour has the PETScWrappers, but you can
* change that.
*/
SolverLAPACK (SolverControl &cn,
const MPI_Comm &mpi_communicator = PETSC_COMM_SELF,
const AdditionalData &data = AdditionalData());
protected:
/**
* Store a copy of the flags for this particular solver.
*/
const AdditionalData additional_data;
/**
* Function that takes a Eigenvalue Problem Solver context object,
* and sets the type of solver that is appropriate for this class.
*/
virtual void set_solver_type (EPS &eps) const;
};
// --------------------------- inline and template functions -----------
/**
* This is declared here to make it possible to take a std::vector
* of different PETScWrappers vector types
*/
// todo: The logic of these functions can be simplified without breaking backward compatibility...
template <typename OutputVector>
void
SolverBase::solve (const PETScWrappers::MatrixBase &A,
std::vector<PetscScalar> &eigenvalues,
std::vector<OutputVector> &eigenvectors,
const unsigned int n_eigenpairs)
{
// Panic if the number of eigenpairs wanted is out of bounds.
AssertThrow ((n_eigenpairs > 0) && (n_eigenpairs <= A.m ()),
ExcSLEPcWrappersUsageError());
// Set the matrices of the problem
set_matrices (A);
// and solve
unsigned int n_converged = 0;
solve (n_eigenpairs, &n_converged);
if (n_converged > n_eigenpairs)
n_converged = n_eigenpairs;
AssertThrow (n_converged == n_eigenpairs,
ExcSLEPcEigenvectorConvergenceMismatchError(n_converged, n_eigenpairs));
AssertThrow (eigenvectors.size() != 0, ExcSLEPcWrappersUsageError());
eigenvectors.resize (n_converged, eigenvectors.front());
eigenvalues.resize (n_converged);
for (unsigned int index=0; index<n_converged; ++index)
get_eigenpair (index, eigenvalues[index], eigenvectors[index]);
}
template <typename OutputVector>
void
SolverBase::solve (const PETScWrappers::MatrixBase &A,
const PETScWrappers::MatrixBase &B,
std::vector<PetscScalar> &eigenvalues,
std::vector<OutputVector> &eigenvectors,
const unsigned int n_eigenpairs)
{
// Guard against incompatible matrix sizes:
AssertThrow (A.m() == B.m (), ExcDimensionMismatch(A.m(), B.m()));
AssertThrow (A.n() == B.n (), ExcDimensionMismatch(A.n(), B.n()));
// Panic if the number of eigenpairs wanted is out of bounds.
AssertThrow ((n_eigenpairs>0) && (n_eigenpairs<=A.m ()),
ExcSLEPcWrappersUsageError());
// Set the matrices of the problem
set_matrices (A, B);
// and solve
unsigned int n_converged = 0;
solve (n_eigenpairs, &n_converged);
if (n_converged>=n_eigenpairs)
n_converged = n_eigenpairs;
AssertThrow (n_converged==n_eigenpairs,
ExcSLEPcEigenvectorConvergenceMismatchError(n_converged, n_eigenpairs));
AssertThrow (eigenvectors.size() != 0, ExcSLEPcWrappersUsageError());
eigenvectors.resize (n_converged, eigenvectors.front());
eigenvalues.resize (n_converged);
for (unsigned int index=0; index<n_converged; ++index)
get_eigenpair (index, eigenvalues[index], eigenvectors[index]);
}
template <typename OutputVector>
void
SolverBase::solve (const PETScWrappers::MatrixBase &A,
const PETScWrappers::MatrixBase &B,
std::vector<double> &real_eigenvalues,
std::vector<double> &imag_eigenvalues,
std::vector<OutputVector> &real_eigenvectors,
std::vector<OutputVector> &imag_eigenvectors,
const unsigned int n_eigenpairs)
{
// Guard against incompatible matrix sizes:
AssertThrow (A.m() == B.m (), ExcDimensionMismatch(A.m(), B.m()));
AssertThrow (A.n() == B.n (), ExcDimensionMismatch(A.n(), B.n()));
// and incompatible eigenvalue/eigenvector sizes
AssertThrow (real_eigenvalues.size() == imag_eigenvalues.size(),
ExcDimensionMismatch(real_eigenvalues.size(), imag_eigenvalues.size()));
AssertThrow (real_eigenvectors.size() == imag_eigenvectors.size (),
ExcDimensionMismatch(real_eigenvectors.size(), imag_eigenvectors.size()));
// Panic if the number of eigenpairs wanted is out of bounds.
AssertThrow ((n_eigenpairs>0) && (n_eigenpairs<=A.m ()),
ExcSLEPcWrappersUsageError());
// Set the matrices of the problem
set_matrices (A, B);
// and solve
unsigned int n_converged = 0;
solve (n_eigenpairs, &n_converged);
if (n_converged>=n_eigenpairs)
n_converged = n_eigenpairs;
AssertThrow (n_converged==n_eigenpairs,
ExcSLEPcEigenvectorConvergenceMismatchError(n_converged, n_eigenpairs));
AssertThrow ((real_eigenvectors.size()!=0) && (imag_eigenvectors.size()!=0),
ExcSLEPcWrappersUsageError());
real_eigenvectors.resize (n_converged, real_eigenvectors.front());
imag_eigenvectors.resize (n_converged, imag_eigenvectors.front());
real_eigenvalues.resize (n_converged);
imag_eigenvalues.resize (n_converged);
for (unsigned int index=0; index<n_converged; ++index)
get_eigenpair (index,
real_eigenvalues[index], imag_eigenvalues[index],
real_eigenvectors[index], imag_eigenvectors[index]);
}
}
DEAL_II_NAMESPACE_CLOSE
#endif // DEAL_II_WITH_SLEPC
/*---------------------------- slepc_solver.h ---------------------------*/
#endif
/*---------------------------- slepc_solver.h ---------------------------*/
|