iterative-solver 0.0
molpro::linalg::itsolv::IterativeSolver< R, Q, P > Class Template Referenceabstract

Base class defining the interface common to all iterative solvers. More...

#include <IterativeSolver.h>

Inheritance diagram for molpro::linalg::itsolv::IterativeSolver< R, Q, P >:

Detailed Description

template<class R, class Q, class P>
class molpro::linalg::itsolv::IterativeSolver< R, Q, P >

Base class defining the interface common to all iterative solvers.

As well as through the interface, some behaviour (profiling, and tuning of BLAS operations) is influenced by the contents of a global molpro::Options object, which defaults to molpro::Options("ITERATIVE-SOLVER", "").

Template Parameters
Rcontainer for "working-set" vectors. These are typically implemented in memory, and are created by the client program. R vectors are never created inside IterativeSolver.
Qcontainer for other vectors. These are typically implemented on backing store and/or distributed across processors. IterativeSolver constructs a number of instances of Q containers to store history.
Pa class that specifies the definition of a single P-space vector, which is a strictly sparse vector in the underlying space.

Public Types

using value_type = typename R::value_type
 The underlying type of elements of vectors. More...
 
using scalar_type = typename array::ArrayHandler< R, Q >::value_type
 
using value_type_abs = typename array::ArrayHandler< R, R >::value_type_abs
 
using VectorP = std::vector< value_type >
 
using fapply_on_p_type = std::function< void(const std::vector< VectorP > &, const CVecRef< P > &, const VecRef< R > &)>
 

Public Member Functions

virtual ~IterativeSolver ()=default
 
 IterativeSolver ()=default
 
 IterativeSolver (const IterativeSolver< R, Q, P > &)=delete
 
IterativeSolver< R, Q, P > & operator= (const IterativeSolver< R, Q, P > &)=delete
 
 IterativeSolver (IterativeSolver< R, Q, P > &&) noexcept=default
 
IterativeSolver< R, Q, P > & operator= (IterativeSolver< R, Q, P > &&) noexcept=default
 
virtual bool solve (const VecRef< R > &parameters, const VecRef< R > &actions, const Problem< R > &problem, bool generate_initial_guess=false)=0
 Simplified one-call solver. More...
 
virtual bool solve (R &parameters, R &actions, const Problem< R > &problem, bool generate_initial_guess=false)=0
 
virtual bool solve (std::vector< R > &parameters, std::vector< R > &actions, const Problem< R > &problem, bool generate_initial_guess=false)=0
 
virtual int add_vector (const VecRef< R > &parameters, const VecRef< R > &actions)=0
 Take, typically, a current solution and residual, and add it to the solution space. More...
 
virtual int add_vector (std::vector< R > &parameters, std::vector< R > &action)=0
 
virtual int add_vector (R &parameters, R &action, value_type value=0)=0
 
virtual size_t add_p (const CVecRef< P > &pparams, const array::Span< value_type > &pp_action_matrix, const VecRef< R > &parameters, const VecRef< R > &action, fapply_on_p_type apply_p)=0
 Add P-space vectors to the expansion set for linear methods. More...
 
virtual void clearP ()=0
 
virtual void solution (const std::vector< int > &roots, const VecRef< R > &parameters, const VecRef< R > &residual)=0
 Construct solution and residual for a given set of roots. More...
 
virtual void solution_params (const std::vector< int > &roots, const VecRef< R > &parameters)=0
 Constructs parameters of selected roots. More...
 
virtual size_t end_iteration (const VecRef< R > &parameters, const VecRef< R > &residual)=0
 Behaviour depends on the solver. More...
 
virtual bool end_iteration_needed ()=0
 signal whether end_iteration should be called More...
 
virtual std::vector< size_t > suggest_p (const CVecRef< R > &solution, const CVecRef< R > &residual, size_t max_number, double threshold)=0
 Get the solver's suggestion of which degrees of freedom would be best to add to the P-space. More...
 
virtual void solution (const std::vector< int > &roots, std::vector< R > &parameters, std::vector< R > &residual)=0
 
virtual void solution (R &parameters, R &residual)=0
 
virtual void solution_params (const std::vector< int > &roots, std::vector< R > &parameters)=0
 
virtual void solution_params (R &parameters)=0
 
virtual size_t end_iteration (std::vector< R > &parameters, std::vector< R > &action)=0
 
virtual size_t end_iteration (R &parameters, R &action)=0
 
virtual const std::vector< int > & working_set () const =0
 Working set of roots that are not yet converged. More...
 
virtual std::vector< scalar_typeworking_set_eigenvalues () const
 
virtual size_t n_roots () const =0
 
virtual void set_n_roots (size_t nroots)=0
 
virtual const std::vector< scalar_type > & errors () const =0
 
virtual const Statisticsstatistics () const =0
 
virtual void report (std::ostream &cout, bool endl=true) const =0
 Writes a report to cout output stream. More...
 
virtual void report () const =0
 Writes a report to std::cout. More...
 
virtual void set_convergence_threshold (double thresh)=0
 Sets the convergence threshold. More...
 
virtual double convergence_threshold () const =0
 Reports the convergence threshold. More...
 
virtual void set_convergence_threshold_value (double thresh)=0
 Sets the value convergence threshold. More...
 
virtual double convergence_threshold_value () const =0
 Reports the value convergence threshold. More...
 
virtual void set_verbosity (Verbosity v)=0
 
virtual void set_verbosity (int v)=0
 
virtual Verbosity get_verbosity () const =0
 
virtual void set_max_iter (int n)=0
 
virtual int get_max_iter () const =0
 
virtual void set_max_p (int n)=0
 
virtual int get_max_p () const =0
 
virtual void set_p_threshold (double thresh)=0
 
virtual double get_p_threshold () const =0
 
virtual const subspace::Dimensionsdimensions () const =0
 
virtual void set_options (const Options &options)=0
 
virtual std::shared_ptr< Optionsget_options () const =0
 
virtual scalar_type value () const =0
 Report the function value for the current optimum solution. More...
 
virtual bool nonlinear () const =0
 Report whether the class is a non-linear solver. More...
 
virtual void set_profiler (molpro::profiler::Profiler &profiler)=0
 Attach a profiler in order to collect performance data. More...
 
virtual const std::shared_ptr< molpro::profiler::Profiler > & profiler () const =0
 
virtual bool test_problem (const Problem< R > &problem, R &v0, R &v1, int verbosity=0, double threshold=1e-5) const =0
 Test a supplied problem class. More...
 

Member Typedef Documentation

◆ fapply_on_p_type

template<class R , class Q , class P >
using molpro::linalg::itsolv::IterativeSolver< R, Q, P >::fapply_on_p_type = std::function<void(const std::vector<VectorP>&, const CVecRef<P>&, const VecRef<R>&)>

Function type for applying matrix to the P space vectors and accumulating result in a residual

◆ scalar_type

template<class R , class Q , class P >
using molpro::linalg::itsolv::IterativeSolver< R, Q, P >::scalar_type = typename array::ArrayHandler<R, Q>::value_type

The type of scalar products of vectors

◆ value_type

template<class R , class Q , class P >
using molpro::linalg::itsolv::IterativeSolver< R, Q, P >::value_type = typename R::value_type

The underlying type of elements of vectors.

◆ value_type_abs

template<class R , class Q , class P >
using molpro::linalg::itsolv::IterativeSolver< R, Q, P >::value_type_abs = typename array::ArrayHandler<R, R>::value_type_abs

◆ VectorP

template<class R , class Q , class P >
using molpro::linalg::itsolv::IterativeSolver< R, Q, P >::VectorP = std::vector<value_type>

type for vectors projected on to P space, each element is a coefficient for the corresponding P space parameter

Constructor & Destructor Documentation

◆ ~IterativeSolver()

template<class R , class Q , class P >
virtual molpro::linalg::itsolv::IterativeSolver< R, Q, P >::~IterativeSolver ( )
virtualdefault

◆ IterativeSolver() [1/3]

template<class R , class Q , class P >
molpro::linalg::itsolv::IterativeSolver< R, Q, P >::IterativeSolver ( )
default

◆ IterativeSolver() [2/3]

template<class R , class Q , class P >
molpro::linalg::itsolv::IterativeSolver< R, Q, P >::IterativeSolver ( const IterativeSolver< R, Q, P > &  )
delete

◆ IterativeSolver() [3/3]

template<class R , class Q , class P >
molpro::linalg::itsolv::IterativeSolver< R, Q, P >::IterativeSolver ( IterativeSolver< R, Q, P > &&  )
defaultnoexcept

Member Function Documentation

◆ add_p()

template<class R , class Q , class P >
virtual size_t molpro::linalg::itsolv::IterativeSolver< R, Q, P >::add_p ( const CVecRef< P > &  pparams,
const array::Span< value_type > &  pp_action_matrix,
const VecRef< R > &  parameters,
const VecRef< R > &  action,
fapply_on_p_type  apply_p 
)
pure virtual

Add P-space vectors to the expansion set for linear methods.

Note
the apply_p function is stored and used by the solver internally.
Parameters
Pparamsthe vectors to add. Each Pvector specifies a sparse vector in the underlying space. The size of the P space must be at least the number of roots to be sought.
pp_action_matrixMatrix projected onto the existing+new, new P space. It should be provided as a 1-dimensional array, with the existing+new index running fastest.
parametersUsed as scratch working space
actionOn exit, the residual of the interpolated solution. The contribution from the new, and any existing, P parameters is missing, and should be added in subsequently.
apply_pA function that evaluates the action of the matrix on vectors in the P space
Returns
The number of vectors contained in parameters, action, parametersP

Implemented in molpro::linalg::itsolv::IterativeSolverTemplate< LinearEigensystem, R, R, std::map< size_t, typename R::value_type > >, molpro::linalg::itsolv::IterativeSolverTemplate< LinearEquations, R, R, std::map< size_t, typename R::value_type > >, and molpro::linalg::itsolv::IterativeSolverTemplate< Optimize, R, R, std::map< size_t, typename R::value_type > >.

◆ add_vector() [1/3]

template<class R , class Q , class P >
virtual int molpro::linalg::itsolv::IterativeSolver< R, Q, P >::add_vector ( const VecRef< R > &  parameters,
const VecRef< R > &  actions 
)
pure virtual

Take, typically, a current solution and residual, and add it to the solution space.

Parameters
parametersOn input, the current solution or expansion vector. On exit, undefined.
actionsOn input, the residual for parameters (non-linear), or action of matrix on parameters (linear). On exit, a vector set that should be preconditioned before returning to end_iteration().
valueThe value of the objective function for parameters. Used only in Optimize classes.
Returns
The size of the new working set. In non-linear optimisation, the special value -1 can also be returned, indicating that preconditioning should not be carried out on action.

Implemented in molpro::linalg::itsolv::IterativeSolverTemplate< LinearEigensystem, R, R, std::map< size_t, typename R::value_type > >, molpro::linalg::itsolv::IterativeSolverTemplate< LinearEquations, R, R, std::map< size_t, typename R::value_type > >, and molpro::linalg::itsolv::IterativeSolverTemplate< Optimize, R, R, std::map< size_t, typename R::value_type > >.

◆ add_vector() [2/3]

◆ add_vector() [3/3]

◆ clearP()

◆ convergence_threshold()

◆ convergence_threshold_value()

◆ dimensions()

◆ end_iteration() [1/3]

template<class R , class Q , class P >
virtual size_t molpro::linalg::itsolv::IterativeSolver< R, Q, P >::end_iteration ( const VecRef< R > &  parameters,
const VecRef< R > &  residual 
)
pure virtual

◆ end_iteration() [2/3]

◆ end_iteration() [3/3]

◆ end_iteration_needed()

◆ errors()

◆ get_max_iter()

◆ get_max_p()

◆ get_options()

◆ get_p_threshold()

◆ get_verbosity()

◆ n_roots()

◆ nonlinear()

◆ operator=() [1/2]

template<class R , class Q , class P >
IterativeSolver< R, Q, P > & molpro::linalg::itsolv::IterativeSolver< R, Q, P >::operator= ( const IterativeSolver< R, Q, P > &  )
delete

◆ operator=() [2/2]

template<class R , class Q , class P >
IterativeSolver< R, Q, P > & molpro::linalg::itsolv::IterativeSolver< R, Q, P >::operator= ( IterativeSolver< R, Q, P > &&  )
defaultnoexcept

◆ profiler()

◆ report() [1/2]

◆ report() [2/2]

◆ set_convergence_threshold()

◆ set_convergence_threshold_value()

◆ set_max_iter()

◆ set_max_p()

◆ set_n_roots()

◆ set_options()

◆ set_p_threshold()

◆ set_profiler()

◆ set_verbosity() [1/2]

◆ set_verbosity() [2/2]

◆ solution() [1/3]

template<class R , class Q , class P >
virtual void molpro::linalg::itsolv::IterativeSolver< R, Q, P >::solution ( const std::vector< int > &  roots,
const VecRef< R > &  parameters,
const VecRef< R > &  residual 
)
pure virtual

◆ solution() [2/3]

◆ solution() [3/3]

◆ solution_params() [1/3]

◆ solution_params() [2/3]

◆ solution_params() [3/3]

◆ solve() [1/3]

template<class R , class Q , class P >
virtual bool molpro::linalg::itsolv::IterativeSolver< R, Q, P >::solve ( const VecRef< R > &  parameters,
const VecRef< R > &  actions,
const Problem< R > &  problem,
bool  generate_initial_guess = false 
)
pure virtual

Simplified one-call solver.

Parameters
parametersA set of scratch vectors. On entry, these vectors should be filled with starting guesses. Where possible, the number of vectors should be equal to the number of solutions sought, but a smaller array is permitted.
actionsA set of scratch vectors. It should have the same size as parameters.
problemA Problem object defining the problem to be solved
generate_initial_guessWhether to start with a guess based on diagonal elements (true) or on the contents of parameters (false)
Returns
true if the solution was found

Implemented in molpro::linalg::itsolv::IterativeSolverTemplate< LinearEigensystem, R, R, std::map< size_t, typename R::value_type > >, molpro::linalg::itsolv::IterativeSolverTemplate< LinearEquations, R, R, std::map< size_t, typename R::value_type > >, molpro::linalg::itsolv::IterativeSolverTemplate< Optimize, R, R, std::map< size_t, typename R::value_type > >, and molpro::linalg::itsolv::LinearEquationsDavidson< R, Q, P >.

◆ solve() [2/3]

◆ solve() [3/3]

◆ statistics()

◆ suggest_p()

template<class R , class Q , class P >
virtual std::vector< size_t > molpro::linalg::itsolv::IterativeSolver< R, Q, P >::suggest_p ( const CVecRef< R > &  solution,
const CVecRef< R > &  residual,
size_t  max_number,
double  threshold 
)
pure virtual

Get the solver's suggestion of which degrees of freedom would be best to add to the P-space.

Parameters
solutionCurrent solution
residualCurrent residual
max_numberSuggest no more than this number
thresholdSuggest only axes for which the current residual and update indicate an energy improvement in the next iteration of this amount or more.
Returns

Implemented in molpro::linalg::itsolv::IterativeSolverTemplate< LinearEigensystem, R, R, std::map< size_t, typename R::value_type > >, molpro::linalg::itsolv::IterativeSolverTemplate< LinearEquations, R, R, std::map< size_t, typename R::value_type > >, and molpro::linalg::itsolv::IterativeSolverTemplate< Optimize, R, R, std::map< size_t, typename R::value_type > >.

◆ test_problem()

template<class R , class Q , class P >
virtual bool molpro::linalg::itsolv::IterativeSolver< R, Q, P >::test_problem ( const Problem< R > &  problem,
R &  v0,
R &  v1,
int  verbosity = 0,
double  threshold = 1e-5 
) const
pure virtual

◆ value()

◆ working_set()

◆ working_set_eigenvalues()

template<class R , class Q , class P >
virtual std::vector< scalar_type > molpro::linalg::itsolv::IterativeSolver< R, Q, P >::working_set_eigenvalues ( ) const
inlinevirtual

The calculated eigenvalues for roots in the working set (eigenvalue problems) or zero (otherwise)

Reimplemented in molpro::linalg::itsolv::LinearEigensystemDavidson< R, Q, P >, and molpro::linalg::itsolv::LinearEigensystemRSPT< R, Q, P >.