iterative-solver 0.0
IterativeSolverTemplate.h
1#ifndef LINEARALGEBRA_SRC_MOLPRO_LINALG_ITSOLV_ITERATIVESOLVERTEMPLATE_H
2#define LINEARALGEBRA_SRC_MOLPRO_LINALG_ITSOLV_ITERATIVESOLVERTEMPLATE_H
3#include <molpro/Profiler.h>
4#include <molpro/iostream.h>
5#include <molpro/linalg/itsolv/IterativeSolver.h>
6#include <molpro/linalg/itsolv/Logger.h>
7#include <molpro/linalg/itsolv/subspace/ISubspaceSolver.h>
8#include <molpro/linalg/itsolv/subspace/IXSpace.h>
9#include <molpro/linalg/itsolv/subspace/Matrix.h>
10#include <molpro/linalg/itsolv/subspace/util.h>
11#include <molpro/linalg/itsolv/util.h>
12#include <molpro/linalg/itsolv/wrap.h>
13#include <molpro/profiler/Profiler.h>
14
15#include <cassert>
16#include <cmath>
17#include <format>
18#include <iostream>
19#include <optional>
20#include <sstream>
21
22namespace molpro::linalg::itsolv {
23namespace detail {
24
25inline std::vector<std::pair<size_t, size_t>> parameter_batches(const size_t nsol, const size_t nparam) {
26 auto batches = std::vector<std::pair<size_t, size_t>>{};
27 if (nparam && nsol) {
28 auto n_batch = nsol / nparam + (nsol % nparam ? 1 : 0);
29 for (size_t ib = 0, start_sol = 0, end_sol = 0; ib < n_batch; ++ib, start_sol = end_sol) {
30 end_sol = std::min(start_sol + nparam, nsol);
31 batches.emplace_back(start_sol, end_sol);
32 }
33 }
34 return batches;
35}
36
37template <class R, class Q, class P>
38void construct_solution(const VecRef<R>& params, const std::vector<int>& roots,
39 const subspace::Matrix<double>& solutions,
40 const std::vector<std::reference_wrapper<P>>& pparams,
41 const std::vector<std::reference_wrapper<Q>>& qparams,
42 const std::vector<std::reference_wrapper<Q>>& dparams, size_t oP, size_t oQ, size_t oD,
43 ArrayHandlers<R, Q, P>& handlers) {
44 auto prof = molpro::Profiler::single();
45 prof->start("get rd_mat");
46 if (roots.empty())
47 return;
48 assert(params.size() >= roots.size());
49 for (size_t i = 0; i < roots.size(); ++i) {
50 handlers.rr().fill(0, params.at(i));
51 }
52 subspace::Matrix<double> rp_mat(std::make_pair(pparams.size(), roots.size())),
53 rq_mat(std::make_pair(qparams.size(), roots.size())), rd_mat(std::make_pair(dparams.size(), roots.size()));
54 for (size_t i = 0; i < roots.size(); ++i) {
55 for (size_t j = 0; j < pparams.size(); ++j) {
56 rp_mat(j, i) = solutions(roots[i], oP + j);
57 }
58 for (size_t j = 0; j < qparams.size(); ++j) {
59 rq_mat(j, i) = solutions(roots[i], oQ + j);
60 }
61 for (size_t j = 0; j < dparams.size(); ++j) {
62 rd_mat(j, i) = solutions(roots[i], oD + j);
63 }
64 }
65 prof->stop();
66 handlers.rp().gemm_outer(rp_mat, cwrap(pparams), params);
67 handlers.rq().gemm_outer(rq_mat, cwrap(qparams), params);
68 handlers.rq().gemm_outer(rd_mat, cwrap(dparams), params);
69}
70
71template <typename T>
72std::vector<std::vector<T>> construct_vectorP(const std::vector<int>& roots, const subspace::Matrix<T>& solutions,
73 const size_t oP, const size_t nP) {
74 auto vectorP = std::vector<std::vector<T>>{};
75 for (auto root : roots) {
76 vectorP.emplace_back();
77 for (size_t j = 0; j < nP; ++j)
78 vectorP.back().push_back(solutions(root, oP + j));
79 }
80 return vectorP;
81}
82
83template <class R>
84void normalise(const size_t n_roots, const VecRef<R>& params, const VecRef<R>& actions,
85 array::ArrayHandler<R, R>& handler, Logger& logger) {
86 assert(params.size() >= n_roots && actions.size() >= n_roots);
87 for (size_t i = 0; i < n_roots; ++i) {
88 auto dot = handler.dot(params.at(i), params.at(i));
89 dot = std::sqrt(std::abs(dot));
90 if (dot > 1.0e-14) {
91 handler.scal(1. / dot, params.at(i));
92 handler.scal(1. / dot, actions.at(i));
93 } else {
94 logger.warn("solution parameter's length is too small, dot = " + std::format("{:.2e}", dot));
95 }
96 }
97}
98
99template <class R, typename T>
100void update_errors(std::vector<T>& errors, const CVecRef<R>& residual, array::ArrayHandler<R, R>& handler) {
101 assert(residual.size() >= errors.size());
102 for (size_t i = 0; i < errors.size(); ++i) {
103 auto a = handler.dot(residual[i], residual[i]);
104 errors[i] = std::sqrt(std::abs(a));
105 }
106}
107
108template <typename T>
109std::vector<int> select_working_set(const size_t nw, const std::vector<T>& errors, const T threshold,
110 const std::vector<T>& value_errors, const T value_threshold) {
111 auto ordered_errors = std::multimap<T, size_t, std::greater<T>>{};
112 for (size_t i = 0; i < errors.size(); ++i) {
113 if (errors[i] > threshold or (i < value_errors.size() and value_errors[i] > value_threshold))
114 ordered_errors.emplace(errors[i], i);
115 }
116 auto working_set = std::vector<int>{};
117 auto end = (ordered_errors.size() < nw ? ordered_errors.end() : next(begin(ordered_errors), nw));
118 std::transform(begin(ordered_errors), end, std::back_inserter(working_set), [](const auto& el) { return el.second; });
119 std::sort(working_set.begin(), working_set.end());
120 return working_set;
121}
122
123} // namespace detail
124
125namespace log {
126
130struct NewIteration : ContextBase<NewIteration, true, int, std::vector<double>> {
131 static constexpr const char *name = "NewIteration";
132 static constexpr std::size_t iter = 0;
133 static constexpr std::size_t errors = 1;
134};
135static_assert(context<NewIteration>);
136
140struct IterationReport : ContextBase<IterationReport, true> {
141 static constexpr const char *name = "IterationReport";
142};
143static_assert(context<IterationReport>);
144
148struct SummaryReport : ContextBase<SummaryReport, true> {
149 static constexpr const char *name = "SummaryReport";
150};
151static_assert(context<SummaryReport>);
152
153}
154
160template <template <class, class, class> class Solver, class R, class Q, class P>
161class IterativeSolverTemplate : public Solver<R, Q, P> {
162public:
163 using typename Solver<R, Q, P>::fapply_on_p_type;
164 using typename Solver<R, Q, P>::scalar_type;
165 using typename Solver<R, Q, P>::value_type;
166 using typename Solver<R, Q, P>::VectorP;
167
171 IterativeSolverTemplate<Solver, R, Q, P>& operator=(const IterativeSolverTemplate<Solver, R, Q, P>&) = delete;
172 IterativeSolverTemplate<Solver, R, Q, P>& operator=(IterativeSolverTemplate<Solver, R, Q, P>&&) noexcept = default;
173
174 void set_logger(std::shared_ptr<Logger> logger) override {
175 assert(logger);
176 m_logger = std::move(logger);
177 m_xspace->set_logger(m_logger);
178 m_subspace_solver->set_logger(m_logger);
179 }
180
181 Logger &logger() override { return *m_logger; }
182
183 int add_vector(const VecRef<R>& parameters, const VecRef<R>& actions) override {
184 auto outer = profiler()->push("itsolv::add_vector");
185 auto prof = molpro::Profiler::single();
186 m_logger->trace("IterativeSolverTemplate::add_vector iteration = ", m_stats->iterations);
187 m_logger->debug("IterativeSolverTemplate::add_vector size of {params, actions, working_set} = ",
188 parameters.size(), actions.size(), m_working_set.size());
189 if (m_xspace->dimensions().nP != 0 && !m_apply_p)
190 throw std::runtime_error(
191 "Solver contains P space but no valid apply_p function. Make sure add_p was called correctly.");
192 auto nW = std::min(m_working_set.size(), parameters.size());
193 auto cwparams = cwrap(begin(parameters), begin(parameters) + nW);
194 auto cwactions = cwrap(begin(actions), begin(actions) + nW);
195 m_stats->r_creations += nW;
196 prof->start("update_qspace");
197 m_xspace->update_qspace(cwparams, cwactions);
198 m_stats->q_creations += 2 * nW;
199 prof->stop();
200 prof->start("solve_and_generate_working_set");
201 auto working_set = solve_and_generate_working_set(parameters, actions);
202 prof->stop();
204 this->m_end_iteration_needed = true;
205 return working_set;
206 }
207
208 int add_vector(std::vector<R>& parameters, std::vector<R>& actions) override {
209 return add_vector(wrap(parameters), wrap(actions));
210 }
211 int add_vector(R& parameters, R& actions, value_type value = 0) override {
212 return add_vector(wrap_arg(parameters), wrap_arg(actions));
213 }
214
215 // FIXME Currently only works if called on an empty subspace. Either enforce it or generalise.
216 size_t add_p(const CVecRef<P>& pparams, const array::Span<value_type>& pp_action_matrix, const VecRef<R>& parameters,
217 const VecRef<R>& actions, fapply_on_p_type apply_p) override {
218 auto prof = profiler()->push("itsolv::add_p");
219 if (not pparams.empty() and pparams.size() < n_roots())
220 throw std::runtime_error("P space must be empty or at least as large as number of roots sought");
221 if (apply_p)
222 m_apply_p = std::move(apply_p);
223 m_xspace->update_pspace(pparams, pp_action_matrix);
224 auto working_set = solve_and_generate_working_set(parameters, actions);
226 return working_set;
227 };
228
229 void clearP() override {}
230
231 void solution(const std::vector<int>& roots, const VecRef<R>& parameters, const VecRef<R>& residual) override {
232 // auto prof = profiler()->push("itsolv::solution"); // FIXME two profilers
233 auto prof = molpro::Profiler::single();
234 check_consistent_number_of_roots_and_solutions(roots, parameters.size());
235 prof->start("construct_solution (parameters)");
236 detail::construct_solution(parameters, roots, m_subspace_solver->solutions(), m_xspace->paramsp(),
237 m_xspace->paramsq(), m_xspace->paramsd(), m_xspace->dimensions().oP,
238 m_xspace->dimensions().oQ, m_xspace->dimensions().oD, *m_handlers);
239 prof->stop();
240 prof->start("construct_solution (residual)");
241 detail::construct_solution(residual, roots, m_subspace_solver->solutions(), {}, m_xspace->actionsq(),
242 m_xspace->actionsd(), m_xspace->dimensions().oP, m_xspace->dimensions().oQ,
243 m_xspace->dimensions().oD, *m_handlers);
244 prof->stop();
245 prof->start("apply");
246 auto pvectors = detail::construct_vectorP(roots, m_subspace_solver->solutions(), m_xspace->dimensions().oP,
247 m_xspace->dimensions().nP);
249 detail::normalise(roots.size(), parameters, residual, m_handlers->rr(), *m_logger);
250 if (m_apply_p)
251 m_apply_p(pvectors, m_xspace->cparamsp(), residual);
252 construct_residual(roots, cwrap(parameters), residual);
254 prof->stop();
255 };
256
257 void solution(const std::vector<int>& roots, std::vector<R>& parameters, std::vector<R>& residual) override {
258 return solution(roots, wrap(parameters), wrap(residual));
259 }
260 void solution(R& parameters, R& residual) override {
261 return solution(std::vector<int>(1, 0), wrap_arg(parameters), wrap_arg(residual));
262 }
263
264 void solution_params(const std::vector<int>& roots, std::vector<R>& parameters) override {
265 return solution_params(roots, wrap(parameters));
266 }
267
268 void solution_params(const std::vector<int>& roots, const VecRef<R>& parameters) override {
269 check_consistent_number_of_roots_and_solutions(roots, parameters.size());
270 detail::construct_solution(parameters, roots, m_subspace_solver->solutions(), m_xspace->paramsp(),
271 m_xspace->paramsq(), m_xspace->paramsd(), m_xspace->dimensions().oP,
272 m_xspace->dimensions().oQ, m_xspace->dimensions().oD, *m_handlers);
273 };
274
275 void solution_params(R& parameters) override { return solution_params(std::vector<int>(1, 0), wrap_arg(parameters)); }
276
277 // TODO Implement this
278 std::vector<size_t> suggest_p(const CVecRef<R>& solution, const CVecRef<R>& residual, size_t max_number,
279 double threshold) override {
280 return {};
281 }
282
283 const std::vector<int>& working_set() const override { return m_working_set; }
284
285 size_t n_roots() const override { return m_nroots; }
286
287 void set_n_roots(size_t roots) override {
288 m_nroots = roots;
289 m_working_set.resize(roots);
290 std::iota(begin(m_working_set), end(m_working_set), (int)0);
291 }
292
293 void set_options(const Options& options) override {
294 if (options.n_roots)
295 set_n_roots(options.n_roots.value());
296 if (options.convergence_threshold)
297 set_convergence_threshold(options.convergence_threshold.value());
298 if (options.verbosity)
299 set_verbosity(options.verbosity.value());
300 if (options.max_iter)
301 set_max_iter(options.max_iter.value());
302 if (options.max_p)
303 set_max_p(options.max_p.value());
304 if (options.p_threshold)
305 set_p_threshold(options.p_threshold.value());
306 }
307
308 std::shared_ptr<Options> get_options() const override {
309 auto options = std::make_shared<Options>();
310 options->n_roots = n_roots();
311 options->convergence_threshold = convergence_threshold();
312 options->verbosity = get_verbosity();
313 options->max_iter = get_max_iter();
314 options->max_p = get_max_p();
315 options->p_threshold = get_p_threshold();
316 return options;
317 }
318
319 const std::vector<scalar_type>& errors() const override { return m_errors; }
320
321 const Statistics& statistics() const override { return *m_stats; }
322
323 void report(std::ostream& cout, bool endl = true) const override {
324 cout << "iteration " << m_stats->iterations;
325 if (not m_errors.empty()) {
326 auto it_max_error = std::max_element(m_errors.cbegin(), m_errors.cend());
327 if (n_roots() > 1)
328 cout << ", |residual[" << std::distance(m_errors.cbegin(), it_max_error) << "]| = ";
329 else
330 cout << ", |residual| = ";
331 cout << std::scientific << *it_max_error << std::defaultfloat;
332 }
333 if (endl)
334 cout << std::endl;
335 }
336
337 void iteration_report() const {
338 std::stringstream sstream;
339 report(sstream, false);
340 this->m_logger->info<log::IterationReport>(sstream.str());
341 }
342
343 void summary_report() const {
344 std::stringstream sstream;
345 report(sstream, false);
346 this->m_logger->info<log::SummaryReport>(sstream.str());
347 }
348
349 void report() const override {
350 // We don't know what kind of report this is supposed to be -> default to summary
352 }
353
354 void set_convergence_threshold(double thresh) override { m_convergence_threshold = thresh; }
355 double convergence_threshold() const override { return m_convergence_threshold; }
356 void set_convergence_threshold_value(double thresh) override { m_convergence_threshold_value = thresh; }
358 void set_verbosity(Verbosity v) override { m_verbosity = v; }
359 void set_verbosity(int v) override {
360 if (v == 0) {
361 m_logger->set_verbosity(log::Verbosity::None);
362 m_logger->set_min_severity(log::Severity::Error);
363 } else if (v == 1) {
364 m_logger->set_verbosity(log::Verbosity::None);
365 m_logger->set_min_severity(log::Severity::Warning);
366 } else if (v == 2) {
367 m_logger->set_verbosity(log::Verbosity::Info);
368 m_logger->set_min_severity(log::Severity::Normal);
369 } else if (v == 3) {
370 m_logger->set_verbosity(log::Verbosity::Debug);
371 m_logger->set_min_severity(log::Severity::Normal);
372 } else {
373 m_logger->set_verbosity(log::Verbosity::Trace);
374 m_logger->set_min_severity(log::Severity::Normal);
375 }
376 }
377 Verbosity get_verbosity() const override { return m_verbosity.value_or(Verbosity::None); }
378 void set_max_iter(int n) override { m_max_iter = n; }
379 int get_max_iter() const override { return m_max_iter; }
380 void set_max_p(int n) override { m_max_p = n; }
381 int get_max_p() const override { return m_max_p; }
382 void set_p_threshold(double threshold) override { m_p_threshold = threshold; }
383 double get_p_threshold() const override { return m_p_threshold; }
385 const subspace::Dimensions& dimensions() const override { return m_xspace->dimensions(); }
386 scalar_type value() const override {
387 return m_xspace->data.count(subspace::EqnData::value) > 0 ? m_xspace->data[subspace::EqnData::value](0, 0)
388 : nan("molpro::linalg::itsolv::IterativeSolver::value");
389 }
390 // void set_profiler(molpro::profiler::Profiler& profiler) override { m_profiler.reset(&profiler); }
392 m_profiler = std::shared_ptr<molpro::profiler::Profiler>(&profiler);
393 }
394 const std::shared_ptr<molpro::profiler::Profiler>& profiler() const override { return m_profiler; }
395
396 bool solve(const VecRef<R>& parameters, const VecRef<R>& actions, const Problem<R>& problem,
397 bool generate_initial_guess = false) override {
398 if (parameters.empty())
399 throw std::runtime_error("Empty container passed to IterativeSolver::solve()");
400 if (parameters.size() != actions.size())
401 throw std::runtime_error("Inconsistent container sizes in IterativeSolver::solve()");
402 if (this->m_verbosity == Verbosity::None) {
403 // Backwards compatibility
404 this->m_logger->set_verbosity(log::Verbosity::None);
405 }
406 if (this->m_verbosity == Verbosity::Detailed) {
407 this->m_logger->set_verbosity(log::Verbosity::Trace);
408 this->m_logger->enable_data_dumps(true);
409 }
410 bool use_diagonals = problem.diagonals(actions.at(0));
411 std::unique_ptr<Q> diagonals;
412 if (use_diagonals)
413 diagonals.reset(new Q{m_handlers->qr().copy(actions.at(0))});
414// std::cout << "solve() generate_initial_guess "<<generate_initial_guess<<", roots "<<this->n_roots()<<std::endl;
415 if (generate_initial_guess) {
416 if (not use_diagonals)
417 throw std::runtime_error("Default initial guess requested, but diagonal elements are not available");
418 auto guess = m_handlers->qq().select(parameters.size(), *diagonals);
419 size_t root = 0;
420 if (!this->m_verbosity.has_value() || this->m_verbosity >= Verbosity::Summary) {
421 this->m_logger->info("Initial guess generated from diagonal elements");
422 }
423 for (const auto& g : guess) {
424 m_handlers->rp().copy(parameters[root], P{{g.first, 1}});
425 if (!this->m_verbosity.has_value() || this->m_verbosity >= Verbosity::Detailed) {
426 m_logger->debug("-> initial guess with index ", g.first);
427 }
428 root++;
429 }
430 }
431 int nwork = parameters.size();
432 std::vector<P> pspace;
433 if (use_diagonals and m_max_p > 0) {
434 auto selectp = m_handlers->qq().select(m_max_p, *diagonals);
435 if (!selectp.empty()) {
436 // selectp is keyed by index, not value; find the smallest selected
437 // diagonal explicitly before applying the threshold.
438 auto min_val = std::min_element(selectp.begin(), selectp.end(),
439 [](const auto& a, const auto& b) { return a.second < b.second; })
440 ->second;
441 for (auto s = selectp.begin(); s != selectp.end();) {
442 if (s->second > min_val + m_p_threshold)
443 s = selectp.erase(s);
444 else
445 ++s;
446 }
447 }
448 if (!selectp.empty() && (!this->m_verbosity.has_value() || this->m_verbosity >= Verbosity::Summary)) {
449 this->m_logger->info("P-space dimension, threshold and limit", selectp.size(), m_p_threshold, m_max_p);
450 }
451 if (!this->m_verbosity.has_value() || this->m_verbosity >= Verbosity::Detailed) {
452 for (const auto& s : selectp) {
453 this->m_logger->debug(std::format("P space element {}: {}", s.first, s.second));
454 }
455 }
456 for (const auto& s : selectp)
457 pspace.emplace_back((P){{s.first, 1}});
458 fapply_on_p_type apply_on_p = [&problem](const std::vector<std::vector<value_type>>& pcoeff,
459 const CVecRef<P>& pparams, const VecRef<R>& actions) {
460 problem.p_action(pcoeff, pparams, actions);
461 };
462 auto action_matrix = problem.pp_action_matrix(pspace);
463 nwork = add_p(cwrap(pspace), array::Span<value_type>(action_matrix.data(), action_matrix.size()), parameters,
464 actions, apply_on_p);
465 }
466 for (auto iter = 0; iter < this->m_max_iter && nwork > 0; iter++) {
467 m_logger->info<log::NewIteration>("Iteration, Current errors", iter, m_errors);
468 value_type value;
469 if (this->nonlinear()) {
470 value = problem.residual(*parameters.begin(), *actions.begin());
471 nwork = this->add_vector(*parameters.begin(), *actions.begin(), value);
472 } else if (iter > 0 or pspace.empty()) {
473 problem.action(cwrap(parameters.begin(), parameters.begin() + nwork),
474 wrap(actions.begin(), actions.begin() + nwork));
475 nwork = this->add_vector(parameters, actions);
476 }
477 // std::cout << "** nwork="<<nwork<<"use_diagonals="<<use_diagonals<<std::endl;
478 while (this->end_iteration_needed()) {
479 if (nwork > 0) {
480 if (use_diagonals) {
481 m_handlers->rq().copy(parameters.at(0), *diagonals);
482 problem.precondition(wrap(actions.begin(), actions.begin() + nwork), this->working_set_eigenvalues(),
483 parameters.at(0));
484 } else
485 problem.precondition(wrap(actions.begin(), actions.begin() + nwork), this->working_set_eigenvalues());
486 }
487 nwork = this->end_iteration(parameters, actions);
488 }
489 if (!this->m_verbosity.has_value() || this->m_verbosity >= Verbosity::Iteration) {
491 }
492 }
493 if (!this->m_verbosity.has_value() || this->m_verbosity >= Verbosity::Summary) {
495 }
496 if (!this->m_verbosity.has_value() || this->m_verbosity >= Verbosity::Summary) {
497 if (*std::max_element(m_errors.begin(), m_errors.end()) > m_convergence_threshold) {
498 this->m_logger->warn("Solver has not converged to threshold ", m_convergence_threshold);
499 } else {
500 this->m_logger->info("Solver converged");
501 }
502 }
503 return nwork == 0 and *std::max_element(m_errors.begin(), m_errors.end()) <= m_convergence_threshold;
504 }
505
506 bool solve(R& parameters, R& actions, const Problem<R>& problem, bool generate_initial_guess = false) override {
507 auto wparams = std::vector<std::reference_wrapper<R>>{std::ref(parameters)};
508 auto wactions = std::vector<std::reference_wrapper<R>>{std::ref(actions)};
509 return solve(wparams, wactions, problem, generate_initial_guess);
510 }
511 bool solve(std::vector<R>& parameters, std::vector<R>& actions, const Problem<R>& problem,
512 bool generate_initial_guess = false) override {
513 return solve(wrap(parameters), wrap(actions), problem, generate_initial_guess);
514 }
515
516 bool test_problem(const Problem<R>& problem, R& v0, R& v1, int verbosity, double threshold) const override {
517 bool success = true;
518 double step=1e-4;
519 if (this->nonlinear()) {
520 if (!problem.test_parameters(0, v0))
521 return true;
522 auto value0 = problem.residual(v0, v1);
523 if (verbosity > 1) {
524 std::cout << "value0 " << value0 << std::endl;
525 // std::cout << "parameters0 " << v0[0] << "," << v0[1] << ",..." << std::endl;
526 // std::cout << "residual0 " << v1[0] << "," << v1[1] << ",..." << std::endl;
527 }
528 Q parameters0 = m_handlers->qr().copy(v0);
529 Q residual0 = m_handlers->qr().copy(v1);
530 for (int instance = 1; problem.test_parameters(instance, v0); ++instance) {
531 m_handlers->rq().axpy(-1.0, parameters0, v0);
532 m_handlers->rr().scal(1 / std::sqrt(m_handlers->rr().dot(v0, v0)), v0);
533 Q step1 = m_handlers->qr().copy(v0);
534 auto residual_analytic = m_handlers->rq().dot(v0, residual0);
535 m_handlers->rq().copy(v0, parameters0);
536 m_handlers->rq().axpy(-2*step, step1, v0);
537 auto valuem2 = problem.residual(v0, v1);
538 m_handlers->rq().copy(v0, parameters0);
539 m_handlers->rq().axpy(-1*step, step1, v0);
540 auto valuem1 = problem.residual(v0, v1);
541 m_handlers->rq().copy(v0, parameters0);
542 m_handlers->rq().axpy(+1*step, step1, v0);
543 auto valuep1 = problem.residual(v0, v1);
544 m_handlers->rq().copy(v0, parameters0);
545 m_handlers->rq().axpy(+2*step, step1, v0);
546 auto valuep2 = problem.residual(v0, v1);
547 auto residual_numerical = (valuem2 - 8*valuem1 + 8*valuep1 - valuep2) / (12*step);
548 if (verbosity > 1)
549 std::cout << "testing problem class, instance: " << instance << ", numerical: " << residual_numerical <<", analytical: "<<residual_analytic <<", difference: "<<residual_numerical-residual_analytic << std::endl;
550 success = success && std::abs(residual_numerical - residual_analytic) < threshold;
551 }
552 } else {
553 for (int instance = 0; problem.test_parameters(instance, v0); ++instance) {
554 problem.action({v0}, {v1});
555 Q residual = m_handlers->qr().copy(v1);
556 auto norm2_residual = std::sqrt(m_handlers->rr().dot(v1, v1));
557 constexpr double scale_factor{10.0};
558 m_handlers->rr().scal(scale_factor, v0);
559 problem.action({v0}, {v1});
560 m_handlers->rq().axpy(-scale_factor, residual, v1);
561 auto norm2 = std::sqrt(m_handlers->rr().dot(v1, v1));
562 success = success && std::abs(norm2 / norm2_residual) < threshold;
563 if (verbosity > 0 or (verbosity > -1 and not success))
564 std::cout << "Length of residual: " << norm2_residual << ", scaling defect: " << norm2 << std::endl;
565 }
566 }
567 return success;
568 }
569
570protected:
572 std::shared_ptr<subspace::ISubspaceSolver<R, Q, P>> solver,
573 std::shared_ptr<ArrayHandlers<R, Q, P>> handlers, std::shared_ptr<Statistics> stats,
574 std::shared_ptr<Logger> logger)
575 : m_handlers(std::move(handlers)), m_xspace(std::move(xspace)), m_subspace_solver(std::move(solver)),
576 m_stats(std::move(stats)), m_logger(std::move(logger)), m_profiler(molpro::Profiler::single()),
577 m_profiler_saved_depth(m_profiler->get_max_depth()) {
578 set_n_roots(1);
579 m_profiler->set_max_depth(options()->parameter("PROFILER_DEPTH", 0));
580 }
581
583 if (molpro::mpi::rank_global() == 0) {
584 auto file = options()->parameter("PROFILER_OUTPUT", "");
585 if (profiler()->get_max_depth() > 0 and
586 std::find_if(file.begin(), file.end(), [](unsigned char ch) { return !std::isspace(ch); }) != file.end())
587 std::ofstream(file) << *profiler() << std::endl;
588 }
589 if (molpro::mpi::rank_global() == 0) {
590 auto file = options()->parameter("PROFILER_DOTGRAPH", "");
591 if (profiler()->get_max_depth() > 0 and
592 std::find_if(file.begin(), file.end(), [](unsigned char ch) { return !std::isspace(ch); }) != file.end())
593 profiler()->dotgraph(file, options()->parameter("PROFILER_THRESHOLD", .01));
594 }
595 molpro::Profiler::single()->set_max_depth(m_profiler_saved_depth);
596 }
597
599 virtual void set_value_errors() {}
601 virtual void construct_residual(const std::vector<int>& roots, const CVecRef<R>& params,
602 const VecRef<R>& actions) = 0;
603 virtual bool linearEigensystem() const { return false; }
612 size_t solve_and_generate_working_set(const VecRef<R>& parameters, const VecRef<R>& action) {
613 auto prof = profiler();
614 auto p = prof->push("itsolv::solve_and_generate_working_set");
615 prof->start("itsolv::solve");
617 prof->stop();
618 prof->start("itsolv::temp_solutions");
619 auto nsol = m_subspace_solver->size();
620 std::vector<std::pair<Q, Q>> temp_solutions{};
621 const auto batches = detail::parameter_batches(nsol, parameters.size());
622 for (const auto& batch : batches) {
623 auto [start_sol, end_sol] = batch;
624 auto roots = std::vector<int>(end_sol - start_sol);
625 std::iota(begin(roots), end(roots), start_sol);
626 solution(roots, parameters, action);
627 auto errors = std::vector<scalar_type>(roots.size(), 0);
629 if (batches.size() > 1) {
630 for (size_t i = 0; i < roots.size(); ++i)
631 temp_solutions.emplace_back(m_handlers->qr().copy(parameters[i]), m_handlers->qr().copy(action[i]));
632 m_stats->q_creations += 2 * roots.size();
633 }
634 m_subspace_solver->set_error(roots, errors);
635 }
636 prof->stop();
638 m_errors = m_subspace_solver->errors();
641 for (size_t i = 0; i < m_working_set.size(); ++i) {
642 size_t root = m_working_set[i];
643 if (batches.size() > 1) {
644 m_handlers->rq().copy(parameters[i], temp_solutions.at(root).first);
645 m_handlers->rq().copy(action[i], temp_solutions.at(root).second);
646 } else {
647 if (root < i)
648 throw std::logic_error("incorrect ordering of roots");
649 if (root > i && root < parameters.size()) {
650 m_handlers->rr().copy(parameters[i], parameters[root]);
651 m_handlers->rr().copy(action[i], action[root]);
652 }
653 }
654 }
655 m_logger->trace("add_vector::errors = ", m_errors);
656 return m_working_set.size();
657 }
658
659 template <typename TTT>
660 void check_consistent_number_of_roots_and_solutions(const std::vector<TTT>& roots, const size_t nparams) {
661 if (roots.size() > nparams)
662 throw std::runtime_error("asking for more roots than parameters");
663 if (!roots.empty() &&
664 size_t(*std::max_element(roots.begin(), roots.end())) >= m_subspace_solver->solutions().size())
665 throw std::runtime_error("asking for more roots than there are solutions");
666 }
667
669
670
671 std::shared_ptr<ArrayHandlers<R, Q, P>> m_handlers;
672 std::shared_ptr<subspace::IXSpace<R, Q, P>> m_xspace;
673 std::shared_ptr<subspace::ISubspaceSolver<R, Q, P>> m_subspace_solver;
674 std::vector<double> m_errors;
675 std::vector<double> m_value_errors;
676 std::vector<int> m_working_set;
677 size_t m_nroots{0};
678 double m_convergence_threshold{1.0e-8};
680 std::numeric_limits<double>::max()};
681 std::shared_ptr<Statistics> m_stats;
682 std::shared_ptr<Logger> m_logger;
683 bool m_normalise_solution = false;
684 fapply_on_p_type m_apply_p = {};
685 std::optional<Verbosity> m_verbosity = {};
686 int m_max_iter = 100;
687 size_t m_max_p = 0;
688 double m_p_threshold = std::numeric_limits<double>::max();
689private:
690 mutable std::shared_ptr<molpro::profiler::Profiler> m_profiler;
691 int m_profiler_saved_depth;
692protected:
694};
695
696} // namespace molpro::linalg::itsolv
697
698#endif // LINEARALGEBRA_SRC_MOLPRO_LINALG_ITSOLV_ITERATIVESOLVERTEMPLATE_H
Enhances various operations between pairs of arrays and allows dynamic code injection with uniform in...
Definition: ArrayHandler.h:162
virtual value_type dot(const AL &x, const AR &y)=0
virtual void scal(value_type alpha, AL &x)=0
Non-owning container taking a pointer to the data buffer and its size and exposing routines for itera...
Definition: Span.h:31
Class, containing a collection of array handlers used in IterativeSolver Provides a Builder sub-class...
Definition: ArrayHandlers.h:25
auto & rp()
Definition: ArrayHandlers.h:44
auto & rr()
Definition: ArrayHandlers.h:39
auto & rq()
Definition: ArrayHandlers.h:42
Implements IterativeSolver interface that is common to all solvers.
Definition: IterativeSolverTemplate.h:161
void iteration_report() const
Definition: IterativeSolverTemplate.h:337
void set_profiler(molpro::profiler::Profiler &profiler) override
Definition: IterativeSolverTemplate.h:391
std::optional< Verbosity > m_verbosity
how much output to print in solve()
Definition: IterativeSolverTemplate.h:685
size_t add_p(const CVecRef< P > &pparams, const array::Span< value_type > &pp_action_matrix, const VecRef< R > &parameters, const VecRef< R > &actions, fapply_on_p_type apply_p) override
Definition: IterativeSolverTemplate.h:216
void set_convergence_threshold(double thresh) override
Definition: IterativeSolverTemplate.h:354
int get_max_iter() const override
Definition: IterativeSolverTemplate.h:379
void solution_params(const std::vector< int > &roots, std::vector< R > &parameters) override
Definition: IterativeSolverTemplate.h:264
void report(std::ostream &cout, bool endl=true) const override
Definition: IterativeSolverTemplate.h:323
double get_p_threshold() const override
Definition: IterativeSolverTemplate.h:383
std::vector< int > m_working_set
indices of roots in the working set
Definition: IterativeSolverTemplate.h:676
void summary_report() const
Definition: IterativeSolverTemplate.h:343
scalar_type value() const override
Definition: IterativeSolverTemplate.h:386
std::shared_ptr< subspace::IXSpace< R, Q, P > > m_xspace
manages the subspace and associated data
Definition: IterativeSolverTemplate.h:672
void set_n_roots(size_t roots) override
Definition: IterativeSolverTemplate.h:287
std::shared_ptr< ArrayHandlers< R, Q, P > > m_handlers
Array handlers.
Definition: IterativeSolverTemplate.h:671
fapply_on_p_type m_apply_p
function that evaluates effect of action on the P space projection
Definition: IterativeSolverTemplate.h:684
void solution(const std::vector< int > &roots, std::vector< R > &parameters, std::vector< R > &residual) override
Definition: IterativeSolverTemplate.h:257
void set_max_iter(int n) override
Definition: IterativeSolverTemplate.h:378
void check_consistent_number_of_roots_and_solutions(const std::vector< TTT > &roots, const size_t nparams)
Definition: IterativeSolverTemplate.h:660
bool test_problem(const Problem< R > &problem, R &v0, R &v1, int verbosity, double threshold) const override
Definition: IterativeSolverTemplate.h:516
void set_max_p(int n) override
Definition: IterativeSolverTemplate.h:380
std::shared_ptr< subspace::ISubspaceSolver< R, Q, P > > m_subspace_solver
solves the subspace problem
Definition: IterativeSolverTemplate.h:673
int get_max_p() const override
Definition: IterativeSolverTemplate.h:381
size_t n_roots() const override
Definition: IterativeSolverTemplate.h:285
double convergence_threshold() const override
Definition: IterativeSolverTemplate.h:355
std::shared_ptr< Options > get_options() const override
Definition: IterativeSolverTemplate.h:308
size_t m_max_p
maximum size of P space
Definition: IterativeSolverTemplate.h:687
double convergence_threshold_value() const override
Definition: IterativeSolverTemplate.h:357
const subspace::Dimensions & dimensions() const override
Access dimensions of the subspace.
Definition: IterativeSolverTemplate.h:385
size_t solve_and_generate_working_set(const VecRef< R > &parameters, const VecRef< R > &action)
Solves the subspace problems and selects the working set of roots, returning their parameters and res...
Definition: IterativeSolverTemplate.h:612
const std::shared_ptr< molpro::profiler::Profiler > & profiler() const override
Definition: IterativeSolverTemplate.h:394
IterativeSolverTemplate(std::shared_ptr< subspace::IXSpace< R, Q, P > > xspace, std::shared_ptr< subspace::ISubspaceSolver< R, Q, P > > solver, std::shared_ptr< ArrayHandlers< R, Q, P > > handlers, std::shared_ptr< Statistics > stats, std::shared_ptr< Logger > logger)
Definition: IterativeSolverTemplate.h:571
virtual ~IterativeSolverTemplate()
Definition: IterativeSolverTemplate.h:582
void set_p_threshold(double threshold) override
Definition: IterativeSolverTemplate.h:382
int add_vector(std::vector< R > &parameters, std::vector< R > &actions) override
Definition: IterativeSolverTemplate.h:208
void solution_params(R &parameters) override
Definition: IterativeSolverTemplate.h:275
virtual bool linearEigensystem() const
Definition: IterativeSolverTemplate.h:603
const Statistics & statistics() const override
Definition: IterativeSolverTemplate.h:321
void solution(R &parameters, R &residual) override
Definition: IterativeSolverTemplate.h:260
size_t m_nroots
number of roots the solver is searching for
Definition: IterativeSolverTemplate.h:677
double m_p_threshold
threshold for selecting P space
Definition: IterativeSolverTemplate.h:688
bool m_end_iteration_needed
whether end_iteration should be called after any preconditioner
Definition: IterativeSolverTemplate.h:693
void set_convergence_threshold_value(double thresh) override
Definition: IterativeSolverTemplate.h:356
void solution_params(const std::vector< int > &roots, const VecRef< R > &parameters) override
Definition: IterativeSolverTemplate.h:268
int add_vector(const VecRef< R > &parameters, const VecRef< R > &actions) override
Definition: IterativeSolverTemplate.h:183
void clearP() override
Definition: IterativeSolverTemplate.h:229
bool solve(R &parameters, R &actions, const Problem< R > &problem, bool generate_initial_guess=false) override
Definition: IterativeSolverTemplate.h:506
bool end_iteration_needed() override
Definition: IterativeSolverTemplate.h:668
IterativeSolverTemplate(const IterativeSolverTemplate< Solver, R, Q, P > &)=delete
const std::vector< scalar_type > & errors() const override
Definition: IterativeSolverTemplate.h:319
void report() const override
Definition: IterativeSolverTemplate.h:349
bool m_normalise_solution
whether to normalise the solutions
Definition: IterativeSolverTemplate.h:683
double m_convergence_threshold
residual norms less than this mark a converged solution
Definition: IterativeSolverTemplate.h:678
std::vector< size_t > suggest_p(const CVecRef< R > &solution, const CVecRef< R > &residual, size_t max_number, double threshold) override
Definition: IterativeSolverTemplate.h:278
const std::vector< int > & working_set() const override
Definition: IterativeSolverTemplate.h:283
std::shared_ptr< Logger > m_logger
logger
Definition: IterativeSolverTemplate.h:682
virtual void set_value_errors()
Implementation class should overload this to set errors in the current values (e.g....
Definition: IterativeSolverTemplate.h:599
bool solve(std::vector< R > &parameters, std::vector< R > &actions, const Problem< R > &problem, bool generate_initial_guess=false) override
Definition: IterativeSolverTemplate.h:511
void set_logger(std::shared_ptr< Logger > logger) override
Definition: IterativeSolverTemplate.h:174
virtual void construct_residual(const std::vector< int > &roots, const CVecRef< R > &params, const VecRef< R > &actions)=0
Constructs residual for given roots provided their parameters and actions.
void set_options(const Options &options) override
Definition: IterativeSolverTemplate.h:293
std::vector< double > m_errors
errors from the most recent solution
Definition: IterativeSolverTemplate.h:674
int add_vector(R &parameters, R &actions, value_type value=0) override
Definition: IterativeSolverTemplate.h:211
void set_verbosity(Verbosity v) override
Definition: IterativeSolverTemplate.h:358
bool solve(const VecRef< R > &parameters, const VecRef< R > &actions, const Problem< R > &problem, bool generate_initial_guess=false) override
Definition: IterativeSolverTemplate.h:396
int m_max_iter
maximum number of iterations in solve()
Definition: IterativeSolverTemplate.h:686
std::shared_ptr< Statistics > m_stats
accumulates statistics of operations performed by the solver
Definition: IterativeSolverTemplate.h:681
double m_convergence_threshold_value
value changes less than this mark a converged solution
Definition: IterativeSolverTemplate.h:679
Logger & logger() override
Definition: IterativeSolverTemplate.h:181
void set_verbosity(int v) override
Definition: IterativeSolverTemplate.h:359
Verbosity get_verbosity() const override
Definition: IterativeSolverTemplate.h:377
void solution(const std::vector< int > &roots, const VecRef< R > &parameters, const VecRef< R > &residual) override
Definition: IterativeSolverTemplate.h:231
std::vector< double > m_value_errors
value errors from the most recent solution
Definition: IterativeSolverTemplate.h:675
IterativeSolverTemplate(IterativeSolverTemplate< Solver, R, Q, P > &&) noexcept=default
Definition: Logger.h:428
void warn(std::string_view message, Ts &&...args) const
Definition: Logger.h:496
Abstract class defining the problem-specific interface for the simplified solver interface to Iterati...
Definition: IterativeSolver.h:78
virtual void action(const CVecRef< R > &parameters, const VecRef< R > &action) const
Calculate the action of the kernel matrix on a set of parameters. Used by linear solvers,...
Definition: IterativeSolver.h:102
virtual std::vector< double > pp_action_matrix(const std::vector< P > &pparams) const
Calculate the kernel matrix in the P space.
Definition: IterativeSolver.h:155
virtual void precondition(const VecRef< R > &residual, const std::vector< value_t > &shift) const
Apply preconditioning to a residual vector in order to predict a step towards the solution.
Definition: IterativeSolver.h:125
virtual bool diagonals(container_t &d) const
Optionally provide the diagonal elements of the underlying kernel. If implemented and returning true,...
Definition: IterativeSolver.h:114
virtual void p_action(const std::vector< std::vector< value_t > > &p_coefficients, const CVecRef< P > &pparams, const VecRef< container_t > &actions) const
Calculate the action of the kernel matrix on a set of vectors in the P space.
Definition: IterativeSolver.h:167
virtual value_t residual(const R &parameters, R &residual) const
Calculate the residual vector. Used by non-linear solvers (NonLinearEquations, Optimize) only.
Definition: IterativeSolver.h:94
virtual bool test_parameters(unsigned int instance, R &parameters) const
Provide values of R vectors for testing the problem class. For use in a non-linear solver,...
Definition: IterativeSolver.h:182
Matrix container that allows simple data access, slicing, copying and resizing without loosing data.
Definition: Matrix.h:30
static std::shared_ptr< Profiler > single()
Definition: Logger.h:123
std::vector< std::vector< T > > construct_vectorP(const std::vector< int > &roots, const subspace::Matrix< T > &solutions, const size_t oP, const size_t nP)
Definition: IterativeSolverTemplate.h:72
void construct_solution(const VecRef< R > &params, const std::vector< int > &roots, const subspace::Matrix< double > &solutions, const std::vector< std::reference_wrapper< P > > &pparams, const std::vector< std::reference_wrapper< Q > > &qparams, const std::vector< std::reference_wrapper< Q > > &dparams, size_t oP, size_t oQ, size_t oD, ArrayHandlers< R, Q, P > &handlers)
Definition: IterativeSolverTemplate.h:38
void normalise(const size_t n_roots, const VecRef< R > &params, const VecRef< R > &actions, array::ArrayHandler< R, R > &handler, Logger &logger)
Definition: IterativeSolverTemplate.h:84
std::vector< std::pair< size_t, size_t > > parameter_batches(const size_t nsol, const size_t nparam)
Definition: IterativeSolverTemplate.h:25
std::vector< int > select_working_set(const size_t nw, const std::vector< T > &errors, const T threshold, const std::vector< T > &value_errors, const T value_threshold)
Definition: IterativeSolverTemplate.h:109
void update_errors(std::vector< T > &errors, const CVecRef< R > &residual, array::ArrayHandler< R, R > &handler)
Definition: IterativeSolverTemplate.h:100
4-parameter interpolation of a 1-dimensional function given two points for which function values and ...
Definition: helper.h:11
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 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.
Definition: wrap.h:135
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
Verbosity
Specifies how much detail IterativeSolver::solve should write to output.
Definition: Options.h:12
@ Summary
only summary at the start and end of solve
@ Detailed
Iteration plus more detail such as operation count.
@ Iteration
Summary plus minimal results at each iteration.
void read_handler_counts(std::shared_ptr< Statistics > stats, std::shared_ptr< ArrayHandlers< R, Q, P > > handlers)
Definition: Statistics.h:30
const std::shared_ptr< const molpro::Options > options()
Get the Options object associated with iterative-solver.
Definition: options.cpp:4
Access point for different options in iterative solvers.
Definition: Options.h:20
Information about performance of IterativeSolver instance.
Definition: Statistics.h:10
Definition: IterativeSolverTemplate.h:140
static constexpr const char * name
Definition: IterativeSolverTemplate.h:141
Definition: IterativeSolverTemplate.h:130
static constexpr std::size_t iter
Definition: IterativeSolverTemplate.h:132
static constexpr std::size_t errors
Definition: IterativeSolverTemplate.h:133
static constexpr const char * name
Definition: IterativeSolverTemplate.h:131
Definition: IterativeSolverTemplate.h:148
static constexpr const char * name
Definition: IterativeSolverTemplate.h:149
Stores partitioning of XSpace into P, Q and R blocks with sizes and offsets for each one.
Definition: Dimensions.h:8
Manages solution of the subspace problem and storage of those solutions.
Definition: ISubspaceSolver.h:21