iterative-solver 0.0
LinearEquationsDavidson.h
1#ifndef LINEARALGEBRA_SRC_MOLPRO_LINALG_ITSOLV_LINEAREQUATIONSDAVIDSON_H
2#define LINEARALGEBRA_SRC_MOLPRO_LINALG_ITSOLV_LINEAREQUATIONSDAVIDSON_H
3#include <molpro/linalg/itsolv/CastOptions.h>
4#include <molpro/linalg/itsolv/DSpaceResetter.h>
5#include <molpro/linalg/itsolv/IterativeSolverTemplate.h>
6#include <molpro/linalg/itsolv/propose_rspace.h>
7
27template <class R, class Q = R, class P = std::map<size_t, typename R::value_type>>
28class LinearEquationsDavidson : public IterativeSolverTemplate<LinearEquations, R, Q, P> {
29public:
32 using typename SolverTemplate::value_type;
33 using typename SolverTemplate::value_type_abs;
34
35 explicit LinearEquationsDavidson(const std::shared_ptr<ArrayHandlers<R, Q, P>>& handlers,
36 const std::shared_ptr<Logger>& logger_ = std::make_shared<Logger>())
37 : SolverTemplate(std::make_shared<subspace::XSpace<R, Q, P>>(handlers, logger_),
38 std::static_pointer_cast<subspace::ISubspaceSolver<R, Q, P>>(
39 std::make_shared<subspace::SubspaceSolverLinEig<R, Q, P>>(logger_)),
40 handlers, std::make_shared<Statistics>(), logger_),
41 logger(logger_) {
43 this->m_normalise_solution = false;
44 }
45
46 bool solve(const VecRef<R>& parameters, const VecRef<R>& actions, const Problem<R>& problem,
47 bool generate_initial_guess = false) override {
48 R& vector = actions[0];
49 for (unsigned int instance = 0; problem.RHS(vector, instance); ++instance)
50 add_equations(vector);
51 return IterativeSolverTemplate<LinearEquations, R, Q, P>::solve(parameters, actions, problem,
52 generate_initial_guess);
53 }
54
55 bool nonlinear() const override { return false; }
56
57 size_t end_iteration(const VecRef<R>& parameters, const VecRef<R>& action) override {
58 auto prof = this->profiler()->push("itsolv::end_iteration");
59 if (m_dspace_resetter.do_reset(this->m_stats->iterations, this->m_xspace->dimensions())) {
60 this->m_working_set = m_dspace_resetter.run(parameters, *this->m_xspace, this->m_subspace_solver->solutions(),
61 m_norm_thresh, m_svd_thresh, *this->m_handlers, *this->m_logger);
62 } else {
63 this->m_working_set = detail::propose_rspace(*this, parameters, action, *this->m_xspace, *this->m_subspace_solver,
65 m_max_size_qspace, *this->profiler());
66 }
67 this->m_stats->iterations++;
68 this->m_end_iteration_needed = false;
69 return this->working_set().size();
70 }
71
72 // FIXME move this to the template
73 size_t end_iteration(std::vector<R>& parameters, std::vector<R>& action) override {
74 return end_iteration(wrap(parameters), wrap(action));
75 }
76 size_t end_iteration(R& parameters, R& actions) override {
77 auto wparams = std::vector<std::reference_wrapper<R>>{std::ref(parameters)};
78 auto wactions = std::vector<std::reference_wrapper<R>>{std::ref(actions)};
79 return end_iteration(wparams, wactions);
80 }
81
82 void add_equations(const CVecRef<R>& rhs) override {
83 auto prof = this->profiler()->push("itsolv::add_equations");
84 auto xspace = std::static_pointer_cast<subspace::XSpace<R, Q, P>>(this->m_xspace);
85 xspace->add_rhs_equations(rhs);
86 this->set_n_roots(xspace->dimensions().nRHS);
87 }
88
89 void add_equations(const R& rhs) override { add_equations(cwrap_arg(rhs)); }
90 void add_equations(const std::vector<R>& rhs) override { add_equations(cwrap(rhs)); }
91
92 CVecRef<Q> rhs() const override {
93 auto xspace = std::static_pointer_cast<subspace::XSpace<R, Q, P>>(this->m_xspace);
94 return xspace->rhs();
95 }
96
98 void set_norm_thresh(double thresh) { m_norm_thresh = thresh; }
99 double get_norm_thresh() const { return m_norm_thresh; }
102 void set_svd_thresh(double thresh) { m_svd_thresh = thresh; }
103 double get_svd_thresh() const { return m_svd_thresh; }
105 void set_reset_D(size_t n) { m_dspace_resetter.set_nreset(n); }
106 size_t get_reset_D() const { return m_dspace_resetter.get_nreset(); }
108 void set_reset_D_maxQ_size(size_t n) { m_dspace_resetter.set_max_Qsize(n); }
109 int get_reset_D_maxQ_size() const { return m_dspace_resetter.get_max_Qsize(); }
114 if (m_dspace_resetter.get_max_Qsize() > m_max_size_qspace)
116 }
118 void set_hermiticity(bool hermitian) override {
119 m_hermiticity = hermitian;
120 auto xspace = std::dynamic_pointer_cast<subspace::XSpace<R, Q, P>>(this->m_xspace);
121 xspace->set_hermiticity(hermitian);
122 auto subspace_solver = std::dynamic_pointer_cast<subspace::SubspaceSolverLinEig<R, Q, P>>(this->m_subspace_solver);
123 subspace_solver->set_hermiticity(hermitian);
124 }
125 bool get_hermiticity() const override { return m_hermiticity; }
127 void set_augmented_hessian(const double parameter) {
128 auto subspace_solver = std::dynamic_pointer_cast<subspace::SubspaceSolverLinEig<R, Q, P>>(this->m_subspace_solver);
129 subspace_solver->set_augmented_hessian(parameter);
130 }
132 double get_augmented_hessian() const {
133 auto subspace_solver = std::dynamic_pointer_cast<subspace::SubspaceSolverLinEig<R, Q, P>>(this->m_subspace_solver);
134 return subspace_solver->get_augmented_hessian();
135 }
136
137 void set_options(const Options& options) override {
140 if (opt.reset_D)
141 set_reset_D(opt.reset_D.value());
142 if (opt.reset_D_max_Q_size)
143 set_reset_D_maxQ_size(opt.reset_D_max_Q_size.value());
144 if (opt.max_size_qspace)
145 set_max_size_qspace(opt.max_size_qspace.value());
146 if (opt.norm_thresh)
147 set_norm_thresh(opt.norm_thresh.value());
148 if (opt.svd_thresh)
149 set_svd_thresh(opt.svd_thresh.value());
150 if (opt.hermiticity)
151 set_hermiticity(opt.hermiticity.value());
152 if (opt.augmented_hessian)
153 set_augmented_hessian(opt.augmented_hessian.value());
154 }
155
156 std::shared_ptr<Options> get_options() const override {
157 auto opt = std::make_shared<LinearEquationsDavidsonOptions>();
158 opt->copy(*SolverTemplate::get_options());
159 opt->reset_D = get_reset_D();
160 opt->reset_D_max_Q_size = get_reset_D_maxQ_size();
161 opt->max_size_qspace = get_max_size_qspace();
162 opt->norm_thresh = get_norm_thresh();
163 opt->svd_thresh = get_svd_thresh();
164 opt->hermiticity = get_hermiticity();
165 opt->augmented_hessian = get_augmented_hessian();
166 return opt;
167 }
168
169 void report(std::ostream& cout, bool endl = true) const override {
170 SolverTemplate::report(cout, false);
171 cout << ", errors " << std::scientific;
172 auto& err = this->m_errors;
173 std::copy(begin(err), end(err), std::ostream_iterator<value_type_abs>(molpro::cout, ", "));
174 cout << std::defaultfloat;
175 if (endl)
176 cout << std::endl;
177 }
178 std::shared_ptr<Logger> logger;
179
180protected:
181 // FIXME The scale is fixed by the norm of RHS, but if RHS=0 there is no reference. We could use the norm of params
182 void construct_residual(const std::vector<int>& roots, const CVecRef<R>& params, const VecRef<R>& actions) override {
183 assert(params.size() >= roots.size());
184 const auto& norm = std::dynamic_pointer_cast<subspace::XSpace<R, Q, P>>(this->m_xspace)->rhs_norm();
185 for (size_t i = 0; i < roots.size(); ++i) {
186 const auto ii = roots[i];
187 this->m_handlers->rq().axpy(-1, rhs().at(ii), actions.at(i));
188 if (norm.at(ii) != 0) {
189 auto scal = 1 / norm[ii];
190 this->m_handlers->rr().scal(scal, actions.at(i));
191 }
192 }
193 }
194
195 double m_norm_thresh = 1e-10;
196 double m_svd_thresh = 1e-12;
197 int m_max_size_qspace = std::numeric_limits<int>::max();
199 bool m_hermiticity = true;
200};
201
202} // namespace molpro::linalg::itsolv
203#endif // LINEARALGEBRA_SRC_MOLPRO_LINALG_ITSOLV_LINEAREQUATIONSDAVIDSON_H
Class, containing a collection of array handlers used in IterativeSolver Provides a Builder sub-class...
Definition: ArrayHandlers.h:25
Implements IterativeSolver interface that is common to all solvers.
Definition: IterativeSolverTemplate.h:127
std::vector< int > m_working_set
indices of roots in the working set
Definition: IterativeSolverTemplate.h:583
std::shared_ptr< subspace::IXSpace< R, R, std::map< size_t, typename R::value_type > > > m_xspace
manages the subspace and associated data
Definition: IterativeSolverTemplate.h:579
std::shared_ptr< ArrayHandlers< R, R, std::map< size_t, typename R::value_type > > > m_handlers
Array handlers.
Definition: IterativeSolverTemplate.h:578
std::shared_ptr< subspace::ISubspaceSolver< R, R, std::map< size_t, typename R::value_type > > > m_subspace_solver
solves the subspace problem
Definition: IterativeSolverTemplate.h:580
std::shared_ptr< Options > get_options() const override
Definition: IterativeSolverTemplate.h:262
const std::shared_ptr< molpro::profiler::Profiler > & profiler() const override
Definition: IterativeSolverTemplate.h:320
bool m_end_iteration_needed
whether end_iteration should be called after any preconditioner
Definition: IterativeSolverTemplate.h:600
void report() const override
Definition: IterativeSolverTemplate.h:287
bool m_normalise_solution
whether to normalise the solutions
Definition: IterativeSolverTemplate.h:590
const std::vector< int > & working_set() const override
Working set of roots that are not yet converged.
Definition: IterativeSolverTemplate.h:243
void set_options(const Options &options) override
Definition: IterativeSolverTemplate.h:253
std::vector< double > m_errors
errors from the most recent solution
Definition: IterativeSolverTemplate.h:581
bool solve(const VecRef< R > &parameters, const VecRef< R > &actions, const Problem< R > &problem, bool generate_initial_guess=false) override
Definition: IterativeSolverTemplate.h:322
std::shared_ptr< Statistics > m_stats
accumulates statistics of operations performed by the solver
Definition: IterativeSolverTemplate.h:588
Solves a system of linear equation, A x = b.
Definition: LinearEquationsDavidson.h:28
size_t get_reset_D() const
Definition: LinearEquationsDavidson.h:106
detail::DSpaceResetter< Q > m_dspace_resetter
resets D space
Definition: LinearEquationsDavidson.h:198
bool solve(const VecRef< R > &parameters, const VecRef< R > &actions, const Problem< R > &problem, bool generate_initial_guess=false) override
Simplified one-call solver.
Definition: LinearEquationsDavidson.h:46
std::shared_ptr< Logger > logger
Definition: LinearEquationsDavidson.h:178
double m_norm_thresh
vectors with norm less than threshold can be considered null.
Definition: LinearEquationsDavidson.h:195
void add_equations(const std::vector< R > &rhs) override
Definition: LinearEquationsDavidson.h:90
void set_hermiticity(bool hermitian) override
Sets hermiticity of kernel.
Definition: LinearEquationsDavidson.h:118
void set_augmented_hessian(const double parameter)
Set value of augmented hessian parameter. If 0, than augmented Hessian is not used.
Definition: LinearEquationsDavidson.h:127
void set_svd_thresh(double thresh)
Definition: LinearEquationsDavidson.h:102
size_t end_iteration(const VecRef< R > &parameters, const VecRef< R > &action) override
Behaviour depends on the solver.
Definition: LinearEquationsDavidson.h:57
double get_norm_thresh() const
Definition: LinearEquationsDavidson.h:99
double m_svd_thresh
svd values smaller than this mark the null space
Definition: LinearEquationsDavidson.h:196
size_t end_iteration(R &parameters, R &actions) override
Definition: LinearEquationsDavidson.h:76
void add_equations(const CVecRef< R > &rhs) override
Definition: LinearEquationsDavidson.h:82
std::shared_ptr< Options > get_options() const override
Definition: LinearEquationsDavidson.h:156
LinearEquationsDavidson(const std::shared_ptr< ArrayHandlers< R, Q, P > > &handlers, const std::shared_ptr< Logger > &logger_=std::make_shared< Logger >())
Definition: LinearEquationsDavidson.h:35
void set_norm_thresh(double thresh)
Set threshold on the norm of parameters that should be considered null.
Definition: LinearEquationsDavidson.h:98
void add_equations(const R &rhs) override
Definition: LinearEquationsDavidson.h:89
int get_reset_D_maxQ_size() const
Definition: LinearEquationsDavidson.h:109
bool nonlinear() const override
Report whether the class is a non-linear solver.
Definition: LinearEquationsDavidson.h:55
size_t end_iteration(std::vector< R > &parameters, std::vector< R > &action) override
Definition: LinearEquationsDavidson.h:73
void set_options(const Options &options) override
Definition: LinearEquationsDavidson.h:137
bool m_hermiticity
whether the problem is hermitian or not
Definition: LinearEquationsDavidson.h:199
int m_max_size_qspace
maximum size of Q space
Definition: LinearEquationsDavidson.h:197
void construct_residual(const std::vector< int > &roots, const CVecRef< R > &params, const VecRef< R > &actions) override
Constructs residual for given roots provided their parameters and actions.
Definition: LinearEquationsDavidson.h:182
double get_augmented_hessian() const
Definition: LinearEquationsDavidson.h:132
CVecRef< Q > rhs() const override
Definition: LinearEquationsDavidson.h:92
void set_reset_D(size_t n)
Set the period in iterations for resetting the D space.
Definition: LinearEquationsDavidson.h:105
void report(std::ostream &cout, bool endl=true) const override
Writes a report to cout output stream.
Definition: LinearEquationsDavidson.h:169
bool get_hermiticity() const override
Gets hermiticity of kernel, if true than it is hermitian, otherwise it is not.
Definition: LinearEquationsDavidson.h:125
void report() const override
Writes a report to std::cout.
Definition: IterativeSolverTemplate.h:287
double get_svd_thresh() const
Definition: LinearEquationsDavidson.h:103
void set_reset_D_maxQ_size(size_t n)
Set the maximum size of Q space after resetting the D space.
Definition: LinearEquationsDavidson.h:108
int get_max_size_qspace() const
Definition: LinearEquationsDavidson.h:117
void set_max_size_qspace(int n)
Definition: LinearEquationsDavidson.h:112
Abstract class defining the problem-specific interface for the simplified solver interface to Iterati...
Definition: IterativeSolver.h:77
virtual bool RHS(R &RHS, unsigned int instance) const
Return the inhomogeneous part of a linear equation system.
Definition: IterativeSolver.h:146
Resets D space constructing full solutions as the new working set, removing instabilities from Q spac...
Definition: DSpaceResetter.h:70
auto propose_rspace(IterativeSolver< R, Q, P > &solver, const VecRef< R > &parameters, const VecRef< R > &residuals, subspace::IXSpace< R, Q, P > &xspace, subspace::ISubspaceSolver< R, Q, P > &subspace_solver, ArrayHandlers< R, Q, P > &handlers, Logger &logger, value_type_abs svd_thresh, value_type_abs res_norm_thresh, int max_size_qspace, molpro::profiler::Profiler &profiler)
Proposes new parameters for the subspace from the preconditioned residuals.
Definition: propose_rspace.h:554
4-parameter interpolation of a 1-dimensional function given two points for which function values and ...
Definition: helper.h:10
auto cwrap(ForwardIt begin, ForwardIt end)
Takes a begin and end iterators and returns a vector of references to each element.
Definition: wrap.h:52
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.
Definition: wrap.h:149
auto wrap(ForwardIt begin, ForwardIt end)
Takes a begin and end iterators and returns a vector of references to each element.
Definition: wrap.h:32
std::vector< std::reference_wrapper< const A > > CVecRef
Definition: wrap.h:14
std::vector< std::reference_wrapper< A > > VecRef
Definition: wrap.h:11
const std::shared_ptr< const molpro::Options > options()
Get the Options object associated with iterative-solver.
Definition: options.cpp:4
static std::shared_ptr< LinearEquationsDavidsonOptions > LinearEquations(const std::shared_ptr< Options > &options)
Definition: CastOptions.h:54
Access point for different options in iterative solvers.
Definition: Options.h:20
Information about performance of IterativeSolver instance.
Definition: Statistics.h:10