CNLMath : Optimization : min_con_polytope
min_con_polytope

   more...
Minimizes a function of n variables subject to bounds on the variables using a direct search complex algorithm.
Synopsis
#include <imsl.h>
float *imsl_f_min_con_polytope (void fcn(), int n, float xlb[], float xub[], , 0)
The typedouble function is imsl_d_min_con_polytope.
Required Arguments
voidfcn (int n, float x[], float *f) (Input/Output)
User-supplied function, f(x), to be minimized.
Arguments
int n (Input)
The number of variables.
float x[] (Input)
Array of size n at which point the function is evaluated.
float *f (Output)
The function value at x.
int n (Input)
The number of variables.
float xlb[] (Input)
Array of size n containing the lower bounds on the variables.
float xub[] (Input)
Array of size n containing the upper bounds on the variables.
Return Value
A pointer to the solution x containing the best estimate for the minimum. To release this space, use imsl_free. If no solution can be computed, then NULL is returned.
Synopsis with Optional Arguments
#include <imsl.h>
float *imsl_f_min_con_polytope (void fcn(), int n, float xlb[], float xub[],
IMSL_XGUESS, float xguess[],
IMSL_XCPLX, float xcplx[],
IMSL_TOLERANCE, float ftol,
IMSL_REFLCOEF, float alpha,
IMSL_EXPNCOEF, float beta,
IMSL_CNTRCOEF, float gamma,
IMSL_MAX_FCN, int *maxfcn,
IMSL_FVALUE, float *fvalue,
IMSL_RETURN_USER, float x[],
IMSL_FCN_W_DATA, void fcn(), void *data,
0)
Optional Arguments
IMSL_XGUESS, float xguess[] (Input)
Array of size n containing an initial guess of the minimum.
Default: xguess = 0.0.
IMSL_XCPLX, float xcplx[] (Input)
Array of size 2n × n containing the 2 * n initial complex points. If the IMSL_XCPLX optional argument is used, then the initial guess must be specified in the first row of xcplx, and there is no need to provide optional argument IMSL_XGUESS. Thus, if both optional arguments IMSL_XCPLX and IMSL_XGUESS are supplied, IMSL_XGUESS will be ignored.
Default: xcplx is not used.
IMSL_TOLERANCE, float ftol (Input)
First convergence criterion: The algorithm stops when the relative error in the function values is less than ftol, i.e. when F(worst) – F(best) < ftol * (1 + |F(best)|) where F(worst) and F(best) are the function values of the current worst and best point, respectively.
Second convergence criterion: The algorithm stops when the standard deviation of the function values at the 2 * n current points is less than ftol. If the function terminates prematurely, try again with a smaller value for ftol.
Default: ftol = 1.0e-4 for single and 1.0e-8 for double precision.
IMSL_REFLCOEF, float alpha (Input)
The reflection coefficient. alpha must be greater than 0.
Default: alpha = 1.0.
IMSL_EXPNCOEF, float beta (Input)
The expansion coefficient. beta must be greater than 1.
Default: beta = 2.0.
IMSL_CNTRCOEF, float gamma (Input)
The contraction coefficient. gamma must be greater than 0 and less than 1.
Default: gamma = 0.5.
IMSL_MAX_FCN, int *maxfcn (Input/Output)
On input, the maximum allowed number of function evaluations. On output, the actual number of function evaluations needed.
Default: maxfcn = 300.
IMSL_FVALUE, float *fvalue (Output)
Function value at the computed solution.
IMSL_RETURN_USER, float x[] (Output)
User-supplied array of size n containing the computed solution.
IMSL_FCN_W_DATA, void fcn (int n, float x[], float *f, void *data), void *data (Input)
User supplied function to evaluate the function to be minimized, which also accepts a pointer to data that is supplied by the user. data is a pointer to the data to be passed to the user-supplied function. See Passing Data to User-Supplied Functionsin the Introduction of this manual for more details.
Description
The function imsl_f_min_con_polytope uses the complex method to find a minimum point of a function of n variables. The method is based on function comparison; no smoothness is assumed. It starts with an initial complex of 2n points , which is either user-supplied (using optional argument IMSL_XCPLX) or is otherwise, by default, randomly initialized. At each iteration, a new point is generated to replace the “worst” point xj, which has the largest function value among these 2n points, as described below.
Each iteration begins by determining the best and two worst points in the present complex, and then constructing a new “reflection” point xr by the formula
xr = c + α(c - xj)
where
is the centroid of the 2n - 1 best points and α (α > 0) is the reflection coefficient. (See optional argument IMSL_REFLCOEF). Depending on how the new reflection point xr compares with the existing complex points, the iteration proceeds as follows:
If xr is neither better than the best point nor worse than the second worst point, then xr replaces the worst point xj and, if neither of the stopping criteria (see below) is satisfied, a new iteration begins.
If xr is a best point, that is, if f(xr) ≤ f(xi) for i = 1, …, 2n, an expansion point xe is computed to see if an even better point can be obtained by moving further in the reflection direction:
xe = c + β(xr  c)
where β (β > 1) is called the expansion coefficient (see optional argument IMSL_EXPNCOEF), and worst point xj is replaced by the better of xe and xr and, if neither of the stopping criteria is satisfied, a new iteration begins.
If xr and xj are the two worst points, then the complex is contracted to try to get a better new contraction point xc:
where (0 <  < 1) is called the contraction coefficient. (See optional argument IMSL_CNTRCOEF. If the contraction step is successful (i.e. if min( f(xr), f(xj)) > f(xc)), then worst point xj is replaced by xc. If the contraction step is unsuccessful, then the complex is shrunk by moving the 2n  1 worst points halfway towards the current best point. Following the contraction step, if neither of the stopping criteria is satisfied, a new iteration begins.
Whenever the new point generated is beyond the bound, it will be projected onto the bound. If, at the end of an iteration, one of the following stopping criteria is satisfied, then the process ends with the best point returned as the optimum.
Criterion 1:
fworst fbest  ɛf (1. + |fbest|)
Criterion 2:
where fi = f(xi), fj = f(xj), and ɛf is a given tolerance. For a complete description, see Nelder and Mead (1965) or Gill et al. (1981).
Remarks
Since imsl_f_min_con_polytope uses only function-value information at each step to determine a new approximate minimum, it could be quite inefficient on smooth problems compared to other methods. Hence function imsl_f_min_con_polytope should be used only as a last resort. Briefly, a set of 2 * n points in ann -dimensional space is called a complex. The minimization process iterates by replacing the point with the largest function value by a new point with a smaller function value. The iteration continues until all the points cluster sufficiently close to a minimum.
Examples
Example 1
The problem
is solved with an initial guess (1.2, 1.0), and the solution is printed.
 
#include <imsl.h>
#include <stdio.h>
#include <math.h>
 
void fcn(int n, float x[], float *f);
 
int main() {
int n = 2;
float xguess[] = { -1.2, 1.0 };
float xlb[] = { -2.0, -1.0 };
float xub[] = { 0.5, 2.0 };
float ftol = 1.0e-7;
float *x = NULL, fvalue;
 
x = imsl_f_min_con_polytope(fcn, n, xlb, xub,
IMSL_XGUESS, xguess,
IMSL_TOLERANCE, ftol,
IMSL_FVALUE, &fvalue,
0);
 
printf("Solution x =(%.2f, %.2f)\n", x[0], x[1]);
printf("Function value at x = %0.2f\n", fvalue);
 
if (x) imsl_free(x);
}
 
void fcn(int n, float x[], float *f) {
*f = 100.0 * pow((x[1] - x[0] * x[0]), 2.0) + pow((1.0 - x[0]), 2.0);
}
Output
 
Solution x =(0.50, 0.25)
Function value at x = 0.25
 
Example 2
This example is intended to illustrate the use of optional arguments available for imsl_f_min_con_polytope, and how their use can affect the number of function calls needed to complete the optimization process. The same problem as in Example 1 is approached in five ways, including the use of a function penalty in an attempt to constrain the solution space. In each case, the number of function evaluations required is output.
Solution 1 uses xlb and xub to impose bounds, as in Example 1. Solution 2 through 5 use a function penalty to impose constraints, and to varying degrees make use of optional arguments IMSL_XCPLX, IMSL_REFLCOEF, IMSL_EXPNCOEF, and IMSL_CNTRCOEF.
Note that the actual number of function calls required to complete the minimization process may vary, depending on the computing platform and precision used.
 
#include <imsl.h>
#include <stdio.h>
#include <math.h>
 
#define N 2
 
void fcn(int n, double x[], double *f);
void print_results(double x[], double fvalue, int max_fcn);
void fcn_penalty(int n, double x[], double *f);
 
int main() {
int i, j, n = N;
double xguess[] = { -1.2, 1.0 };
double xlb[2], xub[2];
double xcplx[2 * N][N], xcplx0[2 * N][N];
int ibtype = 0;
double ftol = 1.0e-15;
double alpha, beta, gamma;
int max_fcn;
double *x = NULL, fvalue;
 
xcplx0[0][0] = xguess[0];
xcplx0[0][1] = xguess[1];
xcplx0[1][0] = 0.5;
xcplx0[1][1] = 2.0;
xcplx0[2][0] = -2.0;
xcplx0[2][1] = -1.0;
xcplx0[3][0] = 0.5;
xcplx0[3][1] = -1.0;
 
/*
* Solution 1:
* Constraints imposed with bounds.
*/
xlb[0] = -2.0;
xlb[1] = -1.0;
xub[0] = 0.5;
xub[1] = 2.0;
max_fcn = 500;
 
x = imsl_d_min_con_polytope(fcn, n, xlb, xub,
IMSL_XGUESS, xguess,
IMSL_TOLERANCE, ftol,
IMSL_FVALUE, &fvalue,
IMSL_MAX_FCN, &max_fcn,
0);
 
printf("Solution 1:\n-----------\n");
printf("Constraints imposed with simple bounds.\n");
printf(" -2.0 <= x[0] <= 0.5\n");
printf(" -1.0 <= x[1] <= 2.0\n");
printf("Using randomly generated initial complex and\n");
printf("default values for step-size coefficients.\n\n");
print_results(x, fvalue, max_fcn);
 
if (x) {
imsl_free(x);
x = NULL;
}
 
printf("\n-----------------------------------------------------\n");
printf("Note: In addition to simple bounds\n");
printf(" -2 <= x[0] <= 2\n");
printf(" -2 <= x[1] <= 2\n");
printf(" Solutions 2-5 include the following constraints\n");
printf(" imposed with function penalties:\n");
printf(" -2.0 <= x[0] <= 0.5\n");
printf(" -1.0 <= x[1] <= 2.0\n");
printf("-----------------------------------------------------\n\n");
 
/*
* Solution 2:
* User-supplied initial complex xcplx
* chosen to be within constraint boundaries
*/
xlb[0] = xlb[1] = -2.0;
xub[0] = xub[1] = 2.0;
 
/* Set the initial complex. */
for (i = 0; i < 2 * n; i++)
for (j = 0; j < n; j++)
xcplx[i][j] = xcplx0[i][j];
max_fcn = 1000;
 
x = imsl_d_min_con_polytope(fcn_penalty, n, xlb, xub,
IMSL_XGUESS, xguess,
IMSL_TOLERANCE, ftol,
IMSL_FVALUE, &fvalue,
IMSL_MAX_FCN, &max_fcn,
IMSL_XCPLX, xcplx,
0);
 
printf("Solution 2:\n-----------\n");
printf("Using user-specified intial complex and\n");
printf("default values for step-size coefficients.\n\n");
print_results(x, fvalue, max_fcn);
 
if (x) {
imsl_free(x);
x = NULL;
}
 
/*
* Solution 3:
* User-supplied initial complex xcplx
* chosen to be within constraint boundaries, and
* user-specified values for step-size coefficients.
*/
/* Reset the initial complex. */
for (i = 0; i < 2 * n; i++)
for (j = 0; j < n; j++)
xcplx[i][j] = xcplx0[i][j];
 
max_fcn = 1000;
alpha = 1.0;
beta = 3.1841776469083554;
gamma = 0.33464404002126491;
 
x = imsl_d_min_con_polytope(fcn_penalty, n, xlb, xub,
IMSL_XGUESS, xguess,
IMSL_TOLERANCE, ftol,
IMSL_FVALUE, &fvalue,
IMSL_MAX_FCN, &max_fcn,
IMSL_XCPLX, xcplx,
IMSL_REFLCOEF, alpha,
IMSL_EXPNCOEF, beta,
IMSL_CNTRCOEF, gamma,
0);
 
printf("Solution 3:\n-----------\n");
printf("Using user-specified intial complex and\n");
printf("user-specified values for step-size coefficients.\n");
printf(" alpha =%25.17e\n", alpha);
printf(" beta =%25.17e\n", beta);
printf(" gamma =%25.17e\n\n", gamma);
print_results(x, fvalue, max_fcn);
 
if (x) {
imsl_free(x);
x = NULL;
}
 
/*
* Solution 4:
* Using randomly generated initial complex and
* default values for step-size coefficients.
*/
max_fcn = 1000;
 
x = imsl_d_min_con_polytope(fcn_penalty, n, xlb, xub,
IMSL_XGUESS, xguess,
IMSL_TOLERANCE, ftol,
IMSL_FVALUE, &fvalue,
IMSL_MAX_FCN, &max_fcn,
0);
 
printf("Solution 4:\n-----------\n");
printf("Using randomly generated initial complex and\n");
printf("default values for step-size coefficients.\n\n");
print_results(x, fvalue, max_fcn);
 
if (x) {
imsl_free(x);
x = NULL;
}
 
/*
* Solution 5:
* Using randomly generated initial complex and
* user-supplied values for step-size coefficients.
*/
max_fcn = 1000;
beta = 18.204845270362373;
gamma = 0.31542073037934792;
 
x = imsl_d_min_con_polytope(fcn_penalty, n, xlb, xub,
IMSL_XGUESS, xguess,
IMSL_TOLERANCE, ftol,
IMSL_FVALUE, &fvalue,
IMSL_MAX_FCN, &max_fcn,
IMSL_REFLCOEF, alpha,
IMSL_EXPNCOEF, beta,
IMSL_CNTRCOEF, gamma,
0);
 
printf("Solution 5:\n-----------\n");
printf("Using randomly generated initial complex and\n");
printf("user-supplied values for step-size coefficients.\n");
printf(" alpha =%25.17e\n", alpha);
printf(" beta =%25.17e\n", beta);
printf(" gamma =%25.17e\n\n", gamma);
print_results(x, fvalue, max_fcn);
 
if (x) {
imsl_free(x);
x = NULL;
}
}
 
 
void fcn(int n, double x[], double *f) {
*f = 100.0 * pow((x[1] - x[0] * x[0]), 2.0) + pow((1.0 - x[0]), 2.0);
}
 
 
void fcn_penalty(int n, double x[], double *f) {
if (((-2.0 <= x[0]) && (x[0] <= 0.5)) &&
((-1.0 <= x[1]) && (x[1] <= 2.0))) {
*f = 100.0 * pow((x[1] - x[0] * x[0]), 2.0) + pow((1.0 - x[0]), 2.0);
}
else {
*f = 1000.0;
}
}
 
 
void print_results(double x[], double fvalue, int max_fcn){
printf("Results:\n");
printf("x = [%f, %f]\n", x[0], x[1]);
printf("fvalue = %f\n", fvalue);
printf("function calls = %d\n\n\n", max_fcn);
}
 
Output
Solution 1:
-----------
Constraints imposed with simple bounds.
-2.0 <= x[0] <= 0.5
-1.0 <= x[1] <= 2.0
Using randomly generated initial complex and
default values for step-size coefficients.
 
Results:
x = [0.500000, 0.250000]
fvalue = 0.250000
function calls = 226
 
 
 
-----------------------------------------------------
Note: In addition to simple bounds
-2 <= x[0] <= 2
-2 <= x[1] <= 2
Solutions 2-5 include the following constraints
imposed with function penalties:
-2.0 <= x[0] <= 0.5
-1.0 <= x[1] <= 2.0
-----------------------------------------------------
 
Solution 2:
-----------
Using user-specified intial complex and
default values for step-size coefficients.
 
Results:
x = [0.500000, 0.250000]
fvalue = 0.250000
function calls = 379
 
 
Solution 3:
-----------
Using user-specified intial complex and
user-specified values for step-size coefficients.
alpha = 1.00000000000000000e+000
beta = 3.18417764690835540e+000
gamma = 3.34644040021264910e-001
 
Results:
x = [0.500000, 0.250000]
fvalue = 0.250000
function calls = 323
 
 
Solution 4:
-----------
Using randomly generated initial complex and
default values for step-size coefficients.
 
Results:
x = [0.500000, 0.250000]
fvalue = 0.250000
function calls = 430
 
 
Solution 5:
-----------
Using randomly generated initial complex and
user-supplied values for step-size coefficients.
alpha = 1.00000000000000000e+000
beta = 1.82048452703623730e+001
gamma = 3.15420730379347920e-001
 
Results:
x = [0.500000, 0.250000]
fvalue = 0.250000
function calls = 294
 
Warning Errors
IMSL_MAX_FCN_EVALS
Maximum number of function evaluations, “max_fcn” = # exceeded.
Fatal Errors
IMSL_STOP_USER_FCN
Request from user supplied function to stop algorithm. User flag = "#".