iterative-solver 0.0
molpro::linalg::itsolv Namespace Reference

4-parameter interpolation of a 1-dimensional function given two points for which function values and first derivatives are known. More...

Namespaces

namespace  detail
 
namespace  subspace
 
namespace  util
 

Classes

class  ArrayHandlers
 Class, containing a collection of array handlers used in IterativeSolver Provides a Builder sub-class, containing member functions to accumulate user-defined handlers, as well as determine and instantiate default handlers, based on container (template parameter) types. More...
 
struct  CastOptions
 Safely down-cast Options to one of the implementations. More...
 
struct  decay
 Decays CV and reference qualifiers same as std::decay, but also decays std::reference_wrapper to its underlying type. More...
 
struct  decay< std::reference_wrapper< T > >
 
struct  has_iterator
 
class  Interpolate
 
struct  is_complex
 
struct  is_complex< std::complex< T > >
 
class  IterativeSolver
 Base class defining the interface common to all iterative solvers. More...
 
class  IterativeSolverTemplate
 Implements IterativeSolver interface that is common to all solvers. More...
 
class  LinearEigensystem
 Interface for a specific iterative solver, it can add special member functions or variables. More...
 
class  LinearEigensystemDavidson
 One specific implementation of LinearEigensystem using Davidson's algorithm with modifications to manage near linear dependencies, and consequent numerical noise, in candidate expansion vectors. More...
 
struct  LinearEigensystemDavidsonOptions
 Allows setting and getting of options for LinearEigensystemDavidson instance via IterativeSolver base class. More...
 
struct  LinearEigensystemOptions
 
class  LinearEigensystemRSPT
 One specific implementation of LinearEigensystem using Davidson's algorithm with modifications to manage near linear dependencies, and consequent numerical noise, in candidate expansion vectors. More...
 
struct  LinearEigensystemRSPTOptions
 Allows setting and getting of options for LinearEigensystemRSPT instance via IterativeSolver base class. More...
 
class  LinearEquations
 
class  LinearEquationsDavidson
 Solves a system of linear equation, A x = b. More...
 
struct  LinearEquationsDavidsonOptions
 Allows setting and getting of options for LinearEquationsDavidson instance via IterativeSolver base class. More...
 
struct  LinearEquationsOptions
 
struct  Logger
 A dummy structured logger. More...
 
class  NonLinearEquations
 Solves non-linear system of equations using methods such as DIIS. More...
 
struct  NonLinearEquationsDIISOptions
 Allows setting and getting of options for NonLinearEquationsDIIS instance via IterativeSolver base class. More...
 
struct  NonLinearEquationsOptions
 
class  Optimize
 Optimises to a stationary point using methods such as L-BFGS. More...
 
class  OptimizeBFGS
 A class that optimises a function using a Quasi-Newton or other method. More...
 
struct  OptimizeBFGSOptions
 Allows setting and getting of options for OptimizeBFGS instance via IterativeSolver base class. More...
 
struct  OptimizeOptions
 
class  OptimizeSD
 A class that optimises a function using a simple steepest-descent method. More...
 
struct  OptimizeSDOptions
 Allows setting and getting of options for OptimizeBFGS instance via IterativeSolver base class. More...
 
struct  Options
 Access point for different options in iterative solvers. More...
 
class  Problem
 Abstract class defining the problem-specific interface for the simplified solver interface to IterativeSolver. More...
 
class  SolverFactory
 Factory for creating instances of specific solver implementations from the corresponding Options object. More...
 
struct  Statistics
 Information about performance of IterativeSolver instance. More...
 
struct  SVD
 Stores a singular value and corresponding left and right singular vectors. More...
 

Typedefs

template<class A >
using VecRef = std::vector< std::reference_wrapper< A > >
 
template<class A >
using CVecRef = std::vector< std::reference_wrapper< const A > >
 
template<class T >
using decay_t = typename decay< T >::type
 
using options_map = std::map< std::string, std::string >
 

Enumerations

enum class  Verbosity {
  None = 0 ,
  Summary = 1 ,
  Iteration = 2 ,
  Detailed = 3
}
 Specifies how much detail IterativeSolver::solve should write to output. More...
 

Functions

int eigensolver_lapacke_dsyev (const std::vector< double > &matrix, std::vector< double > &eigenvectors, std::vector< double > &eigenvalues, const size_t dimension)
 
std::list< SVD< double > > eigensolver_lapacke_dsyev (size_t dimension, std::vector< double > &matrix)
 
std::list< SVD< double > > eigensolver_lapacke_dsyev (size_t dimension, const molpro::linalg::array::span::Span< double > &matrix)
 
template<typename value_type >
size_t get_rank (std::vector< value_type > eigenvalues, value_type threshold)
 
template<typename value_type >
size_t get_rank (std::list< SVD< value_type > > svd_system, value_type threshold)
 
template<typename value_type , typename std::enable_if_t<!is_complex< value_type >{}, std::nullptr_t > = nullptr>
std::list< SVD< value_type > > svd_system (size_t nrows, size_t ncols, const array::Span< value_type > &m, double threshold, bool hermitian=false, bool reduce_to_rank=false)
 Performs singular value decomposition and returns SVD objects for singular values less than threshold, sorted in ascending order. More...
 
template<typename value_type >
void printMatrix (const std::vector< value_type > &, size_t rows, size_t cols, std::string title="", std::ostream &s=molpro::cout)
 
template<typename value_type , typename std::enable_if_t< is_complex< value_type >{}, int > = 0>
void eigenproblem (std::vector< value_type > &eigenvectors, std::vector< value_type > &eigenvalues, const std::vector< value_type > &matrix, const std::vector< value_type > &metric, size_t dimension, bool hermitian, double svdThreshold, int verbosity, bool condone_complex)
 
template<typename value_type , typename std::enable_if_t< is_complex< value_type >{}, int > = 0>
void solve_LinearEquations (std::vector< value_type > &solution, std::vector< value_type > &eigenvalues, const std::vector< value_type > &matrix, const std::vector< value_type > &metric, const std::vector< value_type > &rhs, size_t dimension, size_t nroot, double augmented_hessian, double svdThreshold, int verbosity)
 
template<typename value_type , typename std::enable_if_t< is_complex< value_type >{}, int > = 0>
void solve_DIIS (std::vector< value_type > &solution, const std::vector< value_type > &matrix, size_t dimension, double svdThreshold, int verbosity=0)
 
template void printMatrix< double > (const std::vector< double > &, size_t rows, size_t cols, std::string title, std::ostream &s)
 
template std::list< SVD< double > > svd_system (size_t nrows, size_t ncols, const array::Span< double > &m, double threshold, bool hermitian, bool reduce_to_rank)
 
template void eigenproblem< double > (std::vector< double > &eigenvectors, std::vector< double > &eigenvalues, const std::vector< double > &matrix, const std::vector< double > &metric, const size_t dimension, bool hermitian, double svdThreshold, int verbosity, bool condone_complex)
 
template void solve_LinearEquations< double > (std::vector< double > &solution, std::vector< double > &eigenvalues, const std::vector< double > &matrix, const std::vector< double > &metric, const std::vector< double > &rhs, size_t dimension, size_t nroot, double augmented_hessian, double svdThreshold, int verbosity)
 
template void solve_DIIS< double > (std::vector< double > &solution, const std::vector< double > &matrix, const size_t dimension, double svdThreshold, int verbosity)
 
template void printMatrix< std::complex< double > > (const std::vector< std::complex< double > > &, size_t rows, size_t cols, std::string title, std::ostream &s)
 
template std::list< SVD< std::complex< double > > > svd_system (size_t nrows, size_t ncols, const array::Span< std::complex< double > > &m, double threshold, bool hermitian, bool reduce_to_rank)
 
template void eigenproblem< std::complex< double > > (std::vector< std::complex< double > > &eigenvectors, std::vector< std::complex< double > > &eigenvalues, const std::vector< std::complex< double > > &matrix, const std::vector< std::complex< double > > &metric, const size_t dimension, bool hermitian, double svdThreshold, int verbosity, bool condone_complex)
 
template void solve_LinearEquations< std::complex< double > > (std::vector< std::complex< double > > &solution, std::vector< std::complex< double > > &eigenvalues, const std::vector< std::complex< double > > &matrix, const std::vector< std::complex< double > > &metric, const std::vector< std::complex< double > > &rhs, size_t dimension, size_t nroot, double augmented_hessian, double svdThreshold, int verbosity)
 
template void solve_DIIS< std::complex< double > > (std::vector< std::complex< double > > &solution, const std::vector< std::complex< double > > &matrix, const size_t dimension, double svdThreshold, int verbosity)
 
template<typename R , typename Q , typename P >
void read_handler_counts (std::shared_ptr< Statistics > stats, std::shared_ptr< ArrayHandlers< R, Q, P > > handlers)
 
std::ostream & operator<< (std::ostream &o, const Statistics &statistics)
 
template<typename T , typename = std::enable_if_t<std::is_base_of<molpro::linalg::array::DistrArray, T>::value>>
void precondition_default (const VecRef< T > &action, const std::vector< double > &shift, const T &diagonals)
 
template<class T >
void precondition_default (const VecRef< T > &action, const std::vector< double > &shift, const T &diagonals, typename T::iterator *=nullptr)
 
template<typename T , typename = std::enable_if_t<!std::is_base_of<molpro::linalg::array::DistrArray, T>::value>, class = void>
void precondition_default (const VecRef< T > &action, const std::vector< double > &shift, const T &diagonals, typename std::enable_if<!has_iterator< T >::value, void * >::type=nullptr)
 
template<class ForwardIt >
auto wrap (ForwardIt begin, ForwardIt end)
 Takes a begin and end iterators and returns a vector of references to each element. More...
 
template<class ForwardIt >
auto const_cast_wrap (ForwardIt begin, ForwardIt end)
 Takes a begin and end iterators and returns a vector of const-casted references to each element. More...
 
template<class ForwardIt >
auto cwrap (ForwardIt begin, ForwardIt end)
 Takes a begin and end iterators and returns a vector of references to each element. More...
 
template<class IterableContainer >
auto wrap (const IterableContainer &parameters)
 Wraps references for each parameter in an iterable container that implements begin()/end() More...
 
template<class IterableContainer >
auto wrap (IterableContainer &parameters)
 Wraps references for each parameter in an iterable container that implements begin()/end() More...
 
template<class IterableContainer >
auto const_cast_wrap (IterableContainer &parameters)
 Wraps const-casted references for each parameter in an iterable container that implements begin()/end() More...
 
template<class IterableContainer >
auto cwrap (IterableContainer &parameters)
 Wraps references for each parameter in an iterable container that implements begin()/end() More...
 
template<typename T , typename... S>
auto wrap_arg (T &&arg, S &&... args) -> std::enable_if_t< std::conjunction_v< std::is_same< decay_t< T >, decay_t< S > >... >, VecRef< decay_t< T > > >
 Constructs a vector of reference wrappers with provided arguments. More...
 
template<typename T , typename... S>
auto cwrap_arg (T &&arg, S &&... args) -> std::enable_if_t< std::conjunction_v< std::is_same< decay_t< T >, decay_t< S > >... >, CVecRef< decay_t< T > > >
 Constructs a vector of const reference wrappers with provided arguments. More...
 
template<typename R , typename ForwardIt >
std::vector< size_t > find_ref (const VecRef< R > &wparams, ForwardIt begin, ForwardIt end)
 Given wrapped references in wparams and a range of original parameters [begin, end), returns indices of parameters that are wrapped in this range. More...
 
template<typename R >
std::vector< size_t > find_ref (const CVecRef< R > &wparams, const CVecRef< R > &wparams_ref)
 Given two wrapped references returns indices of parameters that are in both sets. More...
 
template<typename T , typename U >
auto remove_elements (std::vector< T > params, const std::vector< U > &indices)
 Removes indices from a vector. More...
 
template<class R , class Q = R, class P = std::map<size_t, typename R::value_type>>
std::unique_ptr< LinearEigensystem< R, Q, P > > create_LinearEigensystem (const LinearEigensystemOptions &options, const std::shared_ptr< ArrayHandlers< R, Q, P > > &handlers=std::make_shared< molpro::linalg::itsolv::ArrayHandlers< R, Q, P > >())
 
template<class R , class Q = R, class P = std::map<size_t, typename R::value_type>>
std::unique_ptr< LinearEigensystem< R, Q, P > > create_LinearEigensystem (const std::string &method="Davidson", const std::string &options="", const std::shared_ptr< ArrayHandlers< R, Q, P > > &handlers=std::make_shared< molpro::linalg::itsolv::ArrayHandlers< R, Q, P > >())
 
template<class R , class Q = R, class P = std::map<size_t, typename R::value_type>>
std::unique_ptr< LinearEquations< R, Q, P > > create_LinearEquations (const LinearEquationsOptions &options, const std::shared_ptr< ArrayHandlers< R, Q, P > > &handlers=std::make_shared< molpro::linalg::itsolv::ArrayHandlers< R, Q, P > >())
 
template<class R , class Q = R, class P = std::map<size_t, typename R::value_type>>
std::unique_ptr< LinearEquations< R, Q, P > > create_LinearEquations (const std::string &method="Davidson", const std::string &options="", const std::shared_ptr< ArrayHandlers< R, Q, P > > &handlers=std::make_shared< molpro::linalg::itsolv::ArrayHandlers< R, Q, P > >())
 
template<class R , class Q = R, class P = std::map<size_t, typename R::value_type>>
std::tuple< bool, std::unique_ptr< LinearEquations< R, Q, P > > > Solve_LinearEquations (const VecRef< R > &parameters, const VecRef< R > &actions, const Problem< R > &problem, int verbosity=0, bool generate_initial_guess=true, const std::string &method="Davidson", const std::string &options="", const std::shared_ptr< ArrayHandlers< R, Q, P > > &handlers=std::make_shared< molpro::linalg::itsolv::ArrayHandlers< R, Q, P > >())
 
template<class R , class Q = R, class P = std::map<size_t, typename R::value_type>>
std::tuple< bool, std::unique_ptr< LinearEquations< R, Q, P > > > Solve_LinearEquations (std::vector< R > &parameters, std::vector< R > &actions, const Problem< R > &problem, int verbosity=0, bool generate_initial_guess=true, const std::string &method="Davidson", const std::string &options="", const std::shared_ptr< ArrayHandlers< R, Q, P > > &handlers=std::make_shared< molpro::linalg::itsolv::ArrayHandlers< R, Q, P > >())
 
template<class R , class Q = R, class P = std::map<size_t, typename R::value_type>>
auto Solve_LinearEquations (R &parameters, R &actions, const Problem< R > &problem, int verbosity=0, bool generate_initial_guess=true, const std::string &method="Davidson", const std::string &options="", const std::shared_ptr< ArrayHandlers< R, Q, P > > &handlers=std::make_shared< molpro::linalg::itsolv::ArrayHandlers< R, Q, P > >())
 
template<class R , class Q = R, class P = std::map<size_t, typename R::value_type>>
std::unique_ptr< NonLinearEquations< R, Q, P > > create_NonLinearEquations (const NonLinearEquationsOptions &options=NonLinearEquationsOptions{}, const std::shared_ptr< ArrayHandlers< R, Q, P > > &handlers=std::make_shared< molpro::linalg::itsolv::ArrayHandlers< R, Q, P > >())
 
template<class R , class Q = R, class P = std::map<size_t, typename R::value_type>>
std::unique_ptr< NonLinearEquations< R, Q, P > > create_NonLinearEquations (const std::string &method="DIIS", const std::string &options="", const std::shared_ptr< ArrayHandlers< R, Q, P > > &handlers=std::make_shared< molpro::linalg::itsolv::ArrayHandlers< R, Q, P > >())
 
template<class R , class Q = R, class P = std::map<size_t, typename R::value_type>>
std::unique_ptr< Optimize< R, Q, P > > create_Optimize (const OptimizeOptions &options=OptimizeOptions{}, const std::shared_ptr< ArrayHandlers< R, Q, P > > &handlers=std::make_shared< molpro::linalg::itsolv::ArrayHandlers< R, Q, P > >())
 
template<class R , class Q = R, class P = std::map<size_t, typename R::value_type>>
std::unique_ptr< Optimize< R, Q, P > > create_Optimize (const std::string &method="BFGS", const std::string &options="", const std::shared_ptr< ArrayHandlers< R, Q, P > > &handlers=std::make_shared< molpro::linalg::itsolv::ArrayHandlers< R, Q, P > >())
 
bool operator== (const Interpolate::point &lhs, const Interpolate::point &rhs)
 
std::ostream & operator<< (std::ostream &os, const Interpolate &interpolant)
 
std::ostream & operator<< (std::ostream &os, const Interpolate::point &p)
 
template<typename value_type >
std::list< SVD< value_type > > svd_eigen_jacobi (size_t nrows, size_t ncols, const array::Span< value_type > &m, double threshold)
 
template<typename value_type >
std::list< SVD< value_type > > svd_eigen_bdcsvd (size_t nrows, size_t ncols, const array::Span< value_type > &m, double threshold)
 
template void printMatrix< value_type > (const std::vector< value_type > &, size_t rows, size_t cols, std::string title, std::ostream &s)
 
template size_t get_rank< value_type > (std::vector< value_type > eigenvalues, value_type threshold)
 
template std::list< SVD< value_type > > svd_system< value_type > (size_t nrows, size_t ncols, const array::Span< value_type > &m, double threshold, bool hermitian, bool reduce_to_rank)
 
template void eigenproblem< value_type > (std::vector< value_type > &eigenvectors, std::vector< value_type > &eigenvalues, const std::vector< value_type > &matrix, const std::vector< value_type > &metric, size_t dimension, bool hermitian, double svdThreshold, int verbosity, bool condone_complex)
 
template void solve_LinearEquations< value_type > (std::vector< value_type > &solution, std::vector< value_type > &eigenvalues, const std::vector< value_type > &matrix, const std::vector< value_type > &metric, const std::vector< value_type > &rhs, size_t dimension, size_t nroot, double augmented_hessian, double svdThreshold, int verbosity)
 
template void solve_DIIS< value_type > (std::vector< value_type > &solution, const std::vector< value_type > &matrix, size_t dimension, double svdThreshold, int verbosity)
 

Detailed Description

4-parameter interpolation of a 1-dimensional function given two points for which function values and first derivatives are known.

Typedef Documentation

◆ CVecRef

template<class A >
using molpro::linalg::itsolv::CVecRef = typedef std::vector<std::reference_wrapper<const A> >

◆ decay_t

template<class T >
using molpro::linalg::itsolv::decay_t = typedef typename decay<T>::type

◆ options_map

using molpro::linalg::itsolv::options_map = typedef std::map<std::string, std::string>

◆ VecRef

template<class A >
using molpro::linalg::itsolv::VecRef = typedef std::vector<std::reference_wrapper<A> >

Enumeration Type Documentation

◆ Verbosity

Specifies how much detail IterativeSolver::solve should write to output.

Enumerator
None 

no output

Summary 

only summary at the start and end of solve

Iteration 

Summary plus minimal results at each iteration.

Detailed 

Iteration plus more detail such as operation count.

Function Documentation

◆ const_cast_wrap() [1/2]

template<class ForwardIt >
auto molpro::linalg::itsolv::const_cast_wrap ( ForwardIt  begin,
ForwardIt  end 
)

Takes a begin and end iterators and returns a vector of const-casted references to each element.

◆ const_cast_wrap() [2/2]

template<class IterableContainer >
auto molpro::linalg::itsolv::const_cast_wrap ( IterableContainer &  parameters)

Wraps const-casted references for each parameter in an iterable container that implements begin()/end()

Template Parameters
IterableContainershould have begin()/end() either as free or member functions
Relement type
Parameters
parametersparameters to wrap
Returns
vector of references to each element in parameters

◆ create_LinearEigensystem() [1/2]

template<class R , class Q = R, class P = std::map<size_t, typename R::value_type>>
std::unique_ptr< LinearEigensystem< R, Q, P > > molpro::linalg::itsolv::create_LinearEigensystem ( const LinearEigensystemOptions options,
const std::shared_ptr< ArrayHandlers< R, Q, P > > &  handlers = std::make_shared<molpro::linalg::itsolv::ArrayHandlers<R, Q, P>>() 
)

◆ create_LinearEigensystem() [2/2]

template<class R , class Q = R, class P = std::map<size_t, typename R::value_type>>
std::unique_ptr< LinearEigensystem< R, Q, P > > molpro::linalg::itsolv::create_LinearEigensystem ( const std::string &  method = "Davidson",
const std::string &  options = "",
const std::shared_ptr< ArrayHandlers< R, Q, P > > &  handlers = std::make_shared<molpro::linalg::itsolv::ArrayHandlers<R, Q, P>>() 
)

◆ create_LinearEquations() [1/2]

template<class R , class Q = R, class P = std::map<size_t, typename R::value_type>>
std::unique_ptr< LinearEquations< R, Q, P > > molpro::linalg::itsolv::create_LinearEquations ( const LinearEquationsOptions options,
const std::shared_ptr< ArrayHandlers< R, Q, P > > &  handlers = std::make_shared<molpro::linalg::itsolv::ArrayHandlers<R, Q, P>>() 
)

◆ create_LinearEquations() [2/2]

template<class R , class Q = R, class P = std::map<size_t, typename R::value_type>>
std::unique_ptr< LinearEquations< R, Q, P > > molpro::linalg::itsolv::create_LinearEquations ( const std::string &  method = "Davidson",
const std::string &  options = "",
const std::shared_ptr< ArrayHandlers< R, Q, P > > &  handlers = std::make_shared<molpro::linalg::itsolv::ArrayHandlers<R, Q, P>>() 
)

◆ create_NonLinearEquations() [1/2]

template<class R , class Q = R, class P = std::map<size_t, typename R::value_type>>
std::unique_ptr< NonLinearEquations< R, Q, P > > molpro::linalg::itsolv::create_NonLinearEquations ( const NonLinearEquationsOptions options = NonLinearEquationsOptions{},
const std::shared_ptr< ArrayHandlers< R, Q, P > > &  handlers = std::make_shared<molpro::linalg::itsolv::ArrayHandlers<R, Q, P>>() 
)

◆ create_NonLinearEquations() [2/2]

template<class R , class Q = R, class P = std::map<size_t, typename R::value_type>>
std::unique_ptr< NonLinearEquations< R, Q, P > > molpro::linalg::itsolv::create_NonLinearEquations ( const std::string &  method = "DIIS",
const std::string &  options = "",
const std::shared_ptr< ArrayHandlers< R, Q, P > > &  handlers = std::make_shared<molpro::linalg::itsolv::ArrayHandlers<R, Q, P>>() 
)

◆ create_Optimize() [1/2]

template<class R , class Q = R, class P = std::map<size_t, typename R::value_type>>
std::unique_ptr< Optimize< R, Q, P > > molpro::linalg::itsolv::create_Optimize ( const OptimizeOptions options = OptimizeOptions{},
const std::shared_ptr< ArrayHandlers< R, Q, P > > &  handlers = std::make_shared<molpro::linalg::itsolv::ArrayHandlers<R, Q, P>>() 
)

◆ create_Optimize() [2/2]

template<class R , class Q = R, class P = std::map<size_t, typename R::value_type>>
std::unique_ptr< Optimize< R, Q, P > > molpro::linalg::itsolv::create_Optimize ( const std::string &  method = "BFGS",
const std::string &  options = "",
const std::shared_ptr< ArrayHandlers< R, Q, P > > &  handlers = std::make_shared<molpro::linalg::itsolv::ArrayHandlers<R, Q, P>>() 
)

◆ cwrap() [1/2]

template<class ForwardIt >
auto molpro::linalg::itsolv::cwrap ( ForwardIt  begin,
ForwardIt  end 
)

Takes a begin and end iterators and returns a vector of references to each element.

◆ cwrap() [2/2]

template<class IterableContainer >
auto molpro::linalg::itsolv::cwrap ( IterableContainer &  parameters)

Wraps references for each parameter in an iterable container that implements begin()/end()

Template Parameters
IterableContainershould have begin()/end() either as free or member functions
Relement type
Parameters
parametersparameters to wrap
Returns
vector of constant references to each element in parameters

◆ cwrap_arg()

template<typename T , typename... S>
auto molpro::linalg::itsolv::cwrap_arg ( T &&  arg,
S &&...  args 
) -> std::enable_if_t<std::conjunction_v<std::is_same<decay_t<T>, decay_t<S>>...>, CVecRef<decay_t<T>>>

Constructs a vector of const reference wrappers with provided arguments.

◆ eigenproblem()

template<typename value_type , typename std::enable_if_t< is_complex< value_type >{}, int > = 0>
void molpro::linalg::itsolv::eigenproblem ( std::vector< value_type > &  eigenvectors,
std::vector< value_type > &  eigenvalues,
const std::vector< value_type > &  matrix,
const std::vector< value_type > &  metric,
size_t  dimension,
bool  hermitian,
double  svdThreshold,
int  verbosity,
bool  condone_complex 
)

◆ eigenproblem< double >()

template void molpro::linalg::itsolv::eigenproblem< double > ( std::vector< double > &  eigenvectors,
std::vector< double > &  eigenvalues,
const std::vector< double > &  matrix,
const std::vector< double > &  metric,
const size_t  dimension,
bool  hermitian,
double  svdThreshold,
int  verbosity,
bool  condone_complex 
)

◆ eigenproblem< std::complex< double > >()

template void molpro::linalg::itsolv::eigenproblem< std::complex< double > > ( std::vector< std::complex< double > > &  eigenvectors,
std::vector< std::complex< double > > &  eigenvalues,
const std::vector< std::complex< double > > &  matrix,
const std::vector< std::complex< double > > &  metric,
const size_t  dimension,
bool  hermitian,
double  svdThreshold,
int  verbosity,
bool  condone_complex 
)

◆ eigenproblem< value_type >()

template void molpro::linalg::itsolv::eigenproblem< value_type > ( std::vector< value_type > &  eigenvectors,
std::vector< value_type > &  eigenvalues,
const std::vector< value_type > &  matrix,
const std::vector< value_type > &  metric,
size_t  dimension,
bool  hermitian,
double  svdThreshold,
int  verbosity,
bool  condone_complex 
)

◆ eigensolver_lapacke_dsyev() [1/3]

int molpro::linalg::itsolv::eigensolver_lapacke_dsyev ( const std::vector< double > &  matrix,
std::vector< double > &  eigenvectors,
std::vector< double > &  eigenvalues,
const size_t  dimension 
)
inline

A wrapper function for lapacke_dsyev (linear eigensystem solver) from the lapack C interface (lapacke.h).

Parameters
[in]matrixthe input matrix (will not be altered!). Must be square. Must be dimension*dimension elements long.
[out]eigenvectorsthe matrix of eigenvectors, row major ordering. Must be square and the same size as matrix.
[out]eigenvaluesa list of eigenvalues. Must be the length specified by the 'dimension' parameter.
[in]dimensionlength of one axis of the matrix.
Returns
status. If 0, successful exit. If -i, the ith argument had an illegal value. If i, the algorithm failed to converge.

◆ eigensolver_lapacke_dsyev() [2/3]

std::list< SVD< double > > molpro::linalg::itsolv::eigensolver_lapacke_dsyev ( size_t  dimension,
const molpro::linalg::array::span::Span< double > &  matrix 
)
inline

A wrapper function for lapacke_dsyev (linear eigensystem solver) from the lapack C interface (lapacke.h).

Parameters
[in]dimensionlength of one axis of the matrix.
[in]matrixa span wrapping some data structure containing the elements of the matrix. Should be dimension*dimension in length.
Returns
a std::list of instances of SVD, a struct containing one eigenvalue and one eigenvector. For a real, symmetric matrix, these are equivalent to singular values, and S/D (which are both the same).

◆ eigensolver_lapacke_dsyev() [3/3]

std::list< SVD< double > > molpro::linalg::itsolv::eigensolver_lapacke_dsyev ( size_t  dimension,
std::vector< double > &  matrix 
)
inline

A wrapper function for lapacke_dsyev (linear eigensystem solver) from the lapack C interface (lapacke.h).

Parameters
[in]matrixthe input matrix (will not be altered!). Must be square. Must be dimension*dimension elements long.
[in]dimensionlength of one axis of the matrix.
Returns
a std::list of instances of SVD, a struct containing one eigenvalue and one eigenvector. For a real, symmetric matrix, these are equivalent to singular values, and S/D (which are both the same).

◆ find_ref() [1/2]

template<typename R >
std::vector< size_t > molpro::linalg::itsolv::find_ref ( const CVecRef< R > &  wparams,
const CVecRef< R > &  wparams_ref 
)

Given two wrapped references returns indices of parameters that are in both sets.

Parameters
wparamsfirst set of references
wparams_refsecond set of references

◆ find_ref() [2/2]

template<typename R , typename ForwardIt >
std::vector< size_t > molpro::linalg::itsolv::find_ref ( const VecRef< R > &  wparams,
ForwardIt  begin,
ForwardIt  end 
)

Given wrapped references in wparams and a range of original parameters [begin, end), returns indices of parameters that are wrapped in this range.

Parameters
wparamswrapped references to subset of params
beginstart of range of original paramaters
endend of range of original paramaters

◆ get_rank() [1/2]

template<typename value_type >
size_t molpro::linalg::itsolv::get_rank ( std::list< SVD< value_type > >  svd_system,
value_type  threshold 
)

Get the rank of some matrix, given a threshold.

Parameters
[in]svd_systema std::list containing SVD objects, such as those created by itsolv::svd_system
[in]thresholdthe threshold. Note that this is the normalised threshold, a value between 0 and 1, relative to the largest element in the matrix.
Returns
the rank. For an empty matrix, returns 0.

◆ get_rank() [2/2]

template<typename value_type >
size_t molpro::linalg::itsolv::get_rank ( std::vector< value_type >  eigenvalues,
value_type  threshold 
)

Get the rank of some matrix, given a threshold.

Parameters
[in]eigenvaluesthe matrix, as a vector.
[in]thresholdthe threshold. Note that this is the normalised threshold, a value between 0 and 1, relative to the largest element in the matrix.
Returns
the rank. For an empty matrix, returns 0.

◆ get_rank< value_type >()

template size_t molpro::linalg::itsolv::get_rank< value_type > ( std::vector< value_type >  eigenvalues,
value_type  threshold 
)

◆ operator<<() [1/3]

std::ostream & molpro::linalg::itsolv::operator<< ( std::ostream &  o,
const Statistics statistics 
)
inline

◆ operator<<() [2/3]

std::ostream & molpro::linalg::itsolv::operator<< ( std::ostream &  os,
const Interpolate interpolant 
)

◆ operator<<() [3/3]

std::ostream & molpro::linalg::itsolv::operator<< ( std::ostream &  os,
const Interpolate::point p 
)

◆ operator==()

bool molpro::linalg::itsolv::operator== ( const Interpolate::point lhs,
const Interpolate::point rhs 
)
inline

◆ precondition_default() [1/3]

template<typename T , typename = std::enable_if_t<std::is_base_of<molpro::linalg::array::DistrArray, T>::value>>
void molpro::linalg::itsolv::precondition_default ( const VecRef< T > &  action,
const std::vector< double > &  shift,
const T &  diagonals 
)

◆ precondition_default() [2/3]

template<typename T , typename = std::enable_if_t<!std::is_base_of<molpro::linalg::array::DistrArray, T>::value>, class = void>
void molpro::linalg::itsolv::precondition_default ( const VecRef< T > &  action,
const std::vector< double > &  shift,
const T &  diagonals,
typename std::enable_if<!has_iterator< T >::value, void * >::type  = nullptr 
)

◆ precondition_default() [3/3]

template<class T >
void molpro::linalg::itsolv::precondition_default ( const VecRef< T > &  action,
const std::vector< double > &  shift,
const T &  diagonals,
typename T::iterator *  = nullptr 
)

◆ printMatrix()

template<typename value_type >
void molpro::linalg::itsolv::printMatrix ( const std::vector< value_type > &  m,
size_t  rows,
size_t  cols,
std::string  title = "",
std::ostream &  s = molpro::cout 
)

◆ printMatrix< double >()

template void molpro::linalg::itsolv::printMatrix< double > ( const std::vector< double > &  ,
size_t  rows,
size_t  cols,
std::string  title,
std::ostream &  s 
)

◆ printMatrix< std::complex< double > >()

template void molpro::linalg::itsolv::printMatrix< std::complex< double > > ( const std::vector< std::complex< double > > &  ,
size_t  rows,
size_t  cols,
std::string  title,
std::ostream &  s 
)

◆ printMatrix< value_type >()

template void molpro::linalg::itsolv::printMatrix< value_type > ( const std::vector< value_type > &  ,
size_t  rows,
size_t  cols,
std::string  title,
std::ostream &  s 
)

◆ read_handler_counts()

template<typename R , typename Q , typename P >
void molpro::linalg::itsolv::read_handler_counts ( std::shared_ptr< Statistics stats,
std::shared_ptr< ArrayHandlers< R, Q, P > >  handlers 
)

◆ remove_elements()

template<typename T , typename U >
auto molpro::linalg::itsolv::remove_elements ( std::vector< T >  params,
const std::vector< U > &  indices 
)

Removes indices from a vector.

Template Parameters
Tvalue type
Iindex type
Parameters
paramscontainer
indicesindices to remove
Returns

◆ solve_DIIS()

template<typename value_type , typename std::enable_if_t< is_complex< value_type >{}, int > = 0>
void molpro::linalg::itsolv::solve_DIIS ( std::vector< value_type > &  solution,
const std::vector< value_type > &  matrix,
size_t  dimension,
double  svdThreshold,
int  verbosity = 0 
)

◆ solve_DIIS< double >()

template void molpro::linalg::itsolv::solve_DIIS< double > ( std::vector< double > &  solution,
const std::vector< double > &  matrix,
const size_t  dimension,
double  svdThreshold,
int  verbosity 
)

◆ solve_DIIS< std::complex< double > >()

template void molpro::linalg::itsolv::solve_DIIS< std::complex< double > > ( std::vector< std::complex< double > > &  solution,
const std::vector< std::complex< double > > &  matrix,
const size_t  dimension,
double  svdThreshold,
int  verbosity 
)

◆ solve_DIIS< value_type >()

template void molpro::linalg::itsolv::solve_DIIS< value_type > ( std::vector< value_type > &  solution,
const std::vector< value_type > &  matrix,
size_t  dimension,
double  svdThreshold,
int  verbosity 
)

◆ Solve_LinearEquations() [1/3]

template<class R , class Q = R, class P = std::map<size_t, typename R::value_type>>
std::tuple< bool, std::unique_ptr< LinearEquations< R, Q, P > > > molpro::linalg::itsolv::Solve_LinearEquations ( const VecRef< R > &  parameters,
const VecRef< R > &  actions,
const Problem< R > &  problem,
int  verbosity = 0,
bool  generate_initial_guess = true,
const std::string &  method = "Davidson",
const std::string &  options = "",
const std::shared_ptr< ArrayHandlers< R, Q, P > > &  handlers = std::make_shared<molpro::linalg::itsolv::ArrayHandlers<R, Q, P>>() 
)

◆ Solve_LinearEquations() [2/3]

template<class R , class Q = R, class P = std::map<size_t, typename R::value_type>>
auto molpro::linalg::itsolv::Solve_LinearEquations ( R &  parameters,
R &  actions,
const Problem< R > &  problem,
int  verbosity = 0,
bool  generate_initial_guess = true,
const std::string &  method = "Davidson",
const std::string &  options = "",
const std::shared_ptr< ArrayHandlers< R, Q, P > > &  handlers = std::make_shared<molpro::linalg::itsolv::ArrayHandlers<R, Q, P>>() 
)

◆ Solve_LinearEquations() [3/3]

template<class R , class Q = R, class P = std::map<size_t, typename R::value_type>>
std::tuple< bool, std::unique_ptr< LinearEquations< R, Q, P > > > molpro::linalg::itsolv::Solve_LinearEquations ( std::vector< R > &  parameters,
std::vector< R > &  actions,
const Problem< R > &  problem,
int  verbosity = 0,
bool  generate_initial_guess = true,
const std::string &  method = "Davidson",
const std::string &  options = "",
const std::shared_ptr< ArrayHandlers< R, Q, P > > &  handlers = std::make_shared<molpro::linalg::itsolv::ArrayHandlers<R, Q, P>>() 
)

◆ solve_LinearEquations()

template<typename value_type , typename std::enable_if_t< is_complex< value_type >{}, int > = 0>
void molpro::linalg::itsolv::solve_LinearEquations ( std::vector< value_type > &  solution,
std::vector< value_type > &  eigenvalues,
const std::vector< value_type > &  matrix,
const std::vector< value_type > &  metric,
const std::vector< value_type > &  rhs,
size_t  dimension,
size_t  nroot,
double  augmented_hessian,
double  svdThreshold,
int  verbosity 
)

◆ solve_LinearEquations< double >()

template void molpro::linalg::itsolv::solve_LinearEquations< double > ( std::vector< double > &  solution,
std::vector< double > &  eigenvalues,
const std::vector< double > &  matrix,
const std::vector< double > &  metric,
const std::vector< double > &  rhs,
size_t  dimension,
size_t  nroot,
double  augmented_hessian,
double  svdThreshold,
int  verbosity 
)

◆ solve_LinearEquations< std::complex< double > >()

template void molpro::linalg::itsolv::solve_LinearEquations< std::complex< double > > ( std::vector< std::complex< double > > &  solution,
std::vector< std::complex< double > > &  eigenvalues,
const std::vector< std::complex< double > > &  matrix,
const std::vector< std::complex< double > > &  metric,
const std::vector< std::complex< double > > &  rhs,
size_t  dimension,
size_t  nroot,
double  augmented_hessian,
double  svdThreshold,
int  verbosity 
)

◆ solve_LinearEquations< value_type >()

template void molpro::linalg::itsolv::solve_LinearEquations< value_type > ( std::vector< value_type > &  solution,
std::vector< value_type > &  eigenvalues,
const std::vector< value_type > &  matrix,
const std::vector< value_type > &  metric,
const std::vector< value_type > &  rhs,
size_t  dimension,
size_t  nroot,
double  augmented_hessian,
double  svdThreshold,
int  verbosity 
)

◆ svd_eigen_bdcsvd()

template<typename value_type >
std::list< SVD< value_type > > molpro::linalg::itsolv::svd_eigen_bdcsvd ( size_t  nrows,
size_t  ncols,
const array::Span< value_type > &  m,
double  threshold 
)

◆ svd_eigen_jacobi()

template<typename value_type >
std::list< SVD< value_type > > molpro::linalg::itsolv::svd_eigen_jacobi ( size_t  nrows,
size_t  ncols,
const array::Span< value_type > &  m,
double  threshold 
)

◆ svd_system() [1/3]

template std::list< SVD< double > > molpro::linalg::itsolv::svd_system ( size_t  nrows,
size_t  ncols,
const array::Span< double > &  m,
double  threshold,
bool  hermitian,
bool  reduce_to_rank 
)

◆ svd_system() [2/3]

template std::list< SVD< std::complex< double > > > molpro::linalg::itsolv::svd_system ( size_t  nrows,
size_t  ncols,
const array::Span< std::complex< double > > &  m,
double  threshold,
bool  hermitian,
bool  reduce_to_rank 
)

◆ svd_system() [3/3]

template<typename value_type , typename std::enable_if_t<!is_complex< value_type >{}, std::nullptr_t > = nullptr>
std::list< SVD< value_type > > molpro::linalg::itsolv::svd_system ( size_t  nrows,
size_t  ncols,
const array::Span< value_type > &  m,
double  threshold,
bool  hermitian = false,
bool  reduce_to_rank = false 
)

Performs singular value decomposition and returns SVD objects for singular values less than threshold, sorted in ascending order.

Template Parameters
value_type
Parameters
nrowsnumber of rows in the matrix
ncolsnumber of columns in the matrix
mrow-wise data buffer for a matrix
thresholdsingular values less than threshold will be returned

◆ svd_system< value_type >()

template std::list< SVD< value_type > > molpro::linalg::itsolv::svd_system< value_type > ( size_t  nrows,
size_t  ncols,
const array::Span< value_type > &  m,
double  threshold,
bool  hermitian,
bool  reduce_to_rank 
)

◆ wrap() [1/3]

template<class IterableContainer >
auto molpro::linalg::itsolv::wrap ( const IterableContainer &  parameters)

Wraps references for each parameter in an iterable container that implements begin()/end()

Template Parameters
IterableContainershould have begin()/end() either as free or member functions
Relement type
Parameters
parametersparameters to wrap
Returns
vector of constant references to each element in parameters

◆ wrap() [2/3]

template<class ForwardIt >
auto molpro::linalg::itsolv::wrap ( ForwardIt  begin,
ForwardIt  end 
)

Takes a begin and end iterators and returns a vector of references to each element.

◆ wrap() [3/3]

template<class IterableContainer >
auto molpro::linalg::itsolv::wrap ( IterableContainer &  parameters)

Wraps references for each parameter in an iterable container that implements begin()/end()

Template Parameters
IterableContainershould have begin()/end() either as free or member functions
Relement type
Parameters
parametersparameters to wrap
Returns
vector of references to each element in parameters

◆ wrap_arg()

template<typename T , typename... S>
auto molpro::linalg::itsolv::wrap_arg ( T &&  arg,
S &&...  args 
) -> std::enable_if_t<std::conjunction_v<std::is_same<decay_t<T>, decay_t<S>>...>, VecRef<decay_t<T>>>

Constructs a vector of reference wrappers with provided arguments.