This example is intended to demonstrate how to set up and solve a very simple problem. We will show you how to solve (unconstrained) Rosenbrock's function in two dimensions, i.e.,
minimize
For demonstration purposes, we will assume that there are no analytic derivatives available so we will use finite-difference approximations to the gradient and a quasi-Newton algorithm with a BFGS approximation to the Hessian. Recall that it is necessary to write C++ code for the main routine that sets up the problem and the algorithm and for the subroutines that initialize and evaluate the function. We step through the specifics below.
First include the necessary header files. Start with any C++/C header files that are needed. In this case, a bit of I/O done. The other two header files are OPT++ header files. NLF contains objects, data, and methods required for setting up the function/problem. OptQNewton contains the objects, data, and methods required for using an unconstrained quasi-Newton optimization method.
#include <fstream> #include "NLF.h" #include "OptQNewton.h" |
These two lines serve as the declarations of the pointers to the subroutines that initialize the problem and evaluate the objective function, respectively.
void init_rosen(int ndim, ColumnVector& x); void rosen(int ndim, const ColumnVector& x, double& fx, int& result); |
The next few lines complete the setup of the problem. Set the dimension of the problem. Create the nonlinear function object using the dimension of the problem and the pointers to the subroutines declared above. The FDNLF1 object has built-in finite-difference approximations to the gradient. Finally, initialize the problem using its initFcn method.
int main()
{
int ndim = 2;
FDNLF1 nlp(ndim, rosen, init_rosen);
nlp.initFcn(); |
Build a quasi-Newton algorithm object using the nonlinear problem that has just been created. The quasi-Newton algorithm will use BFGS updates to approximate the Hessian. In addition, set any of the algorithmic parameters to desired values. All parameters have default values, so it is not necessary to set them unless you have specific values you wish to use. In this example, we set the globalization strategy, the maximum number of function evaluations allowed, the function tolerance (used as a stopping criterion), and the name of the output file.
OptQNewton objfcn(&nlp); objfcn.setSearchStrategy(TrustRegion); objfcn.setMaxFeval(200); objfcn.setFcnTol(1.e-4); // The "0" in the second argument says to create a new file. A "1" // would signify appending to an existing file. if (!objfcn.setOutputFile("example1.out", 0)) cerr << "main: output file open failed" << endl; |
Now call the algorithm's optimize method to solve the problem.
objfcn.optimize(); |
Print out some summary information and clean up before exiting. The summary information is handy, but not necessary. The cleanup flushes the I/O buffers.
objfcn.printStatus("Solution from quasi-newton");
objfcn.cleanup();
} |
Now that the main routine is in place, we step through the code required for the initialization and evaluation of the function.
This section contains examples of the user-defined functions that are required. The first performs the initialization of the problem. The second performs the evaluation of the function.
First, include the necessary header files. In this case, we need the OPT++ header file, NLP, for some definitions.
#include "NLP.h" |
The subroutine that initializes the problem should perform any one-time tasks that are needed for the problem. One part of that is checking for error conditions in the setup. In this case, the dimension, ndim, can only take on a value of 2. Using "exit" is not the ideal way to deal with error conditions, but it serves well as an example.
void init_rosen (int ndim, ColumnVector& x)
{
if (ndim != 2)
exit (1); |
The initialization is also an ideal place to set the initial values of the optimization parameters, x. This can be hard coded, as done here, or it can be done in some other manner (e.g., reading them in from a file, the code for which should appear here).
// ColumnVectors are indexed from 1, and they use parentheses around // the index. x(1) = -1.2; x(2) = 1.0; } |
The last piece of code is a subroutine that will evaluate the function. In this problem, we are trying to find the minimum value of Rosenbrock's function, so it is necessary to write the code that compute the value of that function given some set of optimization parameters. Mathematically, Rosenbrock's function is:
The following code will compute the value of f(x).
First, some error checking and manipulation of the optimization parameters, x, are done.
void rosen(int ndim, const ColumnVector& x, double& fx, int& result)
{
double f1, f2, x1, x2;
if (ndim != 2)
exit (1);
x1 = x(1);
x2 = x(2);
f1 = (x2 - x1 * x1);
f2 = 1. - x1; |
Then the function value, fx, is computed, and the variable, result, is set to indicate that a function evaluation was performed.
fx = 100.* f1*f1 + f2*f2; result = NLPFunction; } |
On a more general note, this subroutine could serve as a wrapper to a C or Fortran subroutine. Similarly, it could make a system call to a completely independent executable. As long as the values of fx and result are set when all is said and done, it does not matter how the function value is computed.
Now that we have all of the code necessary to set up and solve Rosenbrock's function, give it a try!
Building and Running the Example
If you want to try running this example, the following steps should do the trick.
Next Section: Example 2: Nonlinear Interior-Point Method with General Constraints | Back to Setting up and Solving an Optimization Problem