oomph::ANASAZI Class Reference

Class for the Anasazi eigensolver. More...

#include <trilinos_eigen_solver.h>

+ Inheritance diagram for oomph::ANASAZI:

Public Member Functions

 ANASAZI ()
 Constructor. More...
 
 ANASAZI (const ANASAZI &)
 Empty copy constructor. More...
 
virtual ~ANASAZI ()
 Destructor, delete the linear solver. More...
 
intnarnoldi ()
 Access function for the number of Arnoldi vectors. More...
 
const intnarnoldi () const
 Access function for the number of Arnoldi vectors (const version) More...
 
void enable_compute_eigenvectors ()
 Set to enable the computation of the eigenvectors (default) More...
 
void disable_compute_eigenvectors ()
 Set to disable the computation of the eigenvectors. More...
 
void solve_eigenproblem (Problem *const &problem_pt, const int &n_eval, Vector< std::complex< double >> &eigenvalue, Vector< DoubleVector > &eigenvector)
 Solve the eigen problem. More...
 
void get_eigenvalues_left_of_shift ()
 Set the desired eigenvalues to be left of the shift. More...
 
void get_eigenvalues_right_of_shift ()
 Set the desired eigenvalues to be right of the shift. More...
 
void track_eigenvalue_real_part ()
 Set the real part to be the quantity of interest (default) More...
 
void track_eigenvalue_imaginary_part ()
 Set the imaginary part fo the quantity of interest. More...
 
void track_eigenvalue_magnitude ()
 Set the magnitude to be the quantity of interest. More...
 
LinearSolver *& linear_solver_pt ()
 Return a pointer to the linear solver object. More...
 
LinearSolver *const & linear_solver_pt () const
 Return a pointer to the linear solver object (const version) More...
 
- Public Member Functions inherited from oomph::EigenSolver
 EigenSolver ()
 Empty constructor. More...
 
 EigenSolver (const EigenSolver &)
 Empty copy constructor. More...
 
virtual ~EigenSolver ()
 Empty destructor. More...
 
void set_shift (const double &shift_value)
 Set the value of the shift. More...
 
const doubleget_shift () const
 Return the value of the shift (const version) More...
 
- Public Member Functions inherited from oomph::DistributableLinearAlgebraObject
 DistributableLinearAlgebraObject ()
 Default constructor - create a distribution. More...
 
 DistributableLinearAlgebraObject (const DistributableLinearAlgebraObject &matrix)=delete
 Broken copy constructor. More...
 
void operator= (const DistributableLinearAlgebraObject &)=delete
 Broken assignment operator. More...
 
virtual ~DistributableLinearAlgebraObject ()
 Destructor. More...
 
LinearAlgebraDistributiondistribution_pt () const
 access to the LinearAlgebraDistribution More...
 
unsigned nrow () const
 access function to the number of global rows. More...
 
unsigned nrow_local () const
 access function for the num of local rows on this processor. More...
 
unsigned nrow_local (const unsigned &p) const
 access function for the num of local rows on this processor. More...
 
unsigned first_row () const
 access function for the first row on this processor More...
 
unsigned first_row (const unsigned &p) const
 access function for the first row on this processor More...
 
bool distributed () const
 distribution is serial or distributed More...
 
bool distribution_built () const
 
void build_distribution (const LinearAlgebraDistribution *const dist_pt)
 
void build_distribution (const LinearAlgebraDistribution &dist)
 

Private Types

typedef double ST
 
typedef Teuchos::ScalarTraits< STSCT
 
typedef SCT::magnitudeType MT
 
typedef Anasazi::MultiVec< STMV
 
typedef Anasazi::Operator< STOP
 
typedef Anasazi::MultiVecTraits< ST, MVMVT
 
typedef Anasazi::OperatorTraits< ST, MV, OPOPT
 

Private Attributes

Anasazi::OutputManager< ST > * Output_manager_pt
 Pointer to output manager. More...
 
LinearSolverLinear_solver_pt
 Pointer to a linear solver. More...
 
LinearSolverDefault_linear_solver_pt
 Pointer to a default linear solver. More...
 
int Spectrum
 
int NArnoldi
 Number of Arnoldi vectors to compute. More...
 
double Sigma
 Set the shifted value. More...
 
bool Small
 
bool Compute_eigenvectors
 Boolean to indicate whether or not to compute the eigenvectors. More...
 

Additional Inherited Members

- Protected Member Functions inherited from oomph::DistributableLinearAlgebraObject
void clear_distribution ()
 
- Protected Attributes inherited from oomph::EigenSolver
double Sigma_real
 

Detailed Description

Class for the Anasazi eigensolver.

Member Typedef Documentation

◆ MT

typedef SCT::magnitudeType oomph::ANASAZI::MT
private

◆ MV

typedef Anasazi::MultiVec<ST> oomph::ANASAZI::MV
private

◆ MVT

typedef Anasazi::MultiVecTraits<ST, MV> oomph::ANASAZI::MVT
private

◆ OP

typedef Anasazi::Operator<ST> oomph::ANASAZI::OP
private

◆ OPT

typedef Anasazi::OperatorTraits<ST, MV, OP> oomph::ANASAZI::OPT
private

◆ SCT

typedef Teuchos::ScalarTraits<ST> oomph::ANASAZI::SCT
private

◆ ST

typedef double oomph::ANASAZI::ST
private

Constructor & Destructor Documentation

◆ ANASAZI() [1/2]

oomph::ANASAZI::ANASAZI ( )
inline

Constructor.

612  : Linear_solver_pt(0),
614  Spectrum(0),
615  NArnoldi(10),
616  Sigma(0.0),
618  {
619  Output_manager_pt = new Anasazi::BasicOutputManager<ST>();
620  // Set verbosity level
621  int verbosity =
622  Anasazi::Warnings + Anasazi::FinalSummary + Anasazi::TimingDetails;
623  Output_manager_pt->setVerbosity(verbosity);
624 
625  // print greeting
626  Output_manager_pt->stream(Anasazi::Warnings)
627  << Anasazi::Anasazi_Version() << std::endl
628  << std::endl;
629  }
LinearSolver * Default_linear_solver_pt
Pointer to a default linear solver.
Definition: trilinos_eigen_solver.h:588
Anasazi::OutputManager< ST > * Output_manager_pt
Pointer to output manager.
Definition: trilinos_eigen_solver.h:582
int Spectrum
Definition: trilinos_eigen_solver.h:593
bool Compute_eigenvectors
Boolean to indicate whether or not to compute the eigenvectors.
Definition: trilinos_eigen_solver.h:606
int NArnoldi
Number of Arnoldi vectors to compute.
Definition: trilinos_eigen_solver.h:596
LinearSolver * Linear_solver_pt
Pointer to a linear solver.
Definition: trilinos_eigen_solver.h:585
double Sigma
Set the shifted value.
Definition: trilinos_eigen_solver.h:599

References Output_manager_pt.

◆ ANASAZI() [2/2]

oomph::ANASAZI::ANASAZI ( const ANASAZI )
inline

Empty copy constructor.

632 : Sigma(0.0) {}

◆ ~ANASAZI()

virtual oomph::ANASAZI::~ANASAZI ( )
inlinevirtual

Destructor, delete the linear solver.

635 {}

Member Function Documentation

◆ disable_compute_eigenvectors()

void oomph::ANASAZI::disable_compute_eigenvectors ( )
inline

Set to disable the computation of the eigenvectors.

657  {
658  Compute_eigenvectors = false;
659  }

References Compute_eigenvectors.

◆ enable_compute_eigenvectors()

void oomph::ANASAZI::enable_compute_eigenvectors ( )
inline

Set to enable the computation of the eigenvectors (default)

651  {
652  Compute_eigenvectors = true;
653  }

References Compute_eigenvectors.

◆ get_eigenvalues_left_of_shift()

void oomph::ANASAZI::get_eigenvalues_left_of_shift ( )
inline

Set the desired eigenvalues to be left of the shift.

770  {
771  Small = true;
772  }
bool Small
Definition: trilinos_eigen_solver.h:603

References Small.

◆ get_eigenvalues_right_of_shift()

void oomph::ANASAZI::get_eigenvalues_right_of_shift ( )
inline

Set the desired eigenvalues to be right of the shift.

776  {
777  Small = false;
778  }

References Small.

◆ linear_solver_pt() [1/2]

LinearSolver*& oomph::ANASAZI::linear_solver_pt ( )
inline

Return a pointer to the linear solver object.

800  {
801  return Linear_solver_pt;
802  }

References Linear_solver_pt.

Referenced by solve_eigenproblem().

◆ linear_solver_pt() [2/2]

LinearSolver* const& oomph::ANASAZI::linear_solver_pt ( ) const
inline

Return a pointer to the linear solver object (const version)

806  {
807  return Linear_solver_pt;
808  }

References Linear_solver_pt.

◆ narnoldi() [1/2]

int& oomph::ANASAZI::narnoldi ( )
inline

Access function for the number of Arnoldi vectors.

639  {
640  return NArnoldi;
641  }

References NArnoldi.

◆ narnoldi() [2/2]

const int& oomph::ANASAZI::narnoldi ( ) const
inline

Access function for the number of Arnoldi vectors (const version)

645  {
646  return NArnoldi;
647  }

References NArnoldi.

◆ solve_eigenproblem()

void oomph::ANASAZI::solve_eigenproblem ( Problem *const &  problem_pt,
const int n_eval,
Vector< std::complex< double >> &  eigenvalue,
Vector< DoubleVector > &  eigenvector 
)
inlinevirtual

Solve the eigen problem.

Implements oomph::EigenSolver.

666  {
667  // No access to sigma, so set from sigma real
668  Sigma = Sigma_real;
669 
670  // Initially be dumb here
671  Linear_solver_pt = problem_pt->linear_solver_pt();
672 
673  // Let's make the initial vector
674  Teuchos::RCP<DoubleMultiVector> initial = Teuchos::rcp(
675  new DoubleMultiVector(1, problem_pt->dof_distribution_pt()));
676  Anasazi::MultiVecTraits<double, DoubleMultiVector>::MvRandom(*initial);
677 
678  // Make the operator
679  Teuchos::RCP<DoubleMultiVectorOperator> Op_pt =
680  Teuchos::rcp(new ProblemBasedShiftInvertOperator(
681  problem_pt, this->linear_solver_pt(), Sigma));
682 
683  // Create the basic problem
684  Teuchos::RCP<Anasazi::BasicEigenproblem<double,
685  DoubleMultiVector,
686  DoubleMultiVectorOperator>>
687  anasazi_pt = Teuchos::rcp(
688  new Anasazi::BasicEigenproblem<double,
689  DoubleMultiVector,
690  DoubleMultiVectorOperator>(Op_pt,
691  initial));
692 
693  // Think I have it?
694  anasazi_pt->setHermitian(false);
695 
696  // set the number of eigenvalues requested
697  anasazi_pt->setNEV(n_eval);
698 
699  // Inform the eigenproblem that you are done passing it information
700  if (anasazi_pt->setProblem() != true)
701  {
702  oomph_info << "Anasazi::BasicEigenproblem::setProblem() had an error."
703  << std::endl
704  << "End Result: TEST FAILED" << std::endl;
705  }
706 
707  // Create the solver manager
708  // No need to have ncv specificed, Triliinos has a sensible default
709  // int ncv = 10;
710  MT tol = 1.0e-10;
711  int verbosity =
712  Anasazi::Warnings + Anasazi::FinalSummary + Anasazi::TimingDetails;
713 
714  Teuchos::ParameterList MyPL;
715  MyPL.set("Which", "LM");
716  MyPL.set("Block Size", 1);
717  // MyPL.set( "Num Blocks", ncv );
718  MyPL.set("Maximum Restarts", 500);
719  MyPL.set("Orthogonalization", "DGKS");
720  MyPL.set("Verbosity", verbosity);
721  MyPL.set("Convergence Tolerance", tol);
722  Anasazi::BlockKrylovSchurSolMgr<double,
723  DoubleMultiVector,
724  DoubleMultiVectorOperator>
725  BKS(anasazi_pt, MyPL);
726 
727  // Solve the problem to the specified tolerances or length
728  Anasazi::ReturnType ret = BKS.solve();
729  bool testFailed = false;
730  if (ret != Anasazi::Converged)
731  {
732  testFailed = true;
733  }
734 
735  if (testFailed)
736  {
737  oomph_info << "Eigensolver not converged\n";
738  }
739 
740  // Get the eigenvalues and eigenvectors from the eigenproblem
741  Anasazi::Eigensolution<double, DoubleMultiVector> sol =
742  anasazi_pt->getSolution();
743  std::vector<Anasazi::Value<double>> evals = sol.Evals;
744  Teuchos::RCP<DoubleMultiVector> evecs = sol.Evecs;
745 
746  eigenvalue.resize(evals.size());
747  eigenvector.resize(evals.size());
748 
749  for (unsigned i = 0; i < evals.size(); i++)
750  {
751  // Undo shift and invert
752  double a = evals[i].realpart;
753  double b = evals[i].imagpart;
754  double det = a * a + b * b;
755  eigenvalue[i] = std::complex<double>(a / det + Sigma, -b / det);
756 
757  // Now set the eigenvectors, I hope
758  eigenvector[i].build(evecs->distribution_pt());
759  unsigned nrow_local = evecs->nrow_local();
760  // Would be faster with pointers, but I'll sort that out later!
761  for (unsigned n = 0; n < nrow_local; n++)
762  {
763  eigenvector[i][n] = (*evecs)(i, n);
764  }
765  }
766  }
int i
Definition: BiCGSTAB_step_by_step.cpp:9
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
Scalar * b
Definition: benchVecAdd.cpp:17
SCT::magnitudeType MT
Definition: trilinos_eigen_solver.h:575
LinearSolver *& linear_solver_pt()
Return a pointer to the linear solver object.
Definition: trilinos_eigen_solver.h:799
unsigned nrow_local() const
access function for the num of local rows on this processor.
Definition: linear_algebra_distribution.h:469
double Sigma_real
Definition: eigen_solver.h:65
Eigen::DenseIndex ret
Definition: level1_cplx_impl.h:43
const Scalar * a
Definition: level2_cplx_impl.h:32
OomphInfo oomph_info
Definition: oomph_definitions.cc:319

References a, b, oomph::Problem::dof_distribution_pt(), i, oomph::Problem::linear_solver_pt(), Linear_solver_pt, linear_solver_pt(), n, oomph::DistributableLinearAlgebraObject::nrow_local(), oomph::oomph_info, ret, Sigma, and oomph::EigenSolver::Sigma_real.

◆ track_eigenvalue_imaginary_part()

void oomph::ANASAZI::track_eigenvalue_imaginary_part ( )
inline

Set the imaginary part fo the quantity of interest.

788  {
789  Spectrum = -1;
790  }

References Spectrum.

◆ track_eigenvalue_magnitude()

void oomph::ANASAZI::track_eigenvalue_magnitude ( )
inline

Set the magnitude to be the quantity of interest.

794  {
795  Spectrum = 0;
796  }

References Spectrum.

◆ track_eigenvalue_real_part()

void oomph::ANASAZI::track_eigenvalue_real_part ( )
inline

Set the real part to be the quantity of interest (default)

782  {
783  Spectrum = 1;
784  }

References Spectrum.

Member Data Documentation

◆ Compute_eigenvectors

bool oomph::ANASAZI::Compute_eigenvectors
private

Boolean to indicate whether or not to compute the eigenvectors.

Referenced by disable_compute_eigenvectors(), and enable_compute_eigenvectors().

◆ Default_linear_solver_pt

LinearSolver* oomph::ANASAZI::Default_linear_solver_pt
private

Pointer to a default linear solver.

◆ Linear_solver_pt

LinearSolver* oomph::ANASAZI::Linear_solver_pt
private

Pointer to a linear solver.

Referenced by linear_solver_pt(), and solve_eigenproblem().

◆ NArnoldi

int oomph::ANASAZI::NArnoldi
private

Number of Arnoldi vectors to compute.

Referenced by narnoldi().

◆ Output_manager_pt

Anasazi::OutputManager<ST>* oomph::ANASAZI::Output_manager_pt
private

Pointer to output manager.

Referenced by ANASAZI().

◆ Sigma

double oomph::ANASAZI::Sigma
private

Set the shifted value.

Referenced by solve_eigenproblem().

◆ Small

bool oomph::ANASAZI::Small
private

Boolean to set which part of the spectrum left (default) or right of the shifted value.

Referenced by get_eigenvalues_left_of_shift(), and get_eigenvalues_right_of_shift().

◆ Spectrum

int oomph::ANASAZI::Spectrum
private

Integer to set whether the real, imaginary or magnitude is required to be small or large.

Referenced by track_eigenvalue_imaginary_part(), track_eigenvalue_magnitude(), and track_eigenvalue_real_part().


The documentation for this class was generated from the following file: