On this page, we look at an example application of the sgpp::optimization module.
Versions of the example are given in all languages currently supported by SG++: C++, Python, Java, and MATLAB.
The example interpolates a bivariate test function with B-splines instead of piecewise linear basis functions to obtain a smoother interpolant. The resulting sparse grid function is then minimized with the method of steepest descent. For comparison, we also minimize the objective function with Nelder-Mead's method.
First, we include all the necessary headers, including those of the sgpp::base and sgpp::optimization module.
#include <algorithm>
#include <iostream>
#include <iterator>
The function \(f\colon [0, 1]^d \to \mathbb{R}\) to be minimized is called objective function and has to derive from sgpp::base::ScalarFunction. In the constructor, we give the dimensionality of the domain (in this case \(d = 2\)). The eval method evaluates the objective function and returns the function value \(f(\vec{x})\) for a given point \(\vec{x} \in [0, 1]^d\). The clone method returns a std::unique_ptr to a clone of the object and is used for parallelization (in case eval is not thread-safe).
public:
ExampleFunction() :
sgpp::base::ScalarFunction(2) {}
return std::sin(8.0 * x[0]) + std::sin(7.0 * x[1]);
}
virtual void clone(std::unique_ptr<sgpp::base::ScalarFunction>& clone)
const {
clone = std::unique_ptr<sgpp::base::ScalarFunction>(new ExampleFunction(*this));
}
};
A class to store one-dimensional data.
Definition DataVector.hpp:25
Abstract base class for scalar-valued functions (e.g., objective functions in optimization).
Definition ScalarFunction.hpp:23
virtual void clone(std::unique_ptr< ScalarFunction > &clone) const =0
Pure virtual method for cloning the function.
virtual double eval(const DataVector &x)=0
Pure virtual method for calculating .
Class that stores, generates and manipulates a density function during online phase in on/off learnin...
Definition AlgorithmDGEMV.hpp:22
Now, we can start with the main
function.
void printLine() {
std::cout << "----------------------------------------"
"----------------------------------------\n";
}
int main(
int argc,
const char* argv[]) {
(void)argc;
(void)argv;
std::cout << "sgpp::optimization example program started.\n\n";
void setVerbosity(int level)
Definition Printer.hpp:154
static Printer & getInstance()
Definition Printer.cpp:43
int main()
Definition densityMultiplication.cpp:22
Here, we define some parameters: objective function, dimensionality, B-spline degree, maximal number of grid points, and adaptivity.
const size_t d =
f.getNumberOfParameters();
const size_t p = 3;
const size_t N = 30;
const double gamma = 0.95;
base::ScalarFunction & f
Definition AugmentedLagrangian.cpp:75
First, we define a grid with modified B-spline basis functions and an iterative grid generator, which can generate the grid adaptively.
Grid with modified Bspline basis functions.
Definition ModBsplineGrid.hpp:21
Iterative grid generation based on Ritter/Novak's refinement criterion.
Definition IterativeGridGeneratorRitterNovak.hpp:31
std::unique_ptr< sgpp::base::Grid > grid(nullptr)
With the iterative grid generator, we generate adaptively a sparse grid.
printLine();
std::cout << "Generating grid...\n\n";
if (!gridGen.generate()) {
std::cout << "Grid generation failed, exiting.\n";
return 1;
}
Then, we hierarchize the function values to get hierarchical B-spline coefficients of the B-spline sparse grid interpolant \(\tilde{f}\colon [0, 1]^d \to \mathbb{R}\).
printLine();
std::cout << "Hierarchizing...\n\n";
if (!sleSolver.
solve(hierSLE, functionValues, coeffs)) {
std::cout << "Solving failed, exiting.\n";
return 1;
}
Linear system of the hierarchization in a sparse grid.
Definition HierarchisationSLE.hpp:74
Automatic choice of external linear solver.
Definition Auto.hpp:20
bool solve(SLE &system, DataVector &b, DataVector &x) const override
Definition Auto.cpp:46
We define the interpolant \(\tilde{f}\) and its gradient \(\nabla\tilde{f}\) for use with the gradient method (steepest descent). Of course, one can also use other optimization algorithms from sgpp::optimization::optimizer.
printLine();
std::cout << "Optimizing smooth interpolant...\n\n";
double fX0;
double ftX0;
Sparse grid interpolant gradient of a scalar-valued function.
Definition InterpolantScalarFunctionGradient.hpp:26
Sparse grid interpolant of a scalar-valued function.
Definition InterpolantScalarFunction.hpp:36
Gradient-based method of steepest descent.
Definition GradientDescent.hpp:23
The gradient method needs a starting point. We use a point of our adaptively generated sparse grid as starting point. More specifically, we use the point with the smallest (most promising) function value and save it in x0.
{
size_t x0Index =
std::distance(functionValues.getPointer(),
std::min_element(functionValues.getPointer(),
functionValues.getPointer() + functionValues.getSize()));
x0 = gridStorage.getCoordinates(gridStorage[x0Index]);
fX0 = functionValues[x0Index];
ftX0 = ft.eval(x0);
}
std::cout << "x0 = " << x0.toString() << "\n";
std::cout << "f(x0) = " << fX0 << ", ft(x0) = " << ftX0 << "\n\n";
Generic hash table based storage of grid points.
Definition HashGridStorage.hpp:42
We apply the gradient method and print the results.
gradientDescent.setStartingPoint(x0);
gradientDescent.optimize();
const double ftXOpt = gradientDescent.getOptimalValue();
const double fXOpt =
f.eval(xOpt);
std::cout <<
"\nxOpt = " << xOpt.
toString() <<
"\n";
std::cout << "f(xOpt) = " << fXOpt << ", ft(xOpt) = " << ftXOpt << "\n\n";
void toString(std::string &text) const
Writes the data stored in the DataVector into a string.
Definition DataVector.cpp:415
For comparison, we apply the classical gradient-free Nelder-Mead method directly to the objective function \(f\).
printLine();
std::cout << "Optimizing objective function (for comparison)...\n\n";
nelderMead.optimize();
const double fXOptNM = nelderMead.getOptimalValue();
const double ftXOptNM = ft.eval(xOptNM);
std::cout <<
"\nnxOptNM = " << xOptNM.
toString() <<
"\n";
std::cout << "f(xOptNM) = " << fXOptNM << ", ft(xOptNM) = " << ftXOptNM << "\n\n";
printLine();
std::cout << "\nsgpp::optimization example program terminated.\n";
return 0;
}
Gradient-free Nelder-Mead method.
Definition NelderMead.hpp:20
The example program outputs the following results:
sgpp::optimization example program started.
--------------------------------------------------------------------------------
Generating grid...
Adaptive grid generation (Ritter-Novak)...
100.0% (N = 29, k = 3)
Done in 3ms.
--------------------------------------------------------------------------------
Hierarchizing...
Solving linear system (automatic method)...
estimated nnz ratio: 59.8%
Solving linear system (Armadillo)...
constructing matrix (100.0%)
nnz ratio: 58.0%
solving with Armadillo
Done in 0ms.
Done in 1ms.
--------------------------------------------------------------------------------
Optimizing smooth interpolant...
x0 = [0.625, 0.75]
f(x0) = -1.81786, ft(x0) = -1.81786
Optimizing (gradient method)...
9 steps, f(x) = -2.000780
Done in 1ms.
xOpt = [0.589526, 0.673268]
f(xOpt) = -1.99999, ft(xOpt) = -2.00078
--------------------------------------------------------------------------------
Optimizing objective function (for comparison)...
Optimizing (Nelder-Mead)...
280 steps, f(x) = -2.000000
Done in 2ms.
xOptNM = [0.589049, 0.673198]
f(xOptNM) = -2, ft(xOptNM) = -2.00077
--------------------------------------------------------------------------------
sgpp::optimization example program terminated.
We see that both the gradient-based optimization of the smooth sparse grid interpolant and the gradient-free optimization of the objective function find reasonable approximations of the minimum, which lies at \((3\pi/16, 3\pi/14) \approx (0.58904862, 0.67319843)\).