/usr/include/deal.II/lac/slepc_solver.h is in libdeal.ii-dev 6.3.1-1.1.
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 | //---------------------------------------------------------------------------
// $Id: slepc_solver.h 21132 2010-05-15 14:17:22Z kronbichler $
// Version: $Name$
// Author: Toby D. Young, Polish Academy of Sciences, 2008, 2009
//
// Copyright (C) 2009, 2010 by the deal.II authors
//
// This file is subject to QPL and may not be distributed
// without copyright and license information. Please refer
// to the file deal.II/doc/license.html for the text and
// further information on this license.
//
//---------------------------------------------------------------------------
#ifndef __deal2__slepc_solver_h
#define __deal2__slepc_solver_h
#include <base/config.h>
#ifdef DEAL_II_USE_SLEPC
# include <base/std_cxx1x/shared_ptr.h>
# include <lac/exceptions.h>
# include <lac/solver_control.h>
# include <lac/petsc_matrix_base.h>
# include <lac/slepc_spectral_transformation.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. On the other hand, 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 M)x=0$, for $x\neq0$;
* where $A$ is a system matrix, $M$ 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 and therefore, most of the
* methods offered by the SLEPc library are projection methods or
* other methods with similar properties. SLEPc implements these basic
* methods as well as more sophisticated algorithms. On the other
* hand, SLEPc is a general library in the sense that it covers
* standard and generalized eigenvalue problems, 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:
@verbatim
SolverControl solver_control (1000, 1e-9);
SolverArnoldi system (solver_control,
mpi_communicator);
system.solve (A, B, eigenvalues, eigenvectors,
size_of_spectrum);
@endverbatim
* 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: See
* also step-36 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:
@verbatim
template <typename OutputVector>
void
SolverBase::solve (const PETScWrappers::MatrixBase &A,
const PETScWrappers::MatrixBase &B,
std::vector<double> &kr,
std::vector<OutputVector> &vr,
const unsigned int n_eigenvectors)
{ ... }
@endverbatim
* 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; and Rickard Armiento 2008
*
* Various tweaks to the SLEPcWrappers have been contributed by Eloy
* Romeoro and Jose Roman.
*/
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>vr</tt> is
* being used (and can be more). To
* avoid doing this, the fairly
* standard calling sequence
* excecuted here is used:
* Initialise; Set up matrices for
* solving; Actually solve the
* system; Gather the solution(s);
* and reset.
*
* 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<double> &kr,
std::vector<OutputVector> &vr,
const unsigned int n_eigenvectors);
/**
* Same as above, but here a
* composite method for solving the
* system $A x=\lambda B x$.
*/
template <typename OutputVector>
void
solve (const PETScWrappers::MatrixBase &A,
const PETScWrappers::MatrixBase &B,
std::vector<double> &kr,
std::vector<OutputVector> &vr,
const unsigned int n_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);
/**
* Set the initial vector for the
* solver.
*/
void
set_initial_vector (const PETScWrappers::VectorBase &initial_vec);
/**
* Set the spectral transformation
* to be used.
*/
void
set_transformation (SLEPcWrappers::TransformationBase &trans);
/**
* Indicate which part of the
* spectrum is to be computed. By
* default largest magnitude
* eigenvalues are computed. For
* other allowed values see the
* SLEPc documentation.
*/
void
set_which_eigenpairs (EPSWhich set_which);
/**
* Solve the linear system for
* n_eigenvectors
* eigenstates. Parameter
* <code>n_converged</code>
* contains the actual number of
* eigenstates that have .
* converged; this can be both
* fewer or more than
* n_eigenvectors, depending on the
* SLEPc eigensolver used.
*/
void
solve (const unsigned int n_eigenvectors, unsigned int *n_converged);
/**
* Access the 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,
#ifndef DEAL_II_USE_PETSC_COMPLEX
double &kr,
#else
std::complex<double> &kr,
#endif
PETScWrappers::VectorBase &vr);
/**
* Reset the solver, and return
* memory for eigenvectors
*/
void
reset();
/**
* Retrieve the SLEPc solver object
* that is internally used.
*/
EPS *
get_EPS ();
/**
* Access to the object that
* controls convergence.
*/
SolverControl &control () const;
/**
* Exceptions.
*/
DeclException0 (ExcSLEPcWrappersUsageError);
DeclException1 (ExcSLEPcError,
int,
<< " An error with error number " << arg1
<< " occured while calling a SLEPc function");
DeclException2 (ExcSLEPcEigenvectorConvergenceMismatchError,
int, int,
<< " The number of converged eigenvectors is " << arg1
<< " but " << arg2 << " were requested. ");
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;
/**
* Attributes that store the
* relevant information for the
* eigenproblem solver context.
*/
EPSWhich set_which;
const PETScWrappers::MatrixBase *opA;
const PETScWrappers::MatrixBase *opB;
const PETScWrappers::VectorBase *ini_vec;
SLEPcWrappers::TransformationBase *transform;
private:
/**
* A function that is used in SLEPc
* as a callback to check on
* convergence. It takes the
* information provided from SLEPc
* and checks it against deal.II's
* own SolverControl objects to see
* if convergence has been reached.
*/
static
int
convergence_test (EPS eps,
#ifdef PETSC_USE_64BIT_INDICES
const PetscInt iteration,
#else
const int iteration,
#endif
const PetscReal residual_norm,
EPSConvergedReason *reason,
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;
};
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_WORLD,
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
*/
class SolverArnoldi : 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.
*/
SolverArnoldi (SolverControl &cn,
const MPI_Comm &mpi_communicator = PETSC_COMM_WORLD,
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_WORLD,
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
*/
template <typename OutputVector>
void
SolverBase::solve (const PETScWrappers::MatrixBase &A,
std::vector<double> &kr,
std::vector<OutputVector> &vr,
const unsigned int n_eigenvectors = 1)
{
unsigned int n_converged;
set_matrices (A);
solve (n_eigenvectors,&n_converged);
if (n_converged > n_eigenvectors)
n_converged = n_eigenvectors;
AssertThrow (n_converged == n_eigenvectors,
ExcSLEPcEigenvectorConvergenceMismatchError(n_converged, n_eigenvectors));
AssertThrow (vr.size() >= 1, ExcSLEPcWrappersUsageError());
vr.resize (n_converged, vr.front());
kr.resize (n_converged);
for (unsigned int index=0; index < n_converged;
+index)
get_eigenpair (index, kr[index], vr[index]);
}
template <typename OutputVector>
void
SolverBase::solve (const PETScWrappers::MatrixBase &A,
const PETScWrappers::MatrixBase &B,
std::vector<double> &kr,
std::vector<OutputVector> &vr,
const unsigned int n_eigenvectors = 1)
{
unsigned int n_converged;
set_matrices (A,B);
solve (n_eigenvectors, &n_converged);
if (n_converged >= n_eigenvectors)
n_converged = n_eigenvectors;
AssertThrow (n_converged == n_eigenvectors,
ExcSLEPcEigenvectorConvergenceMismatchError(n_converged, n_eigenvectors));
AssertThrow (vr.size() >= 1, ExcSLEPcWrappersUsageError());
vr.resize (n_converged, vr.front());
kr.resize (n_converged);
for (unsigned int index=0; index < n_converged;
++index)
get_eigenpair (index, kr[index], vr[index]);
}
}
DEAL_II_NAMESPACE_CLOSE
#endif // DEAL_II_USE_SLEPC
/*---------------------------- slepc_solver.h ---------------------------*/
#endif
/*---------------------------- slepc_solver.h ---------------------------*/
|