numerical-collection-cpp 0.10.0
A collection of algorithms in numerical analysis implemented in C++
Loading...
Searching...
No Matches
random_multi_quadratic_function.h
Go to the documentation of this file.
1/*
2 * Copyright 2023 MusicScience37 (Kenta Kabashima)
3 *
4 * Licensed under the Apache License, Version 2.0 (the "License");
5 * you may not use this file except in compliance with the License.
6 * You may obtain a copy of the License at
7 *
8 * http://www.apache.org/licenses/LICENSE-2.0
9 *
10 * Unless required by applicable law or agreed to in writing, software
11 * distributed under the License is distributed on an "AS IS" BASIS,
12 * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
13 * See the License for the specific language governing permissions and
14 * limitations under the License.
15 */
20#pragma once
21
22#include <random>
23
24#include <Eigen/Core>
25
27
28namespace num_prob_collect::opt {
29
35public:
37 using variable_type = Eigen::VectorXd;
39 using value_type = double;
41 using hessian_type = Eigen::MatrixXd;
42
50 Eigen::VectorXd optimal_variable, Eigen::VectorXd coefficients)
53 constexpr double two = 2.0;
54 hessian_ = two * coefficients_.asDiagonal();
55 }
56
62 void evaluate_on(const Eigen::VectorXd& variable) noexcept {
63 diff_ = variable - optimal_variable_;
64 value_ = diff_.array().square().matrix().dot(coefficients_);
65 constexpr double two = 2.0;
66 grad_ = two * diff_.cwiseProduct(coefficients_);
67 }
68
74 [[nodiscard]] auto value() const noexcept -> const double& {
75 return value_;
76 }
77
83 [[nodiscard]] auto gradient() const noexcept -> const Eigen::VectorXd& {
84 return grad_;
85 }
86
92 [[nodiscard]] auto hessian() const noexcept -> const Eigen::MatrixXd& {
93 return hessian_;
94 }
95
101 [[nodiscard]] auto optimal_variable() const noexcept
102 -> const Eigen::VectorXd& {
103 return optimal_variable_;
104 }
105
111 [[nodiscard]] auto coefficients() const noexcept -> const Eigen::VectorXd& {
112 return coefficients_;
113 }
114
115private:
117 Eigen::VectorXd optimal_variable_;
118
120 Eigen::VectorXd coefficients_;
121
123 Eigen::VectorXd diff_{};
124
126 double value_{0.0};
127
129 Eigen::VectorXd grad_{};
130
132 Eigen::MatrixXd hessian_{};
133};
134
139public:
141 static constexpr double min_variable = -10.0;
142
144 static constexpr double max_variable = 10.0;
145
152 num_collect::index_type num_variables)
153 : num_variables_(num_variables) {}
154
160 [[nodiscard]] auto operator()() noexcept
162 Eigen::VectorXd optimal_variable;
163 optimal_variable.resize(num_variables_);
164 std::generate(optimal_variable.begin(), optimal_variable.end(),
165 [this] { return optimal_variable_dist_(generator_); });
166
167 Eigen::VectorXd coefficients;
168 coefficients.resize(num_variables_);
169 std::generate(coefficients.begin(), coefficients.end(),
170 [this] { return coeff_dist_(generator_); });
171
173 std::move(optimal_variable), std::move(coefficients)};
174 }
175
176private:
178 std::mt19937
179 generator_{}; // NOLINT(cert-msc32-c,cert-msc51-cpp): For reproducibility.
180
182 std::uniform_real_distribution<double> coeff_dist_{0.5, 20.0}; // NOLINT
183
185 std::uniform_real_distribution<double> optimal_variable_dist_{
187
190};
191
192} // namespace num_prob_collect::opt
Class to generate random random_multi_quadratic_function objects.
random_multi_quadratic_function_generator(num_collect::index_type num_variables)
Constructor.
auto operator()() noexcept -> random_multi_quadratic_function
Generate a random quadratic function.
std::uniform_real_distribution< double > optimal_variable_dist_
Distribution of optimal variables.
std::uniform_real_distribution< double > coeff_dist_
Distribution of coefficients.
Class of multi-variate quadratic functions with random coefficients and optimal variables.
random_multi_quadratic_function(Eigen::VectorXd optimal_variable, Eigen::VectorXd coefficients)
Constructor.
auto optimal_variable() const noexcept -> const Eigen::VectorXd &
Get the optimal variable.
auto coefficients() const noexcept -> const Eigen::VectorXd &
Get the coefficients.
auto gradient() const noexcept -> const Eigen::VectorXd &
Get gradient.
auto value() const noexcept -> const double &
Get function value.
void evaluate_on(const Eigen::VectorXd &variable) noexcept
Evaluate function value on variable.
auto hessian() const noexcept -> const Eigen::MatrixXd &
Get Hessian.
Definition of index_type type.
Namespace of Eigen library.
Definition variable.h:416
std::ptrdiff_t index_type
Type of indices in this library.
Definition index_type.h:33
Namespace of optimization problems.
Definition namespaces.h:25
STL namespace.