CNLMath : Optimization : min_con_gen_lin
min_con_gen_lin

   more...
Minimizes a general objective function subject to linear equality/inequality constraints.
Synopsis
#include <imsl.h>
float *imsl_f_min_con_gen_lin (void fcn(), int nvar, int ncon,int neq, float a[], float b[], float xlb[], float xub[], ..., 0)
The type double function is imsl_d_min_con_gen_lin.
Required Arguments
void fcn (int n, float x[], float *f ) (Input/Output)
User-supplied function to evaluate the function to be minimized. Argument x is a vector of length n at which point the function is evaluated, and f contains the function value at x.
int nvar (Input)
Number of variables.
int ncon (Input)
Number of linear constraints (excluding simple bounds).
int neq (Input)
Number of linear equality constraints.
float a[] (Input)
Array of size ncon × nvar containing the equality constraint gradients in the first neq rows followed by the inequality constraint gradients.
float b[] (Input)
Array of size ncon containing the right-hand sides of the linear constraints. Specifically, the constraints on the variables xii = 0, nvar -1, are ak,0x0 +  + ak,nvar-1xnvar-1 = bk, k = 0, , neq - 1 and ak,0x0 + ak,nvar-1xnvar-1  bk, k = neq, , ncon ‑ 1. Note that the data that define the equality constraints come before the data of the inequalities.
float xlb[] (Input)
Array of length nvar containing the lower bounds on the variables; choose a very large negative value if a component should be unbounded below or set xub[i] = xub[i] to freeze the i-th variable. Specifically, these simple bounds are xlb[i xi, for i = 1, , nvar.
float xub[] (Input)
Array of length nvar containing the upper bounds on the variables; choose a very large positive value if a component should be unbounded above. Specifically, these simple bounds are xi  xub[i], for i = 1, nvar.
Return Value
A pointer to the solution x. 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_gen_lin (void fcn(), int nvar, int ncon, int a, float b, float xlb[], float xub[],
IMSL_XGUESS, float xguess[],
IMSL_GRADIENT, void gradient(),
IMSL_MAX_FCN, int max_fcn,
IMSL_NUMBER_ACTIVE_CONSTRAINTS, int *nact,
IMSL_ACTIVE_CONSTRAINTS, int **iact,
IMSL_ACTIVE_CONSTRAINTS_USER, int *iact_user,
IMSL_LAGRANGE_MULTIPLIERS, float **lagrange,
IMSL_LAGRANGE_MULTIPLIERS_USER, float *lagrange_user,
IMSL_TOLERANCE, float tolerance,
IMSL_OBJ, float *obj,
IMSL_RETURN_USER, float x[],
IMSL_FCN_W_DATA, void fcn(), void *data,
IMSL_GRADIENT_W_DATA, void gradient(),void *data,
0)
Optional Arguments
IMSL_XGUESS, float xguess[] (Input)
Array with n components containing an initial guess.
Default: xguess = 0
IMSL_GRADIENT, void gradient (int n, float x[], float g[]) (Input)
User-supplied function to compute the gradient at the point x, where x is a vector of length n, and g is the vector of length n containing the values of the gradient of the objective function.
IMSL_MAX_FCN, int max_fcn (Input)
Maximum number of function evaluations.
Default: max_fcn = 400
IMSL_NUMBER_ACTIVE_CONSTRAINTS, int *nact (Output)
Final number of active constraints.
IMSL_ACTIVE_CONSTRAINTS, int **iact (Output)
The address of a pointer to an int, which on exit, points to an array containing the nact indices of the final active constraints.
IMSL_ACTIVE_CONSTRAINTS_USER, int × iact_user (Output)
A user-supplied array of length at least ncon + 2
containing the indices of the final active constraints in the first nact locations.
IMSL_LAGRANGE_MULTIPLIERS, float **lagrange (Output)
The address of a pointer, which on exit, points to an array containing the Lagrange multiplier estimates of the final active constraints in the first nact locations.
IMSL_LAGRANGE_MULTIPLIERS_USER, float *lagrange_user (Output)
A user-supplied array of length at least nvar containing the Lagrange multiplier estimates of the final active constraints in the first nact locations.
IMSL_TOLERANCE, float tolerance (Input)
The nonnegative tolerance on the first order conditions at the calculated solution.
Default: tolerance =  , where ɛ is machine epsilon
IMSL_OBJ, float *obj (Output)
The value of the objective function.
IMSL_RETURN_USER, float x[] (Output)
User-supplied array with nvar components 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 compute the value of 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 Functions in the introduction to this manual for more details.
IMSL_GRADIENT_W_DATA, void gradient (int n, float x[], float g[], void *data), void *data (Input)
User-supplied function to compute the gradient at the point x, 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 Functions in the introduction to this manual for more details.
Description
The function imsl_f_min_con_gen_lin is based on M.J.D. Powell’s TOLMIN, which solves linearly constrained optimization problems, i.e., problems of the form
min f (x)
subject to
given the vectors b1, b2, xl ,and xu and the matrices A1 and A2.
The algorithm starts by checking the equality constraints for inconsistency and redundancy. If the equality constraints are consistent, the method will revise x0, the initial guess, to satisfy
Next, x0 is adjusted to satisfy the simple bounds and inequality constraints. This is done by solving a sequence of quadratic programming subproblems to minimize the sum of the constraint or bound violations.
Now, for each iteration with a feasible xk, let Jk be the set of indices of inequality constraints that have small residuals. Here, the simple bounds are treated as inequality constraints. Let Ik be the set of indices of active constraints. The following quadratic programming problem
subject to
is solved to get (dk, λk) where aj is a row vector representing either a constraint in A1 or A2 or a bound constraint on x. In the latter case, the aj =ei for the bound constraint xi  (xu)i and aj = -ei for the constraint -xi  (xl)i. Here, ei is a vector with 1 as the i-th component, and zeros elsewhere. Variables λk are the Lagrange multipliers, and Bk is a positive definite approximation to the second derivative 2 f(xk).
After the search direction dk is obtained, a line search is performed to locate a better point. The new point xk+1 =xk +αkdk has to satisfy the conditions
and
The main idea in forming the set Jk is that, if any of the equality constraints restricts the step-length αk, then its index is not in Jk. Therefore, small steps are likely to be avoided.
Finally, the second derivative approximation BK, is updated by the BFGS formula, if the condition
holds. Let xk xk+1, and start another iteration.
The iteration repeats until the stopping criterion
is satisfied. Here is the supplied tolerance. For more details, see Powell (1988, 1989).
Since a finite difference method is used to approximate the gradient for some single precision calculations, an inaccurate estimate of the gradient may cause the algorithm to terminate at a non­critical point. In such cases, high precision arithmetic is recommended. Also, if the gradient can be easily provided, the option IMSL_GRADIENT should be used.
On some platforms, imsl_f_min_con_gen_lin can evaluate the user-supplied functions fcn and gradient in parallel. This is done only if the function imsl_omp_options is called to flag user-defined functions as thread-safe. A function is thread-safe if there are no dependencies between calls. Such dependencies are usually the result of writing to global or static variables
Examples
Example 1
In this example, the problem
is solved.
 
#include <imsl.h>
 
int main()
{
void fcn(int, float *, float *);
int neq = 2;
int ncon = 2;
int nvar = 5;
 
float a[] = {1.0, 1.0, 1.0, 1.0, 1.0,
0.0, 0.0, 1.0, -2.0, -2.0};
float b[] = {5.0, -3.0};
float xlb[] = {0.0, 0.0, 0.0, 0.0, 0.0};
float xub[] = {10.0, 10.0, 10.0, 10.0, 10.0};
float *x;
 
imsl_omp_options(IMSL_SET_FUNCTIONS_THREAD_SAFE, 1, 0);
x = imsl_f_min_con_gen_lin(fcn, nvar, ncon, neq, a, b, xlb, xub,
0);
 
imsl_f_write_matrix("Solution", 1, nvar, x, 0);
}
 
void fcn(int n, float *x, float *f)
{
*f = x[0]*x[0] + x[1]*x[1] + x[2]*x[2] + x[3]*x[3] + x[4]*x[4]
- 2.0*x[1]*x[2] - 2.0*x[3] * x[4] - 2.0*x[0];
}
Output
 
Solution
1 2 3 4 5
1 1 1 1 1
Example 2
In this example, the problem from Schittkowski (1987)
is solved with an initial guess of x0 =10, x1 =10 and x2 =10.
 
#include <imsl.h>
 
int main()
{
void fcn(int, float *, float *);
void grad(int, float *, float *);
int neq = 0;
int ncon = 2;
int nvar = 3;
int lda = 2;
float obj, x[3];
float a[] = {-1.0, -2.0, -2.0,
1.0, 2.0, 2.0};
float xlb[] = {0.0, 0.0, 0.0};
float xub[] = {20.0, 11.0, 42.0};
float xguess[] = {10.0, 10.0, 10.0};
float b[] = {0.0, 72.0};
 
 
imsl_omp_options(IMSL_SET_FUNCTIONS_THREAD_SAFE, 1, 0);
 
imsl_f_min_con_gen_lin(fcn, nvar, ncon, neq, a, b, xlb, xub,
IMSL_GRADIENT, grad,
IMSL_XGUESS, xguess,
IMSL_OBJ, &obj,
IMSL_RETURN_USER, x,
0);
 
imsl_f_write_matrix("Solution", 1, nvar, x, 0);
printf("Objective value = %f\n", obj);
}
 
void fcn(int n, float *x, float *f)
{
*f = -x[0] * x[1] * x[2];
}
 
void grad(int n, float *x, float *g)
{
g[0] = -x[1]*x[2];
g[1] = -x[0]*x[2];
g[2] = -x[0]*x[1];
}
Output
 
Solution
1 2 3
20 11 15
Objective value = -3300.000000
Fatal Errors
IMSL_STOP_USER_FCN
Request from user supplied function to stop algorithm. User flag = "#".