deal.II version GIT relicensing-2165-gc91f007519 2024-11-20 01:40:00+00:00
|
#include <deal.II/lac/parpack_solver.h>
Classes | |
struct | AdditionalData |
Public Types | |
enum | WhichEigenvalues { algebraically_largest , algebraically_smallest , largest_magnitude , smallest_magnitude , largest_real_part , smallest_real_part , largest_imaginary_part , smallest_imaginary_part , both_ends } |
using | size_type = types::global_dof_index |
Public Member Functions | |
SolverControl & | control () const |
PArpackSolver (SolverControl &control, const MPI_Comm mpi_communicator, const AdditionalData &data=AdditionalData()) | |
void | reinit (const IndexSet &locally_owned_dofs) |
void | reinit (const IndexSet &locally_owned_dofs, const std::vector< IndexSet > &partitioning) |
void | reinit (const VectorType &distributed_vector) |
void | set_initial_vector (const VectorType &vec) |
void | set_shift (const std::complex< double > sigma) |
template<typename MatrixType1 , typename MatrixType2 , typename INVERSE > | |
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) |
template<typename MatrixType1 , typename MatrixType2 , typename INVERSE > | |
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) |
std::size_t | memory_consumption () const |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
EnableObserverPointer functionality | |
Classes derived from EnableObserverPointer provide a facility to subscribe to this object. This is mostly used by the ObserverPointer class. | |
void | subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
void | unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
unsigned int | n_subscriptions () const |
template<typename StreamType > | |
void | list_subscribers (StreamType &stream) const |
void | list_subscribers () const |
Static Public Member Functions | |
static ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
Protected Attributes | |
SolverControl & | solver_control |
const AdditionalData | additional_data |
MPI_Comm | mpi_communicator |
MPI_Fint | mpi_communicator_fortran |
int | lworkl |
std::vector< double > | workl |
std::vector< double > | workd |
int | nloc |
int | ncv |
int | ldv |
std::vector< double > | v |
bool | initial_vector_provided |
std::vector< double > | resid |
int | ldz |
std::vector< double > | z |
int | lworkev |
std::vector< double > | workev |
std::vector< int > | select |
VectorType | src |
VectorType | dst |
VectorType | tmp |
std::vector< types::global_dof_index > | local_indices |
double | sigmar |
double | sigmai |
Private Types | |
using | map_value_type = decltype(counter_map)::value_type |
using | map_iterator = decltype(counter_map)::iterator |
Private Member Functions | |
void | internal_reinit (const IndexSet &locally_owned_dofs) |
void | check_no_subscribers () const noexcept |
Static Private Member Functions | |
static ::ExceptionBase & | PArpackExcConvergedEigenvectors (int arg1, int arg2) |
static ::ExceptionBase & | PArpackExcInvalidNumberofEigenvalues (int arg1, int arg2) |
static ::ExceptionBase & | PArpackExcInvalidEigenvectorSize (int arg1, int arg2) |
static ::ExceptionBase & | PArpackExcInvalidEigenvectorSizeNonsymmetric (int arg1, int arg2) |
static ::ExceptionBase & | PArpackExcInvalidEigenvalueSize (int arg1, int arg2) |
static ::ExceptionBase & | PArpackExcInvalidNumberofArnoldiVectors (int arg1, int arg2) |
static ::ExceptionBase & | PArpackExcSmallNumberofArnoldiVectors (int arg1, int arg2) |
static ::ExceptionBase & | PArpackExcIdo (int arg1) |
static ::ExceptionBase & | PArpackExcMode (int arg1) |
static ::ExceptionBase & | PArpackExcInfoPdnaupd (int arg1) |
static ::ExceptionBase & | PArpackExcInfoPdneupd (int arg1) |
static ::ExceptionBase & | PArpackExcInfoMaxIt (int arg1) |
static ::ExceptionBase & | PArpackExcNoShifts (int arg1) |
Private Attributes | |
std::atomic< unsigned int > | counter |
std::map< std::string, unsigned int > | counter_map |
std::vector< std::atomic< bool > * > | validity_pointers |
const std::type_info * | object_info |
Static Private Attributes | |
static std::mutex | mutex |
Interface for using PARPACK. PARPACK is a collection of Fortran77 subroutines designed to solve large scale eigenvalue problems. Here we interface to the routines pdneupd
, pdseupd
, pdnaupd
, pdsaupd
of PARPACK. The package is designed to compute a few eigenvalues and corresponding eigenvectors of a general n by n matrix A. It is most appropriate for large sparse matrices A.
In this class we make use of the method applied to 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 ArpackSolver can be used in application codes in the following way:
for the generalized eigenvalue problem \(Ax=B\lambda x\), where the variable size_of_spectrum
tells PARPACK the number of eigenvector/eigenvalue pairs to solve for. Here, lambda
is a vector that will contain the eigenvalues computed, x
a vector of objects of type V
that will contain the eigenvectors computed.
Currently, only three modes of (P)Arpack are implemented. In mode 3 (default), OP
is an inverse operation for the matrix A - sigma * B
, where sigma
is a shift value, set to zero by default. Whereas in mode 2, OP
is an inverse of M
. Finally, mode 1 corresponds to standard eigenvalue problem without spectral transformation \(Ax=\lambda x\). The mode can be specified via AdditionalData object. Note that for shift-and-invert (mode=3), the sought eigenpairs are those after the spectral transformation is applied.
The OP
can be specified by using a LinearOperator:
The class is intended to be used with MPI and can work on arbitrary vector and matrix distributed classes. Both symmetric and non-symmetric A
are supported.
For further information on how the PARPACK routines pdneupd
, pdseupd
, pdnaupd
, pdsaupd
work and also how to set the parameters appropriately please take a look into the PARPACK manual.
Definition at line 210 of file parpack_solver.h.
using PArpackSolver< VectorType >::size_type = types::global_dof_index |
Declare the type for container size.
Definition at line 216 of file parpack_solver.h.
|
privateinherited |
The data type used in counter_map.
Definition at line 238 of file enable_observer_pointer.h.
|
privateinherited |
The iterator type used in counter_map.
Definition at line 243 of file enable_observer_pointer.h.
enum PArpackSolver::WhichEigenvalues |
An enum that lists the possible choices for which eigenvalues to compute in the solve() function. Note, that this corresponds to the problem after shift-and-invert (the only currently supported spectral transformation) is applied.
A particular choice is limited based on symmetric or non-symmetric matrix A
considered.
Definition at line 227 of file parpack_solver.h.
PArpackSolver< VectorType >::PArpackSolver | ( | SolverControl & | control, |
const MPI_Comm | mpi_communicator, | ||
const AdditionalData & | data = AdditionalData() |
||
) |
Constructor.
Definition at line 639 of file parpack_solver.h.
SolverControl & PArpackSolver< VectorType >::control | ( | ) | const |
Access to the object that controls convergence.
Definition at line 1159 of file parpack_solver.h.
void PArpackSolver< VectorType >::reinit | ( | const IndexSet & | locally_owned_dofs | ) |
Initialize internal variables.
Definition at line 725 of file parpack_solver.h.
void PArpackSolver< VectorType >::reinit | ( | const IndexSet & | locally_owned_dofs, |
const std::vector< IndexSet > & | partitioning | ||
) |
Initialize internal variables when working with BlockVectors. locally_owned_dofs
is used to set the dimension of the problem, whereas partitioning
is used for calling the reinit of the deal.II blockvectors used.
Definition at line 753 of file parpack_solver.h.
void PArpackSolver< VectorType >::reinit | ( | const VectorType & | distributed_vector | ) |
Initialize internal variables from the input distributed_vector
.
Definition at line 739 of file parpack_solver.h.
void PArpackSolver< VectorType >::set_initial_vector | ( | const VectorType & | vec | ) |
Set initial vector for building Krylov space.
Definition at line 671 of file parpack_solver.h.
void PArpackSolver< VectorType >::set_shift | ( | const std::complex< double > | sigma | ) |
Set shift sigma
for shift-and-invert spectral transformation.
If this function is not called, the shift is assumed to be zero.
mode=3
(see the general documentation of this class for a definition of what the different modes are). Definition at line 661 of file parpack_solver.h.
void PArpackSolver< VectorType >::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 | ||
) |
Solve the generalized eigensprectrum problem \(A x=\lambda B x\) by calling the pd(n/s)eupd
and pd(n/s)aupd
functions of PARPACK.
In mode=3
, inverse
should correspond to \([A-\sigma B]^{-1}\), whereas in mode=2
it should represent \(B^{-1}\). For mode=1
both B
and inverse
are ignored.
Definition at line 769 of file parpack_solver.h.
void PArpackSolver< VectorType >::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 | ||
) |
Same as above but takes eigenvectors as pointers.
Definition at line 787 of file parpack_solver.h.
std::size_t PArpackSolver< VectorType >::memory_consumption | ( | ) | const |
Return the memory consumption of this class in bytes.
Definition at line 588 of file parpack_solver.h.
|
private |
Initialize internal variables which depend on locally_owned_dofs
.
This function is called inside the reinit() functions
Definition at line 685 of file parpack_solver.h.
|
inherited |
Subscribes a user of the object by storing the pointer validity
. The subscriber may be identified by text supplied as identifier
.
Definition at line 131 of file enable_observer_pointer.cc.
|
inherited |
Unsubscribes a user from the object.
identifier
and the validity
pointer must be the same as the one supplied to subscribe(). Definition at line 151 of file enable_observer_pointer.cc.
|
inlineinherited |
Return the present number of subscriptions to this object. This allows to use this class for reference counted lifetime determination where the last one to unsubscribe also deletes the object.
Definition at line 322 of file enable_observer_pointer.h.
|
inlineinherited |
List the subscribers to the input stream
.
Definition at line 339 of file enable_observer_pointer.h.
|
inherited |
List the subscribers to deallog
.
Definition at line 199 of file enable_observer_pointer.cc.
|
inlineinherited |
Read or write the data of this object to or from a stream for the purpose of serialization using the BOOST serialization library.
This function does not actually serialize any of the member variables of this class. The reason is that what this class stores is only who subscribes to this object, but who does so at the time of storing the contents of this object does not necessarily have anything to do with who subscribes to the object when it is restored. Consequently, we do not want to overwrite the subscribers at the time of restoring, and then there is no reason to write the subscribers out in the first place.
Definition at line 331 of file enable_observer_pointer.h.
|
privatenoexceptinherited |
Check that there are no objects subscribing to this object. If this check passes then it is safe to destroy the current object. It this check fails then this function will either abort or print an error message to deallog (by using the AssertNothrow mechanism), but will not throw an exception.
Definition at line 53 of file enable_observer_pointer.cc.
|
protected |
Reference to the object that controls convergence of the iterative solver.
Definition at line 380 of file parpack_solver.h.
|
protected |
Store a copy of the flags for this particular solver.
Definition at line 385 of file parpack_solver.h.
|
protected |
C++ MPI communicator.
Definition at line 392 of file parpack_solver.h.
|
protected |
Fortran MPI communicator.
Definition at line 397 of file parpack_solver.h.
|
protected |
Length of the work array workl.
Definition at line 404 of file parpack_solver.h.
|
protected |
Double precision work array of length lworkl
Definition at line 409 of file parpack_solver.h.
|
protected |
Double precision work array of length 3*N
Definition at line 414 of file parpack_solver.h.
|
protected |
Number of local degrees of freedom.
Definition at line 419 of file parpack_solver.h.
|
protected |
Number of Arnoldi basis vectors specified in additional_data
Definition at line 424 of file parpack_solver.h.
|
protected |
The leading dimension of the array v
Definition at line 430 of file parpack_solver.h.
|
protected |
Double precision vector of size ldv by NCV. Will contains the final set of Arnoldi basis vectors.
Definition at line 436 of file parpack_solver.h.
|
protected |
An auxiliary flag which is set to true when initial vector is provided.
Definition at line 441 of file parpack_solver.h.
|
protected |
The initial residual vector, possibly from a previous run. On output, it contains the final residual vector.
Definition at line 447 of file parpack_solver.h.
|
protected |
The leading dimension of the array Z equal to nloc.
Definition at line 452 of file parpack_solver.h.
|
protected |
A vector of minimum size of nloc by NEV+1. Z contains the B-orthonormal Ritz vectors of the eigensystem A*z = lambda*B*z corresponding to the Ritz value approximations.
Definition at line 459 of file parpack_solver.h.
|
protected |
The size of the workev array.
Definition at line 464 of file parpack_solver.h.
|
protected |
Double precision work array of dimension 3* NCV.
Definition at line 469 of file parpack_solver.h.
|
protected |
A vector of dimension NCV.
Definition at line 474 of file parpack_solver.h.
|
protected |
Temporary vectors used between Arpack and deal.II
Definition at line 479 of file parpack_solver.h.
|
protected |
Definition at line 479 of file parpack_solver.h.
|
protected |
Definition at line 479 of file parpack_solver.h.
|
protected |
Indices of local degrees of freedom.
Definition at line 484 of file parpack_solver.h.
|
protected |
Real part of the shift
Definition at line 489 of file parpack_solver.h.
|
protected |
Imaginary part of the shift
Definition at line 494 of file parpack_solver.h.
|
mutableprivateinherited |
Store the number of objects which subscribed to this object. Initially, this number is zero, and upon destruction it shall be zero again (i.e. all objects which subscribed should have unsubscribed again).
The creator (and owner) of an object is counted in the map below if HE manages to supply identification.
We use the mutable
keyword in order to allow subscription to constant objects also.
This counter may be read from and written to concurrently in multithreaded code: hence we use the std::atomic
class template.
Definition at line 227 of file enable_observer_pointer.h.
|
mutableprivateinherited |
In this map, we count subscriptions for each different identification string supplied to subscribe().
Definition at line 233 of file enable_observer_pointer.h.
|
mutableprivateinherited |
In this vector, we store pointers to the validity bool in the ObserverPointer objects that subscribe to this class.
Definition at line 249 of file enable_observer_pointer.h.
|
mutableprivateinherited |
Pointer to the typeinfo object of this object, from which we can later deduce the class name. Since this information on the derived class is neither available in the destructor, nor in the constructor, we obtain it in between and store it here.
Definition at line 257 of file enable_observer_pointer.h.
|
staticprivateinherited |
A mutex used to ensure data consistency when accessing the mutable
members of this class. This lock is used in the subscribe() and unsubscribe() functions, as well as in list_subscribers()
.
Definition at line 280 of file enable_observer_pointer.h.