GALAHAD LPA package#
purpose#
The lpa
package uses the simplex method to solve a
given linear program (LP).
The aim is to minimize the linear objective function
See Section 4 of $GALAHAD/doc/lpa.pdf for a brief description of the method employed and other details.
N.B. The package is simply a sophisticated interface to the HSL package
LA04, and requires that a user has obtained the latter. LA04 is not
included in GALAHAD but is available without charge to recognised
academics, see http://www.hsl.rl.ac.uk/catalogue/la04.html. If LA04 is
unavailable, the interior- point linear programming package lpb
is an alternative.
terminology#
Any required solution \(x\) necessarily satisfies the primal optimality conditions
The so-called dual to this problem is another linear program
method#
The bulk of the work is peformed by the HSL package LA04. The main subbroutine from this package requires that the input problem be transformed into the ``standard form’’
Once a selection of \(m'\) independent (non-basic) variables is made, the enlarged constraints determine the remaining \(m'\) dependent ({basic) variables. The simplex method is a scheme for systematically adjusting the choice of basic and non-basic variables until a set which defines an optimal solution of the standard-form LP is obtained. Each iteration of the simplex method requires the solution of a number of sets of linear equations whose coefficient matrix is the basis matrix \(B\), made up of the columns of \([A' \;\; I]\) corresponding to the basic variables, or its transpose \(B^T\). As the basis matrices for consecutive iterations are closely related, it is normally advantageous to update (rather than recompute) their factorizations as the computation proceeds. If an initial basis is not provided by the user, a set of basic variables which provide a (permuted) triangular basis matrix is found by the simple crash algorithm of Gould and Reid (1989), and initial steepest-edge weights are calculated.
Phases one (finding a feasible solution) and two (solving the standard-form LP) of the simplex method are applied, as appropriate, with the choice of entering variable as described by Goldfarb and Reid (1977) and the choice of leaving variable as proposed by Harris (1973). Refactorizations of the basis matrix are performed whenever doing so will reduce the average iteration time or there is insufficient memory for its factors. The reduced cost for the entering variable is computed afresh. If it is found to be of a different sign from the recurred value or more than 10% different in magnitude, a fresh computation of all the reduced costs is performed. Details of the factorization and updating procedures are given by Reid (1982). Iterative refinement is encouraged for the basic solution and for the reduced costs after each factorization of the basis matrix and when they are recomputed at the end of phase 1.
references#
D. Goldfarb and J. K. Reid, ``A practicable steepest-edge simplex algorithm’’. Mathematical Programming 12 (1977) 361-371.
N. I. M. Gould and J. K. Reid, ``New crash procedures for large systems of linear constraints’’. Mathematical Programming 45 (1989) 475-501.
P. M. J. Harris, ``Pivot selection methods of the Devex LP code’’. Mathematical Programming 5 (1973) 1-28.
J. K. Reid, ``A sparsity-exploiting variant of the Bartels-Golub decomposition for linear-programming bases’’. Mathematical Programming 24 (1982) 55-69.
matrix storage#
The unsymmetric \(m\) by \(n\) matrix \(A\) may be presented and stored in a variety of convenient input formats.
Dense storage format: The matrix \(A\) is stored as a compact dense matrix by rows, that is, the values of the entries of each row in turn are stored in order within an appropriate real one-dimensional array. In this case, component \(n \ast i + j\) of the storage array A_val will hold the value \(A_{ij}\) for \(0 \leq i \leq m-1\), \(0 \leq j \leq n-1\). The string A_type = ‘dense’ should be specified.
Dense by columns storage format: The matrix \(A\) is stored as a compact dense matrix by columns, that is, the values of the entries of each column in turn are stored in order within an appropriate real one-dimensional array. In this case, component \(m \ast j + i\) of the storage array A_val will hold the value \(A_{ij}\) for \(0 \leq i \leq m-1\), \(0 \leq j \leq n-1\). The string A_type = ‘dense_by_columns’ should be specified.
Sparse co-ordinate storage format: Only the nonzero entries of the matrices are stored. For the \(l\)-th entry, \(0 \leq l \leq ne-1\), of \(A\), its row index i, column index j and value \(A_{ij}\), \(0 \leq i \leq m-1\), \(0 \leq j \leq n-1\), are stored as the \(l\)-th components of the integer arrays A_row and A_col and real array A_val, respectively, while the number of nonzeros is recorded as A_ne = \(ne\). The string A_type = ‘coordinate’should be specified.
Sparse row-wise storage format: Again only the nonzero entries are stored, but this time they are ordered so that those in row i appear directly before those in row i+1. For the i-th row of \(A\) the i-th component of the integer array A_ptr holds the position of the first entry in this row, while A_ptr(m) holds the total number of entries. The column indices j, \(0 \leq j \leq n-1\), and values \(A_{ij}\) of the nonzero entries in the i-th row are stored in components l = A_ptr(i), \(\ldots\), A_ptr(i+1)-1, \(0 \leq i \leq m-1\), of the integer array A_col, and real array A_val, respectively. For sparse matrices, this scheme almost always requires less storage than its predecessor. The string A_type = ‘sparse_by_rows’ should be specified.
Sparse column-wise storage format: Once again only the nonzero entries are stored, but this time they are ordered so that those in column j appear directly before those in column j+1. For the j-th column of \(A\) the j-th component of the integer array A_ptr holds the position of the first entry in this column, while A_ptr(n) holds the total number of entries. The row indices i, \(0 \leq i \leq m-1\), and values \(A_{ij}\) of the nonzero entries in the j-th columns are stored in components l = A_ptr(j), \(\ldots\), A_ptr(j+1)-1, \(0 \leq j \leq n-1\), of the integer array A_row, and real array A_val, respectively. As before, for sparse matrices, this scheme almost always requires less storage than the co-ordinate format. The string A_type = ‘sparse_by_columns’ should be specified.
introduction to function calls#
To solve a given problem, functions from the lpa package must be called in the following order:
To solve a given problem, functions from the lpa package must be called in the following order:
lpa_initialize - provide default control parameters and set up initial data structures
lpa_read_specfile (optional) - override control values by reading replacement values from a file
lpa_import - set up problem data structures and fixed values
lpa_reset_control (optional) - possibly change control parameters if a sequence of problems are being solved
lpa_solve_lp - solve the linear program
lpa_information (optional) - recover information about the solution and solution process
lpa_terminate - deallocate data structures
See the examples section for illustrations of use.
callable functions#
overview of functions provided#
// typedefs typedef float spc_; typedef double rpc_; typedef int ipc_; // structs struct lpa_control_type; struct lpa_inform_type; struct lpa_time_type; // function calls void lpa_initialize(void **data, struct lpa_control_type* control, ipc_ *status); void lpa_read_specfile(struct lpa_control_type* control, const char specfile[]); void lpa_import( struct lpa_control_type* control, void **data, ipc_ *status, ipc_ n, ipc_ m, const char A_type[], ipc_ A_ne, const ipc_ A_row[], const ipc_ A_col[], const ipc_ A_ptr[] ); void lpa_reset_control( struct lpa_control_type* control, void **data, ipc_ *status ); void lpa_solve_lp( void **data, ipc_ *status, ipc_ n, ipc_ m, const rpc_ g[], const rpc_ f, ipc_ a_ne, const rpc_ A_val[], const rpc_ c_l[], const rpc_ c_u[], const rpc_ x_l[], const rpc_ x_u[], rpc_ x[], rpc_ c[], rpc_ y[], rpc_ z[], ipc_ x_stat[], ipc_ c_stat[] ); void lpa_information(void **data, struct lpa_inform_type* inform, ipc_ *status); void lpa_terminate( void **data, struct lpa_control_type* control, struct lpa_inform_type* inform );
typedefs#
typedef float spc_
spc_
is real single precision
typedef double rpc_
rpc_
is the real working precision used, but may be changed to float
by
defining the preprocessor variable SINGLE
.
typedef int ipc_
ipc_
is the default integer word length used, but may be changed to
int64_t
by defining the preprocessor variable INTEGER_64
.
function calls#
void lpa_initialize(void **data, struct lpa_control_type* control, ipc_ *status)
Set default control values and initialize private data
Parameters:
data |
holds private internal data |
control |
is a struct containing control information (see lpa_control_type) |
status |
is a scalar variable of type ipc_, that gives the exit status from the package. Possible values are (currently):
|
void lpa_read_specfile(struct lpa_control_type* control, const char specfile[])
Read the content of a specification file, and assign values associated with given keywords to the corresponding control parameters. An in-depth discussion of specification files is available, and a detailed list of keywords with associated default values is provided in $GALAHAD/src/lpa/LPA.template. See also Table 2.1 in the Fortran documentation provided in $GALAHAD/doc/lpa.pdf for a list of how these keywords relate to the components of the control structure.
Parameters:
control |
is a struct containing control information (see lpa_control_type) |
specfile |
is a character string containing the name of the specification file |
void lpa_import( struct lpa_control_type* control, void **data, ipc_ *status, ipc_ n, ipc_ m, const char A_type[], ipc_ A_ne, const ipc_ A_row[], const ipc_ A_col[], const ipc_ A_ptr[] )
Import problem data into internal storage prior to solution.
Parameters:
control |
is a struct whose members provide control paramters for the remaining prcedures (see lpa_control_type) |
data |
holds private internal data |
status |
is a scalar variable of type ipc_, that gives the exit status from the package. Possible values are:
|
n |
is a scalar variable of type ipc_, that holds the number of variables. |
m |
is a scalar variable of type ipc_, that holds the number of general linear constraints. |
A_type |
is a one-dimensional array of type char that specifies the unsymmetric storage scheme used for the constraint Jacobian, \(A\). It should be one of ‘coordinate’, ‘sparse_by_rows’ or ‘dense; lower or upper case variants are allowed. |
A_ne |
is a scalar variable of type ipc_, that holds the number of entries in \(A\) in the sparse co-ordinate storage scheme. It need not be set for any of the other schemes. |
A_row |
is a one-dimensional array of size A_ne and type ipc_, that holds the row indices of \(A\) in the sparse co-ordinate storage scheme. It need not be set for any of the other schemes, and in this case can be NULL. |
A_col |
is a one-dimensional array of size A_ne and type ipc_, that holds the column indices of \(A\) in either the sparse co-ordinate, or the sparse row-wise storage scheme. It need not be set when the dense or diagonal storage schemes are used, and in this case can be NULL. |
A_ptr |
is a one-dimensional array of size n+1 and type ipc_, that holds the starting position of each row of \(A\), as well as the total number of entries, in the sparse row-wise storage scheme. It need not be set when the other schemes are used, and in this case can be NULL. |
void lpa_reset_control( struct lpa_control_type* control, void **data, ipc_ *status )
Reset control parameters after import if required.
Parameters:
control |
is a struct whose members provide control paramters for the remaining prcedures (see lpa_control_type) |
data |
holds private internal data |
status |
is a scalar variable of type ipc_, that gives the exit status from the package. Possible values are:
|
void lpa_solve_lp( void **data, ipc_ *status, ipc_ n, ipc_ m, const rpc_ g[], const rpc_ f, ipc_ a_ne, const rpc_ A_val[], const rpc_ c_l[], const rpc_ c_u[], const rpc_ x_l[], const rpc_ x_u[], rpc_ x[], rpc_ c[], rpc_ y[], rpc_ z[], ipc_ x_stat[], ipc_ c_stat[] )
Solve the linear program.
Parameters:
data |
holds private internal data |
status |
is a scalar variable of type ipc_, that gives the entry and exit status from the package. Possible exit values are:
|
n |
is a scalar variable of type ipc_, that holds the number of variables |
m |
is a scalar variable of type ipc_, that holds the number of general linear constraints. |
g |
is a one-dimensional array of size n and type rpc_, that holds the linear term \(g\) of the objective function. The j-th component of g, j = 0, … , n-1, contains \(g_j\). |
f |
is a scalar of type rpc_, that holds the constant term \(f\) of the objective function. |
a_ne |
is a scalar variable of type ipc_, that holds the number of entries in the constraint Jacobian matrix \(A\). |
A_val |
is a one-dimensional array of size a_ne and type rpc_, that holds the values of the entries of the constraint Jacobian matrix \(A\) in any of the available storage schemes. |
c_l |
is a one-dimensional array of size m and type rpc_, that holds the lower bounds \(c^l\) on the constraints \(A x\). The i-th component of c_l, i = 0, … , m-1, contains \(c^l_i\). |
c_u |
is a one-dimensional array of size m and type rpc_, that holds the upper bounds \(c^l\) on the constraints \(A x\). The i-th component of c_u, i = 0, … , m-1, contains \(c^u_i\). |
x_l |
is a one-dimensional array of size n and type rpc_, that holds the lower bounds \(x^l\) on the variables \(x\). The j-th component of x_l, j = 0, … , n-1, contains \(x^l_j\). |
x_u |
is a one-dimensional array of size n and type rpc_, that holds the upper bounds \(x^l\) on the variables \(x\). The j-th component of x_u, j = 0, … , n-1, contains \(x^l_j\). |
x |
is a one-dimensional array of size n and type rpc_, that holds the values \(x\) of the optimization variables. The j-th component of x, j = 0, … , n-1, contains \(x_j\). |
c |
is a one-dimensional array of size m and type rpc_, that holds the residual \(c(x)\). The i-th component of c, i = 0, … , m-1, contains \(c_i(x)\). |
y |
is a one-dimensional array of size n and type rpc_, that holds the values \(y\) of the Lagrange multipliers for the general linear constraints. The j-th component of y, i = 0, … , m-1, contains \(y_i\). |
z |
is a one-dimensional array of size n and type rpc_, that holds the values \(z\) of the dual variables. The j-th component of z, j = 0, … , n-1, contains \(z_j\). |
x_stat |
is a one-dimensional array of size n and type ipc_, that gives the optimal status of the problem variables. If x_stat(j) is negative, the variable \(x_j\) most likely lies on its lower bound, if it is positive, it lies on its upper bound, and if it is zero, it lies between its bounds. |
c_stat |
is a one-dimensional array of size m and type ipc_, that gives the optimal status of the general linear constraints. If c_stat(i) is negative, the constraint value \(a_i^Tx\) most likely lies on its lower bound, if it is positive, it lies on its upper bound, and if it is zero, it lies between its bounds. |
void lpa_information(void **data, struct lpa_inform_type* inform, ipc_ *status)
Provides output information
Parameters:
data |
holds private internal data |
inform |
is a struct containing output information (see lpa_inform_type) |
status |
is a scalar variable of type ipc_, that gives the exit status from the package. Possible values are (currently):
|
void lpa_terminate( void **data, struct lpa_control_type* control, struct lpa_inform_type* inform )
Deallocate all internal private storage
Parameters:
data |
holds private internal data |
control |
is a struct containing control information (see lpa_control_type) |
inform |
is a struct containing output information (see lpa_inform_type) |
available structures#
lpa_control_type structure#
#include <galahad_lpa.h> struct lpa_control_type { // components bool f_indexing; ipc_ error; ipc_ out; ipc_ print_level; ipc_ start_print; ipc_ stop_print; ipc_ maxit; ipc_ max_iterative_refinements; ipc_ min_real_factor_size; ipc_ min_integer_factor_size; ipc_ random_number_seed; ipc_ sif_file_device; ipc_ qplib_file_device; rpc_ infinity; rpc_ tol_data; rpc_ feas_tol; rpc_ relative_pivot_tolerance; rpc_ growth_limit; rpc_ zero_tolerance; rpc_ change_tolerance; rpc_ identical_bounds_tol; rpc_ cpu_time_limit; rpc_ clock_time_limit; bool scale; bool dual; bool warm_start; bool steepest_edge; bool space_critical; bool deallocate_error_fatal; bool generate_sif_file; bool generate_qplib_file; char sif_file_name[31]; char qplib_file_name[31]; char prefix[31]; };
detailed documentation#
control derived type as a C struct
components#
bool f_indexing
use C or Fortran sparse matrix indexing
ipc_ error
error and warning diagnostics occur on stream error
ipc_ out
general output occurs on stream out
ipc_ print_level
the level of output required is specified by print_level (>= 2 turns on LA04 output)
ipc_ start_print
any printing will start on this iteration
ipc_ stop_print
any printing will stop on this iteration
ipc_ maxit
at most maxit inner iterations are allowed
ipc_ max_iterative_refinements
maximum number of iterative refinements allowed
ipc_ min_real_factor_size
initial size for real array for the factors and other data
ipc_ min_integer_factor_size
initial size for integer array for the factors and other data
ipc_ random_number_seed
the initial seed used when generating random numbers
ipc_ sif_file_device
specifies the unit number to write generated SIF file describing the current problem
ipc_ qplib_file_device
specifies the unit number to write generated QPLIB file describing the current problem
rpc_ infinity
any bound larger than infinity in modulus will be regarded as infinite
rpc_ tol_data
the tolerable relative perturbation of the data (A,g,..) defining the problem
rpc_ feas_tol
any constraint violated by less than feas_tol will be considered to be satisfied
rpc_ relative_pivot_tolerance
pivot threshold used to control the selection of pivot elements in the matrix factorization. Any potential pivot which is less than the largest entry in its row times the threshold is excluded as a candidate
rpc_ growth_limit
limit to control growth in the upated basis factors. A refactorization occurs if the growth exceeds this limit
rpc_ zero_tolerance
any entry in the basis smaller than this is considered zero
rpc_ change_tolerance
any solution component whose change is smaller than a tolerence times the largest change may be considered to be zero
rpc_ identical_bounds_tol
any pair of constraint bounds (c_l,c_u) or (x_l,x_u) that are closer than identical_bounds_tol will be reset to the average of their values
rpc_ cpu_time_limit
the maximum CPU time allowed (-ve means infinite)
rpc_ clock_time_limit
the maximum elapsed clock time allowed (-ve means infinite)
bool scale
if .scale is true, the problem will be automatically scaled prior to solution. This may improve computation time and accuracy
bool dual
should the dual problem be solved rather than the primal?
bool warm_start
should a warm start using the data in C_stat and X_stat be attempted?
bool steepest_edge
should steepest-edge weights be used to detetrmine the variable leaving the basis?
bool space_critical
if .space_critical is true, every effort will be made to use as little space as possible. This may result in longer computation time
bool deallocate_error_fatal
if .deallocate_error_fatal is true, any array/pointer deallocation error will terminate execution. Otherwise, computation will continue
bool generate_sif_file
if .generate_sif_file is .true. if a SIF file describing the current problem is to be generated
bool generate_qplib_file
if .generate_qplib_file is .true. if a QPLIB file describing the current problem is to be generated
char sif_file_name[31]
name of generated SIF file containing input problem
char qplib_file_name[31]
name of generated QPLIB file containing input problem
char prefix[31]
all output lines will be prefixed by .prefix(2:LEN(TRIM(.prefix))-1) where .prefix contains the required string enclosed in quotes, e.g. “string” or ‘string’
lpa_time_type structure#
#include <galahad_lpa.h> struct lpa_time_type { // components rpc_ total; rpc_ preprocess; rpc_ clock_total; rpc_ clock_preprocess; };
detailed documentation#
time derived type as a C struct
components#
rpc_ total
the total CPU time spent in the package
rpc_ preprocess
the CPU time spent preprocessing the problem
rpc_ clock_total
the total clock time spent in the package
rpc_ clock_preprocess
the clock time spent preprocessing the problem
lpa_inform_type structure#
#include <galahad_lpa.h> struct lpa_inform_type { // components ipc_ status; ipc_ alloc_status; char bad_alloc[81]; ipc_ iter; ipc_ la04_job; ipc_ la04_job_info; rpc_ obj; rpc_ primal_infeasibility; bool feasible; rpc_ RINFO[40]; struct lpa_time_type time; struct rpd_inform_type rpd_inform; };
detailed documentation#
inform derived type as a C struct
components#
ipc_ status
return status. See LPA_solve for details
ipc_ alloc_status
the status of the last attempted allocation/deallocation
char bad_alloc[81]
the name of the array for which an allocation/deallocation error occurred
ipc_ iter
the total number of iterations required
ipc_ la04_job
the final value of la04’s job argument
ipc_ la04_job_info
any extra information from an unsuccessfull call to LA04 (LA04’s RINFO(35)
rpc_ obj
the value of the objective function at the best estimate of the solution determined by LPA_solve
rpc_ primal_infeasibility
the value of the primal infeasibility
bool feasible
is the returned “solution” feasible?
rpc_ RINFO[40]
the information array from LA04
struct lpa_time_type time
timings (see above)
struct rpd_inform_type rpd_inform
inform parameters for RPD
example calls#
This is an example of how to use the package to solve a linear program; the code is available in $GALAHAD/src/lpa/C/lpat.c . A variety of supported constraint matrix storage formats are shown.
Notice that C-style indexing is used, and that this is flagged by setting
control.f_indexing
to false
. The floating-point type rpc_
is set in galahad_precision.h
to double
by default, but to float
if the preprocessor variable SINGLE
is defined. Similarly, the integer
type ipc_
from galahad_precision.h
is set to int
by default,
but to int64_t
if the preprocessor variable INTEGER_64
is defined.
/* lpat.c */
/* Full test for the LPA C interface using C sparse matrix indexing */
#include <stdio.h>
#include <math.h>
#include "galahad_precision.h"
#include "galahad_cfunctions.h"
#include "galahad_lpa.h"
int main(void) {
// Derived types
void *data;
struct lpa_control_type control;
struct lpa_inform_type inform;
// Set problem data
ipc_ n = 3; // dimension
ipc_ m = 2; // number of general constraints
rpc_ g[] = {0.0, 2.0, 0.0}; // linear term in the objective
rpc_ f = 1.0; // constant term in the objective
ipc_ A_ne = 4; // Jacobian elements
ipc_ A_row[] = {0, 0, 1, 1}; // row indices
ipc_ A_col[] = {0, 1, 1, 2}; // column indices
ipc_ A_ptr[] = {0, 2, 4}; // row pointers
rpc_ A_val[] = {2.0, 1.0, 1.0, 1.0 }; // values
rpc_ c_l[] = {1.0, 2.0}; // constraint lower bound
rpc_ c_u[] = {2.0, 2.0}; // constraint upper bound
rpc_ x_l[] = {-1.0, - INFINITY, - INFINITY}; // variable lower bound
rpc_ x_u[] = {1.0, INFINITY, 2.0}; // variable upper bound
// Set output storage
rpc_ c[m]; // constraint values
ipc_ x_stat[n]; // variable status
ipc_ c_stat[m]; // constraint status
char st = ' ';
ipc_ status;
printf(" C sparse matrix indexing\n\n");
printf(" basic tests of lp storage formats\n\n");
for( ipc_ d=1; d <= 3; d++){
// Initialize LPA
lpa_initialize( &data, &control, &status );
// Set user-defined control options
control.f_indexing = false; // C sparse matrix indexing
// Start from 0
rpc_ x[] = {0.0,0.0,0.0};
rpc_ y[] = {0.0,0.0};
rpc_ z[] = {0.0,0.0,0.0};
switch(d){
case 1: // sparse co-ordinate storage
st = 'C';
lpa_import( &control, &data, &status, n, m,
"coordinate", A_ne, A_row, A_col, NULL );
lpa_solve_lp( &data, &status, n, m, g, f,
A_ne, A_val, c_l, c_u, x_l, x_u, x, c, y, z,
x_stat, c_stat );
break;
printf(" case %1" i_ipc_ " break\n",d);
case 2: // sparse by rows
st = 'R';
lpa_import( &control, &data, &status, n, m,
"sparse_by_rows", A_ne, NULL, A_col, A_ptr );
lpa_solve_lp( &data, &status, n, m, g, f,
A_ne, A_val, c_l, c_u, x_l, x_u, x, c, y, z,
x_stat, c_stat );
break;
case 3: // dense
st = 'D';
ipc_ A_dense_ne = 6; // number of elements of A
rpc_ A_dense[] = {2.0, 1.0, 0.0, 0.0, 1.0, 1.0};
lpa_import( &control, &data, &status, n, m,
"dense", A_ne, NULL, NULL, NULL );
lpa_solve_lp( &data, &status, n, m, g, f,
A_dense_ne, A_dense, c_l, c_u, x_l, x_u,
x, c, y, z, x_stat, c_stat );
break;
}
lpa_information( &data, &inform, &status );
if(inform.status == 0){
printf("%c:%6" i_ipc_ " iterations. Optimal objective value = %5.2f status = %1" i_ipc_ "\n",
st, inform.iter, inform.obj, inform.status);
}else{
printf("%c: LPA_solve exit status = %1" i_ipc_ "\n", st, inform.status);
}
//printf("x: ");
//for( ipc_ i = 0; i < n; i++) printf("%f ", x[i]);
//printf("\n");
//printf("gradient: ");
//for( ipc_ i = 0; i < n; i++) printf("%f ", g[i]);
//printf("\n");
// Delete internal workspace
lpa_terminate( &data, &control, &inform );
}
}
This is the same example, but now fortran-style indexing is used; the code is available in $GALAHAD/src/lpa/C/lpatf.c .
/* lpatf.c */
/* Full test for the LPA C interface using Fortran sparse matrix indexing */
#include <stdio.h>
#include <math.h>
#include "galahad_precision.h"
#include "galahad_cfunctions.h"
#include "galahad_lpa.h"
int main(void) {
// Derived types
void *data;
struct lpa_control_type control;
struct lpa_inform_type inform;
// Set problem data
ipc_ n = 3; // dimension
ipc_ m = 2; // number of general constraints
rpc_ g[] = {0.0, 2.0, 0.0}; // linear term in the objective
rpc_ f = 1.0; // constant term in the objective
ipc_ A_ne = 4; // Jacobian elements
ipc_ A_row[] = {1, 1, 2, 2}; // row indices
ipc_ A_col[] = {1, 2, 2, 3}; // column indices
ipc_ A_ptr[] = {1, 3, 5}; // row pointers
rpc_ A_val[] = {2.0, 1.0, 1.0, 1.0 }; // values
rpc_ c_l[] = {1.0, 2.0}; // constraint lower bound
rpc_ c_u[] = {2.0, 2.0}; // constraint upper bound
rpc_ x_l[] = {-1.0, - INFINITY, - INFINITY}; // variable lower bound
rpc_ x_u[] = {1.0, INFINITY, 2.0}; // variable upper bound
// Set output storage
rpc_ c[m]; // constraint values
ipc_ x_stat[n]; // variable status
ipc_ c_stat[m]; // constraint status
char st = ' ';
ipc_ status;
printf(" Fortran sparse matrix indexing\n\n");
printf(" basic tests of lp storage formats\n\n");
for( ipc_ d=1; d <= 3; d++){
// Initialize LPA
lpa_initialize( &data, &control, &status );
// Set user-defined control options
control.f_indexing = true; // Fortran sparse matrix indexing
// Start from 0
rpc_ x[] = {0.0,0.0,0.0};
rpc_ y[] = {0.0,0.0};
rpc_ z[] = {0.0,0.0,0.0};
switch(d){
case 1: // sparse co-ordinate storage
st = 'C';
lpa_import( &control, &data, &status, n, m,
"coordinate", A_ne, A_row, A_col, NULL );
lpa_solve_lp( &data, &status, n, m, g, f,
A_ne, A_val, c_l, c_u, x_l, x_u, x, c, y, z,
x_stat, c_stat );
break;
printf(" case %1" i_ipc_ " break\n",d);
case 2: // sparse by rows
st = 'R';
lpa_import( &control, &data, &status, n, m,
"sparse_by_rows", A_ne, NULL, A_col, A_ptr );
lpa_solve_lp( &data, &status, n, m, g, f,
A_ne, A_val, c_l, c_u, x_l, x_u, x, c, y, z,
x_stat, c_stat );
break;
case 3: // dense
st = 'D';
ipc_ A_dense_ne = 6; // number of elements of A
rpc_ A_dense[] = {2.0, 1.0, 0.0, 0.0, 1.0, 1.0};
lpa_import( &control, &data, &status, n, m,
"dense", A_ne, NULL, NULL, NULL );
lpa_solve_lp( &data, &status, n, m, g, f,
A_dense_ne, A_dense, c_l, c_u, x_l, x_u,
x, c, y, z, x_stat, c_stat );
break;
}
lpa_information( &data, &inform, &status );
if(inform.status == 0){
printf("%c:%6" i_ipc_ " iterations. Optimal objective value = %5.2f status = %1" i_ipc_ "\n",
st, inform.iter, inform.obj, inform.status);
}else{
printf("%c: LPA_solve exit status = %1" i_ipc_ "\n", st, inform.status);
}
//printf("x: ");
//for( ipc_ i = 0; i < n; i++) printf("%f ", x[i]);
//printf("\n");
//printf("gradient: ");
//for( ipc_ i = 0; i < n; i++) printf("%f ", g[i]);
//printf("\n");
// Delete internal workspace
lpa_terminate( &data, &control, &inform );
}
}