iterative-solver 0.0
OptimizeSD.h
1#ifndef LINEARALGEBRA_SRC_MOLPRO_LINALG_ITSOLV_OPTIMIZESD_H
2#define LINEARALGEBRA_SRC_MOLPRO_LINALG_ITSOLV_OPTIMIZESD_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#include <molpro/linalg/itsolv/subspace/SubspaceSolverOptSD.h>
8#include <molpro/linalg/itsolv/subspace/XSpace.h>
9
10namespace molpro::linalg::itsolv {
19template <class R, class Q = R, class P = std::map<size_t, typename R::value_type>>
20class OptimizeSD : public IterativeSolverTemplate<Optimize, R, Q, P> {
21public:
24 using typename SolverTemplate::scalar_type;
25 using typename SolverTemplate::value_type;
26 using typename SolverTemplate::value_type_abs;
28
29 explicit OptimizeSD(const std::shared_ptr<ArrayHandlers<R, Q, P>>& handlers,
30 const std::shared_ptr<Logger>& logger_ = std::make_shared<Logger>())
31 : SolverTemplate(std::make_shared<subspace::XSpace<R, Q, P>>(handlers, logger_),
32 std::static_pointer_cast<subspace::ISubspaceSolver<R, Q, P>>(
33 std::make_shared<subspace::SubspaceSolverOptSD<R, Q, P>>(logger_)),
34 handlers, std::make_shared<Statistics>(), logger_),
35 logger(logger_) {}
36
37 bool nonlinear() const override { return true; }
38
39 size_t end_iteration(const VecRef<R>& parameters, const VecRef<R>& action) override {
40 this->solution_params(this->m_working_set, parameters);
41 if (this->m_errors.front() < this->m_convergence_threshold) {
42 this->m_working_set.clear();
43 return 0;
44 }
45 this->m_working_set.assign(1, 0);
46 this->m_handlers->rr().axpy(-1, action.front(), parameters.front());
47 this->m_stats->iterations++;
48 return 1;
49 }
50
51 size_t end_iteration(std::vector<R>& parameters, std::vector<R>& action) override {
52 auto result = end_iteration(wrap(parameters), wrap(action));
53 return result;
54 }
55
56 void set_value_errors() override {
57 auto& Value = this->m_xspace->data[subspace::EqnData::value];
58 this->m_value_errors.assign(1, std::numeric_limits<double>::max());
59 if (this->m_xspace->size() > 1 and Value(0, 0) < Value(1, 0))
60 this->m_value_errors.front() = Value(1, 0) - Value(0, 0);
61 }
62
63 void set_options(const Options& options) override {
65 if (auto opt2 = dynamic_cast<const molpro::linalg::itsolv::OptimizeSDOptions*>(&options)) {
67 }
68 }
69
70 std::shared_ptr<Options> get_options() const override {
71 auto opt = std::make_shared<OptimizeSDOptions>();
72 opt->copy(*SolverTemplate::get_options());
73 return opt;
74 }
75
76 void report(std::ostream& cout, bool endl = true) const override {
77 SolverTemplate::report(cout, false);
78 cout << ", value " << this->value();
79 if (endl)
80 cout << std::endl;
81 }
82 std::shared_ptr<Logger> logger;
83
84 int add_vector(R& parameters, R& residual, value_type value) override {
85 using namespace subspace;
86 auto& xspace = this->m_xspace;
87 auto& xdata = xspace->data;
88 const auto n = this->m_xspace->dimensions().nX;
89 xdata[EqnData::value].resize({n + 1, 1});
90 xdata[EqnData::value](0, 0) = value;
91 auto nwork = IterativeSolverTemplate<Optimize, R, Q, P>::add_vector(parameters, residual);
92 return nwork;
93 }
94
95 size_t end_iteration(R& parameters, R& actions) override {
96 auto wparams = std::vector<std::reference_wrapper<R>>{std::ref(parameters)};
97 auto wactions = std::vector<std::reference_wrapper<R>>{std::ref(actions)};
98 return end_iteration(wparams, wactions);
99 }
100
101protected:
102 // for non-linear problems, actions already contains the residual
103 void construct_residual(const std::vector<int>& roots, const CVecRef<R>& params, const VecRef<R>& actions) override {}
104};
105
106} // namespace molpro::linalg::itsolv
107#endif // LINEARALGEBRA_SRC_MOLPRO_LINALG_ITSOLV_OPTIMIZESD_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
void solution_params(const std::vector< int > &roots, std::vector< R > &parameters) override
Definition: IterativeSolverTemplate.h:224
std::vector< int > m_working_set
indices of roots in the working set
Definition: IterativeSolverTemplate.h:583
scalar_type value() const override
Report the function value for the current optimum solution.
Definition: IterativeSolverTemplate.h:312
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< Options > get_options() const override
Definition: IterativeSolverTemplate.h:262
int add_vector(const VecRef< R > &parameters, const VecRef< R > &actions) override
Definition: IterativeSolverTemplate.h:140
void report() const override
Definition: IterativeSolverTemplate.h:287
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
std::shared_ptr< Statistics > m_stats
accumulates statistics of operations performed by the solver
Definition: IterativeSolverTemplate.h:588
std::vector< double > m_value_errors
value errors from the most recent solution
Definition: IterativeSolverTemplate.h:582
typename R::value_type value_type
The underlying type of elements of vectors.
Definition: IterativeSolver.h:202
A class that optimises a function using a simple steepest-descent method.
Definition: OptimizeSD.h:20
size_t end_iteration(R &parameters, R &actions) override
Definition: OptimizeSD.h:95
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: OptimizeSD.h:103
void report(std::ostream &cout, bool endl=true) const override
Writes a report to cout output stream.
Definition: OptimizeSD.h:76
bool nonlinear() const override
Report whether the class is a non-linear solver.
Definition: OptimizeSD.h:37
size_t end_iteration(const VecRef< R > &parameters, const VecRef< R > &action) override
Behaviour depends on the solver.
Definition: OptimizeSD.h:39
void set_value_errors() override
Implementation class should overload this to set errors in the current values (e.g....
Definition: OptimizeSD.h:56
std::shared_ptr< Logger > logger
Definition: OptimizeSD.h:82
size_t end_iteration(std::vector< R > &parameters, std::vector< R > &action) override
Definition: OptimizeSD.h:51
int add_vector(R &parameters, R &residual, value_type value) override
Definition: OptimizeSD.h:84
std::shared_ptr< Options > get_options() const override
Definition: OptimizeSD.h:70
OptimizeSD(const std::shared_ptr< ArrayHandlers< R, Q, P > > &handlers, const std::shared_ptr< Logger > &logger_=std::make_shared< Logger >())
Definition: OptimizeSD.h:29
void report() const override
Writes a report to std::cout.
Definition: IterativeSolverTemplate.h:287
void set_options(const Options &options) override
Definition: OptimizeSD.h:63
Solves subspace problem for minimisation using the Steepest Descent algorithm.
Definition: SubspaceSolverOptSD.h:13
4-parameter interpolation of a 1-dimensional function given two points for which function values and ...
Definition: helper.h:10
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< OptimizeSDOptions > OptimizeSD(const std::shared_ptr< Options > &options)
Definition: CastOptions.h:90
Allows setting and getting of options for OptimizeBFGS instance via IterativeSolver base class.
Definition: OptimizeSDOptions.h:9
Access point for different options in iterative solvers.
Definition: Options.h:20
Information about performance of IterativeSolver instance.
Definition: Statistics.h:10