alpar@1
|
1 |
%* glpk04.tex *%
|
alpar@1
|
2 |
|
alpar@1
|
3 |
\chapter{Advanced API Routines}
|
alpar@1
|
4 |
|
alpar@1
|
5 |
\section{Background}
|
alpar@1
|
6 |
\label{basbgd}
|
alpar@1
|
7 |
|
alpar@1
|
8 |
Using vector and matrix notations LP problem (1.1)---(1.3) (see Section
|
alpar@1
|
9 |
\ref{seclp}, page \pageref{seclp}) can be stated as follows:
|
alpar@1
|
10 |
|
alpar@1
|
11 |
\medskip
|
alpar@1
|
12 |
|
alpar@1
|
13 |
\noindent
|
alpar@1
|
14 |
\hspace{.5in} minimize (or maximize)
|
alpar@1
|
15 |
$$z=c^Tx_S+c_0\eqno(3.1)$$
|
alpar@1
|
16 |
\hspace{.5in} subject to linear constraints
|
alpar@1
|
17 |
$$x_R=Ax_S\eqno(3.2)$$
|
alpar@1
|
18 |
\hspace{.5in} and bounds of variables
|
alpar@1
|
19 |
$$
|
alpar@1
|
20 |
\begin{array}{l@{\ }c@{\ }l@{\ }c@{\ }l}
|
alpar@1
|
21 |
l_R&\leq&x_R&\leq&u_R\\
|
alpar@1
|
22 |
l_S&\leq&x_S&\leq&u_S\\
|
alpar@1
|
23 |
\end{array}\eqno(3.3)
|
alpar@1
|
24 |
$$
|
alpar@1
|
25 |
where:
|
alpar@1
|
26 |
|
alpar@1
|
27 |
\noindent
|
alpar@1
|
28 |
$x_R=(x_1,\dots,x_m)$ is the vector of auxiliary variables;
|
alpar@1
|
29 |
|
alpar@1
|
30 |
\noindent
|
alpar@1
|
31 |
$x_S=(x_{m+1},\dots,x_{m+n})$ is the vector of structural
|
alpar@1
|
32 |
variables;
|
alpar@1
|
33 |
|
alpar@1
|
34 |
\noindent
|
alpar@1
|
35 |
$z$ is the objective function;
|
alpar@1
|
36 |
|
alpar@1
|
37 |
\noindent
|
alpar@1
|
38 |
$c=(c_1,\dots,c_n)$ is the vector of objective coefficients;
|
alpar@1
|
39 |
|
alpar@1
|
40 |
\noindent
|
alpar@1
|
41 |
$c_0$ is the constant term (``shift'') of the objective function;
|
alpar@1
|
42 |
|
alpar@1
|
43 |
\noindent
|
alpar@1
|
44 |
$A=(a_{11},\dots,a_{mn})$ is the constraint matrix;
|
alpar@1
|
45 |
|
alpar@1
|
46 |
\noindent
|
alpar@1
|
47 |
$l_R=(l_1,\dots,l_m)$ is the vector of lower bounds of auxiliary
|
alpar@1
|
48 |
variables;
|
alpar@1
|
49 |
|
alpar@1
|
50 |
\noindent
|
alpar@1
|
51 |
$u_R=(u_1,\dots,u_m)$ is the vector of upper bounds of auxiliary
|
alpar@1
|
52 |
variables;
|
alpar@1
|
53 |
|
alpar@1
|
54 |
\noindent
|
alpar@1
|
55 |
$l_S=(l_{m+1},\dots,l_{m+n})$ is the vector of lower bounds of
|
alpar@1
|
56 |
structural variables;
|
alpar@1
|
57 |
|
alpar@1
|
58 |
\noindent
|
alpar@1
|
59 |
$u_S=(u_{m+1},\dots,u_{m+n})$ is the vector of upper bounds of
|
alpar@1
|
60 |
structural variables.
|
alpar@1
|
61 |
|
alpar@1
|
62 |
\medskip
|
alpar@1
|
63 |
|
alpar@1
|
64 |
From the simplex method's standpoint there is no difference between
|
alpar@1
|
65 |
auxiliary and structural variables. This allows combining all these
|
alpar@1
|
66 |
variables into one vector that leads to the following problem statement:
|
alpar@1
|
67 |
|
alpar@1
|
68 |
\medskip
|
alpar@1
|
69 |
|
alpar@1
|
70 |
\noindent
|
alpar@1
|
71 |
\hspace{.5in} minimize (or maximize)
|
alpar@1
|
72 |
$$z=(0\ |\ c)^Tx+c_0\eqno(3.4)$$
|
alpar@1
|
73 |
\hspace{.5in} subject to linear constraints
|
alpar@1
|
74 |
$$(I\ |-\!A)x=0\eqno(3.5)$$
|
alpar@1
|
75 |
\hspace{.5in} and bounds of variables
|
alpar@1
|
76 |
$$l\leq x\leq u\eqno(3.6)$$
|
alpar@1
|
77 |
where:
|
alpar@1
|
78 |
|
alpar@1
|
79 |
\noindent
|
alpar@1
|
80 |
$x=(x_R\ |\ x_S)$ is the $(m+n)$-vector of (all) variables;
|
alpar@1
|
81 |
|
alpar@1
|
82 |
\noindent
|
alpar@1
|
83 |
$(0\ |\ c)$ is the $(m+n)$-vector of objective
|
alpar@1
|
84 |
coefficients;\footnote{Subvector 0 corresponds to objective coefficients
|
alpar@1
|
85 |
at auxiliary variables.}
|
alpar@1
|
86 |
|
alpar@1
|
87 |
\noindent
|
alpar@1
|
88 |
$(I\ |-\!A)$ is the {\it augmented} constraint
|
alpar@1
|
89 |
$m\times(m+n)$-matrix;\footnote{Note that due to auxiliary variables
|
alpar@1
|
90 |
matrix $(I\ |-\!A)$ contains the unity submatrix and therefore has full
|
alpar@1
|
91 |
rank. This means, in particular, that the system (3.5) has no linearly
|
alpar@1
|
92 |
dependent constraints.}
|
alpar@1
|
93 |
|
alpar@1
|
94 |
\noindent
|
alpar@1
|
95 |
$l=(l_R\ |\ l_S)$ is the $(m+n)$-vector of lower bounds of (all)
|
alpar@1
|
96 |
variables;
|
alpar@1
|
97 |
|
alpar@1
|
98 |
\noindent
|
alpar@1
|
99 |
$u=(u_R\ |\ u_S)$ is the $(m+n)$-vector of upper bounds of (all)
|
alpar@1
|
100 |
variables.
|
alpar@1
|
101 |
|
alpar@1
|
102 |
\medskip
|
alpar@1
|
103 |
|
alpar@1
|
104 |
By definition an {\it LP basic solution} geometrically is a point in
|
alpar@1
|
105 |
the space of all variables, which is the intersection of planes
|
alpar@1
|
106 |
corresponding to active constraints\footnote{A constraint is called
|
alpar@1
|
107 |
{\it active} if in a given point it is satisfied as equality, otherwise
|
alpar@1
|
108 |
it is called {\it inactive}.}. The space of all variables has the
|
alpar@1
|
109 |
dimension $m+n$, therefore, to define some basic solution we have to
|
alpar@1
|
110 |
define $m+n$ active constraints. Note that $m$ constraints (3.5) being
|
alpar@1
|
111 |
linearly independent equalities are always active, so remaining $n$
|
alpar@1
|
112 |
active constraints can be chosen only from bound constraints (3.6).
|
alpar@1
|
113 |
|
alpar@1
|
114 |
A variable is called {\it non-basic}, if its (lower or upper) bound is
|
alpar@1
|
115 |
active, otherwise it is called {\it basic}. Since, as was said above,
|
alpar@1
|
116 |
exactly $n$ bound constraints must be active, in any basic solution
|
alpar@1
|
117 |
there are always $n$ non-basic variables and $m$ basic variables.
|
alpar@1
|
118 |
(Note that a free variable also can be non-basic. Although such
|
alpar@1
|
119 |
variable has no bounds, we can think it as the difference between two
|
alpar@1
|
120 |
non-negative variables, which both are non-basic in this case.)
|
alpar@1
|
121 |
|
alpar@1
|
122 |
Now consider how to determine numeric values of all variables for a
|
alpar@1
|
123 |
given basic solution.
|
alpar@1
|
124 |
|
alpar@1
|
125 |
Let $\Pi$ be an appropriate permutation matrix of the order $(m+n)$.
|
alpar@1
|
126 |
Then we can write:
|
alpar@1
|
127 |
$$\left(\begin{array}{@{}c@{}}x_B\\x_N\\\end{array}\right)=
|
alpar@1
|
128 |
\Pi\left(\begin{array}{@{}c@{}}x_R\\x_S\\\end{array}\right)=\Pi x,
|
alpar@1
|
129 |
\eqno(3.7)$$
|
alpar@1
|
130 |
where $x_B$ is the vector of basic variables, $x_N$ is the vector of
|
alpar@1
|
131 |
non-basic variables, $x=(x_R\ |\ x_S)$ is the vector of all variables
|
alpar@1
|
132 |
in the original order. In this case the system of linear constraints
|
alpar@1
|
133 |
(3.5) can be rewritten as follows:
|
alpar@1
|
134 |
$$(I\ |-\!A)\Pi^T\Pi x=0\ \ \ \Rightarrow\ \ \ (B\ |\ N)
|
alpar@1
|
135 |
\left(\begin{array}{@{}c@{}}x_B\\x_N\\\end{array}\right)=0,\eqno(3.8)$$
|
alpar@1
|
136 |
where
|
alpar@1
|
137 |
$$(B\ |\ N)=(I\ |-\!A)\Pi^T.\eqno(3.9)$$
|
alpar@1
|
138 |
Matrix $B$ is a square non-singular $m\times m$-matrix, which is
|
alpar@1
|
139 |
composed from columns of the augmented constraint matrix corresponding
|
alpar@1
|
140 |
to basic variables. It is called the {\it basis matrix} or simply the
|
alpar@1
|
141 |
{\it basis}. Matrix $N$ is a rectangular $m\times n$-matrix, which is
|
alpar@1
|
142 |
composed from columns of the augmented constraint matrix corresponding
|
alpar@1
|
143 |
to non-basic variables.
|
alpar@1
|
144 |
|
alpar@1
|
145 |
From (3.8) it follows that:
|
alpar@1
|
146 |
$$Bx_B+Nx_N=0,\eqno(3.10)$$
|
alpar@1
|
147 |
therefore,
|
alpar@1
|
148 |
$$x_B=-B^{-1}Nx_N.\eqno(3.11)$$
|
alpar@1
|
149 |
Thus, the formula (3.11) shows how to determine numeric values of basic
|
alpar@1
|
150 |
variables $x_B$ assuming that non-basic variables $x_N$ are fixed on
|
alpar@1
|
151 |
their active bounds.
|
alpar@1
|
152 |
|
alpar@1
|
153 |
The $m\times n$-matrix
|
alpar@1
|
154 |
$$\Xi=-B^{-1}N,\eqno(3.12)$$
|
alpar@1
|
155 |
which appears in (3.11), is called the {\it simplex
|
alpar@1
|
156 |
tableau}.\footnote{This definition corresponds to the GLPK
|
alpar@1
|
157 |
implementation.} It shows how basic variables depend on non-basic
|
alpar@1
|
158 |
variables:
|
alpar@1
|
159 |
$$x_B=\Xi x_N.\eqno(3.13)$$
|
alpar@1
|
160 |
|
alpar@1
|
161 |
The system (3.13) is equivalent to the system (3.5) in the sense that
|
alpar@1
|
162 |
they both define the same set of points in the space of (primal)
|
alpar@1
|
163 |
variables, which satisfy to these systems. If, moreover, values of all
|
alpar@1
|
164 |
basic variables satisfy to their bound constraints (3.3), the
|
alpar@1
|
165 |
corresponding basic solution is called {\it (primal) feasible},
|
alpar@1
|
166 |
otherwise {\it (primal) infeasible}. It is understood that any (primal)
|
alpar@1
|
167 |
feasible basic solution satisfy to all constraints (3.2) and (3.3).
|
alpar@1
|
168 |
|
alpar@1
|
169 |
The LP theory says that if LP has optimal solution, it has (at least
|
alpar@1
|
170 |
one) basic feasible solution, which corresponds to the optimum. And the
|
alpar@1
|
171 |
most natural way to determine whether a given basic solution is optimal
|
alpar@1
|
172 |
or not is to use the Karush---Kuhn---Tucker optimality conditions.
|
alpar@1
|
173 |
|
alpar@1
|
174 |
\def\arraystretch{1.5}
|
alpar@1
|
175 |
|
alpar@1
|
176 |
For the problem statement (3.4)---(3.6) the optimality conditions are
|
alpar@1
|
177 |
the following:\footnote{These conditions can be appiled to any solution,
|
alpar@1
|
178 |
not only to a basic solution.}
|
alpar@1
|
179 |
$$(I\ |-\!A)x=0\eqno(3.14)$$
|
alpar@1
|
180 |
$$(I\ |-\!A)^T\pi+\lambda_l+\lambda_u=\nabla z=(0\ |\ c)^T\eqno(3.15)$$
|
alpar@1
|
181 |
$$l\leq x\leq u\eqno(3.16)$$
|
alpar@1
|
182 |
$$\lambda_l\geq 0,\ \ \lambda_u\leq 0\ \ \mbox{(minimization)}
|
alpar@1
|
183 |
\eqno(3.17)$$
|
alpar@1
|
184 |
$$\lambda_l\leq 0,\ \ \lambda_u\geq 0\ \ \mbox{(maximization)}
|
alpar@1
|
185 |
\eqno(3.18)$$
|
alpar@1
|
186 |
$$(\lambda_l)_k(x_k-l_k)=0,\ \ (\lambda_u)_k(x_k-u_k)=0,\ \ k=1,2,\dots,
|
alpar@1
|
187 |
m+n\eqno(3.19)$$
|
alpar@1
|
188 |
where:
|
alpar@1
|
189 |
$\pi=(\pi_1,\pi_2,\dots,\pi_m)$ is a $m$-vector of Lagrange
|
alpar@1
|
190 |
multipliers for equality constraints (3.5);
|
alpar@1
|
191 |
$\lambda_l=[(\lambda_l)_1,(\lambda_l)_2,\dots,(\lambda_l)_n]$ is a
|
alpar@1
|
192 |
$n$-vector of Lagrange multipliers for lower bound constraints (3.6);
|
alpar@1
|
193 |
$\lambda_u=[(\lambda_u)_1,(\lambda_u)_2,\dots,(\lambda_u)_n]$ is a
|
alpar@1
|
194 |
$n$-vector of Lagrange multipliers for upper bound constraints (3.6).
|
alpar@1
|
195 |
|
alpar@1
|
196 |
Condition (3.14) is the {\it primal} (original) system of equality
|
alpar@1
|
197 |
constraints (3.5).
|
alpar@1
|
198 |
|
alpar@1
|
199 |
Condition (3.15) is the {\it dual} system of equality constraints.
|
alpar@1
|
200 |
It requires the gradient of the objective function to be a linear
|
alpar@1
|
201 |
combination of normals to the planes defined by constraints of the
|
alpar@1
|
202 |
original problem.
|
alpar@1
|
203 |
|
alpar@1
|
204 |
Condition (3.16) is the primal (original) system of bound constraints
|
alpar@1
|
205 |
(3.6).
|
alpar@1
|
206 |
|
alpar@1
|
207 |
Condition (3.17) (or (3.18) in case of maximization) is the dual system
|
alpar@1
|
208 |
of bound constraints.
|
alpar@1
|
209 |
|
alpar@1
|
210 |
Condition (3.19) is the {\it complementary slackness condition}. It
|
alpar@1
|
211 |
requires, for each original (auxiliary or structural) variable $x_k$,
|
alpar@1
|
212 |
that either its (lower or upper) bound must be active, or zero bound of
|
alpar@1
|
213 |
the corresponding Lagrange multiplier ($(\lambda_l)_k$ or
|
alpar@1
|
214 |
$(\lambda_u)_k$) must be active.
|
alpar@1
|
215 |
|
alpar@1
|
216 |
In GLPK two multipliers $(\lambda_l)_k$ and $(\lambda_u)_k$ for each
|
alpar@1
|
217 |
primal (original) variable $x_k$, $k=1,2,\dots,m+n$, are combined into
|
alpar@1
|
218 |
one multiplier:
|
alpar@1
|
219 |
$$\lambda_k=(\lambda_l)_k+(\lambda_u)_k,\eqno(3.20)$$
|
alpar@1
|
220 |
which is called a {\it dual variable} for $x_k$. This {\it cannot} lead
|
alpar@1
|
221 |
to the ambiguity, because both lower and upper bounds of $x_k$ cannot be
|
alpar@1
|
222 |
active at the same time,\footnote{If $x_k$ is a fixed variable, we can
|
alpar@1
|
223 |
think it as double-bounded variable $l_k\leq x_k\leq u_k$, where
|
alpar@1
|
224 |
$l_k=u_k.$} so at least one of $(\lambda_l)_k$ and $(\lambda_u)_k$ must
|
alpar@1
|
225 |
be equal to zero, and because these multipliers have different signs,
|
alpar@1
|
226 |
the combined multiplier, which is their sum, uniquely defines each of
|
alpar@1
|
227 |
them.
|
alpar@1
|
228 |
|
alpar@1
|
229 |
\def\arraystretch{1}
|
alpar@1
|
230 |
|
alpar@1
|
231 |
Using dual variables $\lambda_k$ the dual system of bound constraints
|
alpar@1
|
232 |
(3.17) and (3.18) can be written in the form of so called {\it ``rule of
|
alpar@1
|
233 |
signs''} as follows:
|
alpar@1
|
234 |
|
alpar@1
|
235 |
\begin{center}
|
alpar@1
|
236 |
\begin{tabular}{|@{\,}c@{$\,$}|@{$\,$}c@{$\,$}|@{$\,$}c@{$\,$}|
|
alpar@1
|
237 |
@{$\,$}c|c@{$\,$}|@{$\,$}c@{$\,$}|@{$\,$}c@{$\,$}|}
|
alpar@1
|
238 |
\hline
|
alpar@1
|
239 |
Original bound&\multicolumn{3}{c|}{Minimization}&\multicolumn{3}{c|}
|
alpar@1
|
240 |
{Maximization}\\
|
alpar@1
|
241 |
\cline{2-7}
|
alpar@1
|
242 |
constraint&$(\lambda_l)_k$&$(\lambda_u)_k$&$(\lambda_l)_k+
|
alpar@1
|
243 |
(\lambda_u)_k$&$(\lambda_l)_k$&$(\lambda_u)_k$&$(\lambda_l)_k+
|
alpar@1
|
244 |
(\lambda_u)_k$\\
|
alpar@1
|
245 |
\hline
|
alpar@1
|
246 |
$-\infty<x_k<+\infty$&$=0$&$=0$&$\lambda_k=0$&$=0$&$=0$&$\lambda_k=0$\\
|
alpar@1
|
247 |
$x_k\geq l_k$&$\geq 0$&$=0$&$\lambda_k\geq 0$&$\leq 0$&$=0$&$\lambda_k
|
alpar@1
|
248 |
\leq0$\\
|
alpar@1
|
249 |
$x_k\leq u_k$&$=0$&$\leq 0$&$\lambda_k\leq 0$&$=0$&$\geq 0$&$\lambda_k
|
alpar@1
|
250 |
\geq0$\\
|
alpar@1
|
251 |
$l_k\leq x_k\leq u_k$&$\geq 0$& $\leq 0$& $-\infty\!<\!\lambda_k\!<
|
alpar@1
|
252 |
\!+\infty$
|
alpar@1
|
253 |
&$\leq 0$& $\geq 0$& $-\infty\!<\!\lambda_k\!<\!+\infty$\\
|
alpar@1
|
254 |
$x_k=l_k=u_k$&$\geq 0$& $\leq 0$& $-\infty\!<\!\lambda_k\!<\!+\infty$&
|
alpar@1
|
255 |
$\leq 0$&
|
alpar@1
|
256 |
$\geq 0$& $-\infty\!<\!\lambda_k\!<\!+\infty$\\
|
alpar@1
|
257 |
\hline
|
alpar@1
|
258 |
\end{tabular}
|
alpar@1
|
259 |
\end{center}
|
alpar@1
|
260 |
|
alpar@1
|
261 |
May note that each primal variable $x_k$ has its dual counterpart
|
alpar@1
|
262 |
$\lambda_k$ and vice versa. This allows applying the same partition for
|
alpar@1
|
263 |
the vector of dual variables as (3.7):
|
alpar@1
|
264 |
$$\left(\begin{array}{@{}c@{}}\lambda_B\\\lambda_N\\\end{array}\right)=
|
alpar@1
|
265 |
\Pi\lambda,\eqno(3.21)$$
|
alpar@1
|
266 |
where $\lambda_B$ is a vector of dual variables for basic variables
|
alpar@1
|
267 |
$x_B$, $\lambda_N$ is a vector of dual variables for non-basic variables
|
alpar@1
|
268 |
$x_N$.
|
alpar@1
|
269 |
|
alpar@1
|
270 |
By definition, bounds of basic variables are inactive constraints, so in
|
alpar@1
|
271 |
any basic solution $\lambda_B=0$. Corresponding values of dual variables
|
alpar@1
|
272 |
$\lambda_N$ for non-basic variables $x_N$ can be determined in the
|
alpar@1
|
273 |
following way. From the dual system (3.15) we have:
|
alpar@1
|
274 |
$$(I\ |-\!A)^T\pi+\lambda=(0\ |\ c)^T,\eqno(3.22)$$
|
alpar@1
|
275 |
so multiplying both sides of (3.22) by matrix $\Pi$ gives:
|
alpar@1
|
276 |
$$\Pi(I\ |-\!A)^T\pi+\Pi\lambda=\Pi(0\ |\ c)^T.\eqno(3.23)$$
|
alpar@1
|
277 |
From (3.9) it follows that
|
alpar@1
|
278 |
$$\Pi(I\ |-\!A)^T=[(I\ |-\!A)\Pi^T]^T=(B\ |\ N)^T.\eqno(3.24)$$
|
alpar@1
|
279 |
Further, we can apply the partition (3.7) also to the vector of
|
alpar@1
|
280 |
objective coefficients (see (3.4)):
|
alpar@1
|
281 |
$$\left(\begin{array}{@{}c@{}}c_B\\c_N\\\end{array}\right)=
|
alpar@1
|
282 |
\Pi\left(\begin{array}{@{}c@{}}0\\c\\\end{array}\right),\eqno(3.25)$$
|
alpar@1
|
283 |
where $c_B$ is a vector of objective coefficients at basic variables,
|
alpar@1
|
284 |
$c_N$ is a vector of objective coefficients at non-basic variables.
|
alpar@1
|
285 |
Now, substituting (3.24), (3.21), and (3.25) into (3.23), leads to:
|
alpar@1
|
286 |
$$(B\ |\ N)^T\pi+(\lambda_B\ |\ \lambda_N)^T=(c_B\ |\ c_N)^T,
|
alpar@1
|
287 |
\eqno(3.26)$$
|
alpar@1
|
288 |
and transposing both sides of (3.26) gives the system:
|
alpar@1
|
289 |
$$\left(\begin{array}{@{}c@{}}B^T\\N^T\\\end{array}\right)\pi+
|
alpar@1
|
290 |
\left(\begin{array}{@{}c@{}}\lambda_B\\\lambda_N\\\end{array}\right)=
|
alpar@1
|
291 |
\left(\begin{array}{@{}c@{}}c_B\\c_T\\\end{array}\right),\eqno(3.27)$$
|
alpar@1
|
292 |
which can be written as follows:
|
alpar@1
|
293 |
$$\left\{
|
alpar@1
|
294 |
\begin{array}{@{\ }r@{\ }c@{\ }r@{\ }c@{\ }l@{\ }}
|
alpar@1
|
295 |
B^T\pi&+&\lambda_B&=&c_B\\
|
alpar@1
|
296 |
N^T\pi&+&\lambda_N&=&c_N\\
|
alpar@1
|
297 |
\end{array}
|
alpar@1
|
298 |
\right.\eqno(3.28)
|
alpar@1
|
299 |
$$
|
alpar@1
|
300 |
Lagrange multipliers $\pi=(\pi_i)$ correspond to equality constraints
|
alpar@1
|
301 |
(3.5) and therefore can have any sign. This allows resolving the first
|
alpar@1
|
302 |
subsystem of (3.28) as follows:\footnote{$B^{-T}$ means $(B^T)^{-1}=
|
alpar@1
|
303 |
(B^{-1})^T$.}
|
alpar@1
|
304 |
$$\pi=B^{-T}(c_B-\lambda_B)=-B^{-T}\lambda_B+B^{-T}c_B,\eqno(3.29)$$
|
alpar@1
|
305 |
and substitution of $\pi$ from (3.29) into the second subsystem of
|
alpar@1
|
306 |
(3.28) gives:
|
alpar@1
|
307 |
$$\lambda_N=-N^T\pi+c_N=N^TB^{-T}\lambda_B+(c_N-N^TB^{-T}c_B).
|
alpar@1
|
308 |
\eqno(3.30)$$
|
alpar@1
|
309 |
The latter system can be written in the following final form:
|
alpar@1
|
310 |
$$\lambda_N=-\Xi^T\lambda_B+d,\eqno(3.31)$$
|
alpar@1
|
311 |
where $\Xi$ is the simplex tableau (see (3.12)), and
|
alpar@1
|
312 |
$$d=c_N-N^TB^{-T}c_B=c_N+\Xi^Tc_B\eqno(3.32)$$
|
alpar@1
|
313 |
is the vector of so called {\it reduced costs} of non-basic variables.
|
alpar@1
|
314 |
|
alpar@1
|
315 |
\pagebreak
|
alpar@1
|
316 |
|
alpar@1
|
317 |
Above it was said that in any basic solution $\lambda_B=0$, so
|
alpar@1
|
318 |
$\lambda_N=d$ as it follows from (3.31).
|
alpar@1
|
319 |
|
alpar@1
|
320 |
The system (3.31) is equivalent to the system (3.15) in the sense that
|
alpar@1
|
321 |
they both define the same set of points in the space of dual variables
|
alpar@1
|
322 |
$\lambda$, which satisfy to these systems. If, moreover, values of all
|
alpar@1
|
323 |
dual variables $\lambda_N$ (i.e. reduced costs $d$) satisfy to their
|
alpar@1
|
324 |
bound constraints (i.e. to the ``rule of signs''; see the table above),
|
alpar@1
|
325 |
the corresponding basic solution is called {\it dual feasible},
|
alpar@1
|
326 |
otherwise {\it dual infeasible}. It is understood that any dual feasible
|
alpar@1
|
327 |
solution satisfy to all constraints (3.15) and (3.17) (or (3.18) in case
|
alpar@1
|
328 |
of maximization).
|
alpar@1
|
329 |
|
alpar@1
|
330 |
It can be easily shown that the complementary slackness condition
|
alpar@1
|
331 |
(3.19) is always satisfied for {\it any} basic solution. Therefore,
|
alpar@1
|
332 |
a basic solution\footnote{It is assumed that a complete basic solution
|
alpar@1
|
333 |
has the form $(x,\lambda)$, i.e. it includes primal as well as dual
|
alpar@1
|
334 |
variables.} is {\it optimal} if and only if it is primal and dual
|
alpar@1
|
335 |
feasible, because in this case it satifies to all the optimality
|
alpar@1
|
336 |
conditions (3.14)---(3.19).
|
alpar@1
|
337 |
|
alpar@1
|
338 |
\def\arraystretch{1.5}
|
alpar@1
|
339 |
|
alpar@1
|
340 |
The meaning of reduced costs $d=(d_j)$ of non-basic variables can be
|
alpar@1
|
341 |
explained in the following way. From (3.4), (3.7), and (3.25) it follows
|
alpar@1
|
342 |
that:
|
alpar@1
|
343 |
$$z=c_B^Tx_B+c_N^Tx_N+c_0.\eqno(3.33)$$
|
alpar@1
|
344 |
Substituting $x_B$ from (3.11) into (3.33) we can eliminate basic
|
alpar@1
|
345 |
variables and express the objective only through non-basic variables:
|
alpar@1
|
346 |
$$
|
alpar@1
|
347 |
\begin{array}{r@{\ }c@{\ }l}
|
alpar@1
|
348 |
z&=&c_B^T(-B^{-1}Nx_N)+c_N^Tx_N+c_0=\\
|
alpar@1
|
349 |
&=&(c_N^T-c_B^TB^{-1}N)x_N+c_0=\\
|
alpar@1
|
350 |
&=&(c_N-N^TB^{-T}c_B)^Tx_N+c_0=\\
|
alpar@1
|
351 |
&=&d^Tx_N+c_0.
|
alpar@1
|
352 |
\end{array}\eqno(3.34)
|
alpar@1
|
353 |
$$
|
alpar@1
|
354 |
From (3.34) it is seen that reduced cost $d_j$ shows how the objective
|
alpar@1
|
355 |
function $z$ depends on non-basic variable $(x_N)_j$ in the neighborhood
|
alpar@1
|
356 |
of the current basic solution, i.e. while the current basis remains
|
alpar@1
|
357 |
unchanged.
|
alpar@1
|
358 |
|
alpar@1
|
359 |
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
|
alpar@1
|
360 |
|
alpar@1
|
361 |
\newpage
|
alpar@1
|
362 |
|
alpar@1
|
363 |
\section{LP basis routines}
|
alpar@1
|
364 |
\label{lpbasis}
|
alpar@1
|
365 |
|
alpar@1
|
366 |
\subsection{glp\_bf\_exists---check if the basis factorization exists}
|
alpar@1
|
367 |
|
alpar@1
|
368 |
\subsubsection*{Synopsis}
|
alpar@1
|
369 |
|
alpar@1
|
370 |
\begin{verbatim}
|
alpar@1
|
371 |
int glp_bf_exists(glp_prob *lp);
|
alpar@1
|
372 |
\end{verbatim}
|
alpar@1
|
373 |
|
alpar@1
|
374 |
\subsubsection*{Returns}
|
alpar@1
|
375 |
|
alpar@1
|
376 |
If the basis factorization for the current basis associated with the
|
alpar@1
|
377 |
specified problem object exists and therefore is available for
|
alpar@1
|
378 |
computations, the routine \verb|glp_bf_exists| returns non-zero.
|
alpar@1
|
379 |
Otherwise the routine returns zero.
|
alpar@1
|
380 |
|
alpar@1
|
381 |
\subsubsection*{Comments}
|
alpar@1
|
382 |
|
alpar@1
|
383 |
Let the problem object have $m$ rows and $n$ columns. In GLPK the
|
alpar@1
|
384 |
{\it basis matrix} $B$ is a square non-singular matrix of the order $m$,
|
alpar@1
|
385 |
whose columns correspond to basic (auxiliary and/or structural)
|
alpar@1
|
386 |
variables. It is defined by the following main
|
alpar@1
|
387 |
equality:\footnote{For more details see Subsection \ref{basbgd},
|
alpar@1
|
388 |
page \pageref{basbgd}.}
|
alpar@1
|
389 |
$$(B\ |\ N)=(I\ |-\!A)\Pi^T,$$
|
alpar@1
|
390 |
where $I$ is the unity matrix of the order $m$, whose columns correspond
|
alpar@1
|
391 |
to auxiliary variables; $A$ is the original constraint
|
alpar@1
|
392 |
$m\times n$-matrix, whose columns correspond to structural variables;
|
alpar@1
|
393 |
$(I\ |-\!A)$ is the augmented constraint\linebreak
|
alpar@1
|
394 |
$m\times(m+n)$-matrix, whose columns correspond to all (auxiliary and
|
alpar@1
|
395 |
structural) variables following in the original order; $\Pi$ is a
|
alpar@1
|
396 |
permutation matrix of the order $m+n$; and $N$ is a rectangular
|
alpar@1
|
397 |
$m\times n$-matrix, whose columns correspond to non-basic (auxiliary
|
alpar@1
|
398 |
and/or structural) variables.
|
alpar@1
|
399 |
|
alpar@1
|
400 |
For various reasons it may be necessary to solve linear systems with
|
alpar@1
|
401 |
matrix $B$. To provide this possibility the GLPK implementation
|
alpar@1
|
402 |
maintains an invertable form of $B$ (that is, some representation of
|
alpar@1
|
403 |
$B^{-1}$) called the {\it basis factorization}, which is an internal
|
alpar@1
|
404 |
component of the problem object. Typically, the basis factorization is
|
alpar@1
|
405 |
computed by the simplex solver, which keeps it in the problem object
|
alpar@1
|
406 |
to be available for other computations.
|
alpar@1
|
407 |
|
alpar@1
|
408 |
Should note that any changes in the problem object, which affects the
|
alpar@1
|
409 |
basis matrix (e.g. changing the status of a row or column, changing
|
alpar@1
|
410 |
a basic column of the constraint matrix, removing an active constraint,
|
alpar@1
|
411 |
etc.), invalidates the basis factorization. So before calling any API
|
alpar@1
|
412 |
routine, which uses the basis factorization, the application program
|
alpar@1
|
413 |
must make sure (using the routine \verb|glp_bf_exists|) that the
|
alpar@1
|
414 |
factorization exists and therefore available for computations.
|
alpar@1
|
415 |
|
alpar@1
|
416 |
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
|
alpar@1
|
417 |
|
alpar@1
|
418 |
\subsection{glp\_factorize---compute the basis factorization}
|
alpar@1
|
419 |
|
alpar@1
|
420 |
\subsubsection*{Synopsis}
|
alpar@1
|
421 |
|
alpar@1
|
422 |
\begin{verbatim}
|
alpar@1
|
423 |
int glp_factorize(glp_prob *lp);
|
alpar@1
|
424 |
\end{verbatim}
|
alpar@1
|
425 |
|
alpar@1
|
426 |
\subsubsection*{Description}
|
alpar@1
|
427 |
|
alpar@1
|
428 |
The routine \verb|glp_factorize| computes the basis factorization for
|
alpar@1
|
429 |
the current basis associated with the specified problem
|
alpar@1
|
430 |
object.\footnote{The current basis is defined by the current statuses
|
alpar@1
|
431 |
of rows (auxiliary variables) and columns (structural variables).}
|
alpar@1
|
432 |
|
alpar@1
|
433 |
The basis factorization is computed from ``scratch'' even if it exists,
|
alpar@1
|
434 |
so the application program may use the routine \verb|glp_bf_exists|,
|
alpar@1
|
435 |
and, if the basis factorization already exists, not to call the routine
|
alpar@1
|
436 |
\verb|glp_factorize| to prevent an extra work.
|
alpar@1
|
437 |
|
alpar@1
|
438 |
The routine \verb|glp_factorize| {\it does not} compute components of
|
alpar@1
|
439 |
the basic solution (i.e. primal and dual values).
|
alpar@1
|
440 |
|
alpar@1
|
441 |
\subsubsection*{Returns}
|
alpar@1
|
442 |
|
alpar@1
|
443 |
\begin{tabular}{@{}p{25mm}p{97.3mm}@{}}
|
alpar@1
|
444 |
0 & The basis factorization has been successfully computed.\\
|
alpar@1
|
445 |
\verb|GLP_EBADB| & The basis matrix is invalid, because the number of
|
alpar@1
|
446 |
basic (auxiliary and structural) variables is not the same as the number
|
alpar@1
|
447 |
of rows in the problem object.\\
|
alpar@1
|
448 |
\verb|GLP_ESING| & The basis matrix is singular within the working
|
alpar@1
|
449 |
precision.\\
|
alpar@1
|
450 |
\verb|GLP_ECOND| & The basis matrix is ill-conditioned, i.e. its
|
alpar@1
|
451 |
condition number is too large.\\
|
alpar@1
|
452 |
\end{tabular}
|
alpar@1
|
453 |
|
alpar@1
|
454 |
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
|
alpar@1
|
455 |
|
alpar@1
|
456 |
\newpage
|
alpar@1
|
457 |
|
alpar@1
|
458 |
\subsection{glp\_bf\_updated---check if the basis factorization has\\
|
alpar@1
|
459 |
been updated}
|
alpar@1
|
460 |
|
alpar@1
|
461 |
\subsubsection*{Synopsis}
|
alpar@1
|
462 |
|
alpar@1
|
463 |
\begin{verbatim}
|
alpar@1
|
464 |
int glp_bf_updated(glp_prob *lp);
|
alpar@1
|
465 |
\end{verbatim}
|
alpar@1
|
466 |
|
alpar@1
|
467 |
\subsubsection*{Returns}
|
alpar@1
|
468 |
|
alpar@1
|
469 |
If the basis factorization has been just computed from ``scratch'', the
|
alpar@1
|
470 |
routine \verb|glp_bf_updated| returns zero. Otherwise, if the
|
alpar@1
|
471 |
factorization has been updated at least once, the routine returns
|
alpar@1
|
472 |
non-zero.
|
alpar@1
|
473 |
|
alpar@1
|
474 |
\subsubsection*{Comments}
|
alpar@1
|
475 |
|
alpar@1
|
476 |
{\it Updating} the basis factorization means recomputing it to reflect
|
alpar@1
|
477 |
changes in the basis matrix. For example, on every iteration of the
|
alpar@1
|
478 |
simplex method some column of the current basis matrix is replaced by a
|
alpar@1
|
479 |
new column that gives a new basis matrix corresponding to the adjacent
|
alpar@1
|
480 |
basis. In this case computing the basis factorization for the adjacent
|
alpar@1
|
481 |
basis from ``scratch'' (as the routine \verb|glp_factorize| does) would
|
alpar@1
|
482 |
be too time-consuming.
|
alpar@1
|
483 |
|
alpar@1
|
484 |
On the other hand, since the basis factorization update is a numeric
|
alpar@1
|
485 |
computational procedure, applying it many times may lead to accumulating
|
alpar@1
|
486 |
round-off errors. Therefore the basis is periodically refactorized
|
alpar@1
|
487 |
(reinverted) from ``scratch'' (with the routine \verb|glp_factorize|)
|
alpar@1
|
488 |
that allows improving its numerical properties.
|
alpar@1
|
489 |
|
alpar@1
|
490 |
The routine \verb|glp_bf_updated| allows determining if the basis
|
alpar@1
|
491 |
factorization has been updated at least once since it was computed from
|
alpar@1
|
492 |
``scratch''.
|
alpar@1
|
493 |
|
alpar@1
|
494 |
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
|
alpar@1
|
495 |
|
alpar@1
|
496 |
\newpage
|
alpar@1
|
497 |
|
alpar@1
|
498 |
\subsection{glp\_get\_bfcp---retrieve basis factorization control
|
alpar@1
|
499 |
parameters}
|
alpar@1
|
500 |
|
alpar@1
|
501 |
\subsubsection*{Synopsis}
|
alpar@1
|
502 |
|
alpar@1
|
503 |
\begin{verbatim}
|
alpar@1
|
504 |
void glp_get_bfcp(glp_prob *lp, glp_bfcp *parm);
|
alpar@1
|
505 |
\end{verbatim}
|
alpar@1
|
506 |
|
alpar@1
|
507 |
\subsubsection*{Description}
|
alpar@1
|
508 |
|
alpar@1
|
509 |
The routine \verb|glp_get_bfcp| retrieves control parameters, which are
|
alpar@1
|
510 |
used on computing and updating the basis factorization associated with
|
alpar@1
|
511 |
the specified problem object.
|
alpar@1
|
512 |
|
alpar@1
|
513 |
Current values of the control parameters are stored in a \verb|glp_bfcp|
|
alpar@1
|
514 |
structure, which the parameter \verb|parm| points to. For a detailed
|
alpar@1
|
515 |
description of the structure \verb|glp_bfcp| see comments to the routine
|
alpar@1
|
516 |
\verb|glp_set_bfcp| in the next subsection.
|
alpar@1
|
517 |
|
alpar@1
|
518 |
\subsubsection*{Comments}
|
alpar@1
|
519 |
|
alpar@1
|
520 |
The purpose of the routine \verb|glp_get_bfcp| is two-fold. First, it
|
alpar@1
|
521 |
allows the application program obtaining current values of control
|
alpar@1
|
522 |
parameters used by internal GLPK routines, which compute and update the
|
alpar@1
|
523 |
basis factorization.
|
alpar@1
|
524 |
|
alpar@1
|
525 |
The second purpose of this routine is to provide proper values for all
|
alpar@1
|
526 |
fields of the structure \verb|glp_bfcp| in the case when the application
|
alpar@1
|
527 |
program needs to change some control parameters.
|
alpar@1
|
528 |
|
alpar@1
|
529 |
\subsection{glp\_set\_bfcp---change basis factorization control
|
alpar@1
|
530 |
parameters}
|
alpar@1
|
531 |
|
alpar@1
|
532 |
\subsubsection*{Synopsis}
|
alpar@1
|
533 |
|
alpar@1
|
534 |
\begin{verbatim}
|
alpar@1
|
535 |
void glp_set_bfcp(glp_prob *lp, const glp_bfcp *parm);
|
alpar@1
|
536 |
\end{verbatim}
|
alpar@1
|
537 |
|
alpar@1
|
538 |
\subsubsection*{Description}
|
alpar@1
|
539 |
|
alpar@1
|
540 |
The routine \verb|glp_set_bfcp| changes control parameters, which are
|
alpar@1
|
541 |
used by internal GLPK routines on computing and updating the basis
|
alpar@1
|
542 |
factorization associated with the specified problem object.
|
alpar@1
|
543 |
|
alpar@1
|
544 |
New values of the control parameters should be passed in a structure
|
alpar@1
|
545 |
\verb|glp_bfcp|, which the parameter \verb|parm| points to. For a
|
alpar@1
|
546 |
detailed description of the structure \verb|glp_bfcp| see paragraph
|
alpar@1
|
547 |
``Control parameters'' below.
|
alpar@1
|
548 |
|
alpar@1
|
549 |
The parameter \verb|parm| can be specified as \verb|NULL|, in which case
|
alpar@1
|
550 |
all control parameters are reset to their default values.
|
alpar@1
|
551 |
|
alpar@1
|
552 |
\subsubsection*{Comments}
|
alpar@1
|
553 |
|
alpar@1
|
554 |
Before changing some control parameters with the routine
|
alpar@1
|
555 |
\verb|glp_set_bfcp| the application program should retrieve current
|
alpar@1
|
556 |
values of all control parameters with the routine \verb|glp_get_bfcp|.
|
alpar@1
|
557 |
This is needed for backward compatibility, because in the future there
|
alpar@1
|
558 |
may appear new members in the structure \verb|glp_bfcp|.
|
alpar@1
|
559 |
|
alpar@1
|
560 |
Note that new values of control parameters come into effect on a next
|
alpar@1
|
561 |
computation of the basis factorization, not immediately.
|
alpar@1
|
562 |
|
alpar@1
|
563 |
\subsubsection*{Example}
|
alpar@1
|
564 |
|
alpar@1
|
565 |
\begin{verbatim}
|
alpar@1
|
566 |
glp_prob *lp;
|
alpar@1
|
567 |
glp_bfcp parm;
|
alpar@1
|
568 |
. . .
|
alpar@1
|
569 |
/* retrieve current values of control parameters */
|
alpar@1
|
570 |
glp_get_bfcp(lp, &parm);
|
alpar@1
|
571 |
/* change the threshold pivoting tolerance */
|
alpar@1
|
572 |
parm.piv_tol = 0.05;
|
alpar@1
|
573 |
/* set new values of control parameters */
|
alpar@1
|
574 |
glp_set_bfcp(lp, &parm);
|
alpar@1
|
575 |
. . .
|
alpar@1
|
576 |
\end{verbatim}
|
alpar@1
|
577 |
|
alpar@1
|
578 |
\subsubsection*{Control parameters}
|
alpar@1
|
579 |
|
alpar@1
|
580 |
This paragraph describes all basis factorization control parameters
|
alpar@1
|
581 |
currently used in the package. Symbolic names of control parameters are
|
alpar@1
|
582 |
names of corresponding members in the structure \verb|glp_bfcp|.
|
alpar@1
|
583 |
|
alpar@1
|
584 |
\def\arraystretch{1}
|
alpar@1
|
585 |
|
alpar@1
|
586 |
\medskip
|
alpar@1
|
587 |
|
alpar@1
|
588 |
\noindent\begin{tabular}{@{}p{17pt}@{}p{120.5mm}@{}}
|
alpar@1
|
589 |
\multicolumn{2}{@{}l}{{\tt int type} (default: {\tt GLP\_BF\_FT})} \\
|
alpar@1
|
590 |
&Basis factorization type:\\
|
alpar@1
|
591 |
&\verb|GLP_BF_FT|---$LU$ + Forrest--Tomlin update;\\
|
alpar@1
|
592 |
&\verb|GLP_BF_BG|---$LU$ + Schur complement + Bartels--Golub update;\\
|
alpar@1
|
593 |
&\verb|GLP_BF_GR|---$LU$ + Schur complement + Givens rotation update.
|
alpar@1
|
594 |
\\
|
alpar@1
|
595 |
&In case of \verb|GLP_BF_FT| the update is applied to matrix $U$, while
|
alpar@1
|
596 |
in cases of \verb|GLP_BF_BG| and \verb|GLP_BF_GR| the update is applied
|
alpar@1
|
597 |
to the Schur complement.
|
alpar@1
|
598 |
\end{tabular}
|
alpar@1
|
599 |
|
alpar@1
|
600 |
\medskip
|
alpar@1
|
601 |
|
alpar@1
|
602 |
\noindent\begin{tabular}{@{}p{17pt}@{}p{120.5mm}@{}}
|
alpar@1
|
603 |
\multicolumn{2}{@{}l}{{\tt int lu\_size} (default: {\tt 0})} \\
|
alpar@1
|
604 |
&The initial size of the Sparse Vector Area, in non-zeros, used on
|
alpar@1
|
605 |
computing $LU$-factorization of the basis matrix for the first time.
|
alpar@1
|
606 |
If this parameter is set to 0, the initial SVA size is determined
|
alpar@1
|
607 |
automatically.\\
|
alpar@1
|
608 |
\end{tabular}
|
alpar@1
|
609 |
|
alpar@1
|
610 |
\medskip
|
alpar@1
|
611 |
|
alpar@1
|
612 |
\noindent\begin{tabular}{@{}p{17pt}@{}p{120.5mm}@{}}
|
alpar@1
|
613 |
\multicolumn{2}{@{}l}{{\tt double piv\_tol} (default: {\tt 0.10})} \\
|
alpar@1
|
614 |
&Threshold pivoting (Markowitz) tolerance, 0 $<$ \verb|piv_tol| $<$ 1,
|
alpar@1
|
615 |
used on computing $LU$-factorization of the basis matrix. Element
|
alpar@1
|
616 |
$u_{ij}$ of the active submatrix of factor $U$ fits to be pivot if it
|
alpar@1
|
617 |
satisfies to the stability criterion
|
alpar@1
|
618 |
$|u_{ij}| >= {\tt piv\_tol}\cdot\max|u_{i*}|$, i.e. if it is not very
|
alpar@1
|
619 |
small in the magnitude among other elements in the same row. Decreasing
|
alpar@1
|
620 |
this parameter may lead to better sparsity at the expense of numerical
|
alpar@1
|
621 |
accuracy, and vice versa.\\
|
alpar@1
|
622 |
\end{tabular}
|
alpar@1
|
623 |
|
alpar@1
|
624 |
\medskip
|
alpar@1
|
625 |
|
alpar@1
|
626 |
\noindent\begin{tabular}{@{}p{17pt}@{}p{120.5mm}@{}}
|
alpar@1
|
627 |
\multicolumn{2}{@{}l}{{\tt int piv\_lim} (default: {\tt 4})} \\
|
alpar@1
|
628 |
&This parameter is used on computing $LU$-factorization of the basis
|
alpar@1
|
629 |
matrix and specifies how many pivot candidates needs to be considered
|
alpar@1
|
630 |
on choosing a pivot element, \verb|piv_lim| $\geq$ 1. If \verb|piv_lim|
|
alpar@1
|
631 |
candidates have been considered, the pivoting routine prematurely
|
alpar@1
|
632 |
terminates the search with the best candidate found.\\
|
alpar@1
|
633 |
\end{tabular}
|
alpar@1
|
634 |
|
alpar@1
|
635 |
\medskip
|
alpar@1
|
636 |
|
alpar@1
|
637 |
\noindent\begin{tabular}{@{}p{17pt}@{}p{120.5mm}@{}}
|
alpar@1
|
638 |
\multicolumn{2}{@{}l}{{\tt int suhl} (default: {\tt GLP\_ON})} \\
|
alpar@1
|
639 |
&This parameter is used on computing $LU$-factorization of the basis
|
alpar@1
|
640 |
matrix. Being set to {\tt GLP\_ON} it enables applying the following
|
alpar@1
|
641 |
heuristic proposed by Uwe Suhl: if a column of the active submatrix has
|
alpar@1
|
642 |
no eligible pivot candidates, it is no more considered until it becomes
|
alpar@1
|
643 |
a column singleton. In many cases this allows reducing the time needed
|
alpar@1
|
644 |
for pivot searching. To disable this heuristic the parameter \verb|suhl|
|
alpar@1
|
645 |
should be set to {\tt GLP\_OFF}.\\
|
alpar@1
|
646 |
\end{tabular}
|
alpar@1
|
647 |
|
alpar@1
|
648 |
\medskip
|
alpar@1
|
649 |
|
alpar@1
|
650 |
\noindent\begin{tabular}{@{}p{17pt}@{}p{120.5mm}@{}}
|
alpar@1
|
651 |
\multicolumn{2}{@{}l}{{\tt double eps\_tol} (default: {\tt 1e-15})} \\
|
alpar@1
|
652 |
&Epsilon tolerance, \verb|eps_tol| $\geq$ 0, used on computing
|
alpar@1
|
653 |
$LU$-factorization of the basis matrix. If an element of the active
|
alpar@1
|
654 |
submatrix of factor $U$ is less than \verb|eps_tol| in the magnitude,
|
alpar@1
|
655 |
it is replaced by exact zero.\\
|
alpar@1
|
656 |
\end{tabular}
|
alpar@1
|
657 |
|
alpar@1
|
658 |
\medskip
|
alpar@1
|
659 |
|
alpar@1
|
660 |
\noindent\begin{tabular}{@{}p{17pt}@{}p{120.5mm}@{}}
|
alpar@1
|
661 |
\multicolumn{2}{@{}l}{{\tt double max\_gro} (default: {\tt 1e+10})} \\
|
alpar@1
|
662 |
&Maximal growth of elements of factor $U$, \verb|max_gro| $\geq$ 1,
|
alpar@1
|
663 |
allowable on computing $LU$-factorization of the basis matrix. If on
|
alpar@1
|
664 |
some elimination step the ratio $u_{big}/b_{max}$ (where $u_{big}$ is
|
alpar@1
|
665 |
the largest magnitude of elements of factor $U$ appeared in its active
|
alpar@1
|
666 |
submatrix during all the factorization process, $b_{max}$ is the largest
|
alpar@1
|
667 |
magnitude of elements of the basis matrix to be factorized), the basis
|
alpar@1
|
668 |
matrix is considered as ill-conditioned.\\
|
alpar@1
|
669 |
\end{tabular}
|
alpar@1
|
670 |
|
alpar@1
|
671 |
\medskip
|
alpar@1
|
672 |
|
alpar@1
|
673 |
\noindent\begin{tabular}{@{}p{17pt}@{}p{120.5mm}@{}}
|
alpar@1
|
674 |
\multicolumn{2}{@{}l}{{\tt int nfs\_max} (default: {\tt 100})} \\
|
alpar@1
|
675 |
&Maximal number of additional row-like factors (entries of the eta
|
alpar@1
|
676 |
file), \verb|nfs_max| $\geq$ 1, which can be added to $LU$-factorization
|
alpar@1
|
677 |
of the basis matrix on updating it with the Forrest--Tomlin technique.
|
alpar@1
|
678 |
This parameter is used only once, before $LU$-factorization is computed
|
alpar@1
|
679 |
for the first time, to allocate working arrays. As a rule, each update
|
alpar@1
|
680 |
adds one new factor (however, some updates may need no addition), so
|
alpar@1
|
681 |
this parameter limits the number of updates between refactorizations.\\
|
alpar@1
|
682 |
\end{tabular}
|
alpar@1
|
683 |
|
alpar@1
|
684 |
\medskip
|
alpar@1
|
685 |
|
alpar@1
|
686 |
\noindent\begin{tabular}{@{}p{17pt}@{}p{120.5mm}@{}}
|
alpar@1
|
687 |
\multicolumn{2}{@{}l}{{\tt double upd\_tol} (default: {\tt 1e-6})} \\
|
alpar@1
|
688 |
&Update tolerance, 0 $<$ \verb|upd_tol| $<$ 1, used on updating
|
alpar@1
|
689 |
$LU$-factorization of the basis matrix with the Forrest--Tomlin
|
alpar@1
|
690 |
technique. If after updating the magnitude of some diagonal element
|
alpar@1
|
691 |
$u_{kk}$ of factor $U$ becomes less than
|
alpar@1
|
692 |
${\tt upd\_tol}\cdot\max(|u_{k*}|, |u_{*k}|)$, the factorization is
|
alpar@1
|
693 |
considered as inaccurate.\\
|
alpar@1
|
694 |
\end{tabular}
|
alpar@1
|
695 |
|
alpar@1
|
696 |
\medskip
|
alpar@1
|
697 |
|
alpar@1
|
698 |
\noindent\begin{tabular}{@{}p{17pt}@{}p{120.5mm}@{}}
|
alpar@1
|
699 |
\multicolumn{2}{@{}l}{{\tt int nrs\_max} (default: {\tt 100})} \\
|
alpar@1
|
700 |
&Maximal number of additional rows and columns, \verb|nrs_max| $\geq$ 1,
|
alpar@1
|
701 |
which can be added to $LU$-factorization of the basis matrix on updating
|
alpar@1
|
702 |
it with the Schur complement technique. This parameter is used only
|
alpar@1
|
703 |
once, before $LU$-factorization is computed for the first time, to
|
alpar@1
|
704 |
allocate working arrays. As a rule, each update adds one new row and
|
alpar@1
|
705 |
column (however, some updates may need no addition), so this parameter
|
alpar@1
|
706 |
limits the number of updates between refactorizations.\\
|
alpar@1
|
707 |
\end{tabular}
|
alpar@1
|
708 |
|
alpar@1
|
709 |
\medskip
|
alpar@1
|
710 |
|
alpar@1
|
711 |
\noindent\begin{tabular}{@{}p{17pt}@{}p{120.5mm}@{}}
|
alpar@1
|
712 |
\multicolumn{2}{@{}l}{{\tt int rs\_size} (default: {\tt 0})} \\
|
alpar@1
|
713 |
&The initial size of the Sparse Vector Area, in non-zeros, used to
|
alpar@1
|
714 |
store non-zero elements of additional rows and columns introduced on
|
alpar@1
|
715 |
updating $LU$-factorization of the basis matrix with the Schur
|
alpar@1
|
716 |
complement technique. If this parameter is set to 0, the initial SVA
|
alpar@1
|
717 |
size is determined automatically.\\
|
alpar@1
|
718 |
\end{tabular}
|
alpar@1
|
719 |
|
alpar@1
|
720 |
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
|
alpar@1
|
721 |
|
alpar@1
|
722 |
\newpage
|
alpar@1
|
723 |
|
alpar@1
|
724 |
\subsection{glp\_get\_bhead---retrieve the basis header information}
|
alpar@1
|
725 |
|
alpar@1
|
726 |
\subsubsection*{Synopsis}
|
alpar@1
|
727 |
|
alpar@1
|
728 |
\begin{verbatim}
|
alpar@1
|
729 |
int glp_get_bhead(glp_prob *lp, int k);
|
alpar@1
|
730 |
\end{verbatim}
|
alpar@1
|
731 |
|
alpar@1
|
732 |
\subsubsection*{Description}
|
alpar@1
|
733 |
|
alpar@1
|
734 |
The routine \verb|glp_get_bhead| returns the basis header information
|
alpar@1
|
735 |
for the current basis associated with the specified problem object.
|
alpar@1
|
736 |
|
alpar@1
|
737 |
\subsubsection*{Returns}
|
alpar@1
|
738 |
|
alpar@1
|
739 |
If basic variable $(x_B)_k$, $1\leq k\leq m$, is $i$-th auxiliary
|
alpar@1
|
740 |
variable ($1\leq i\leq m$), the routine returns $i$. Otherwise, if
|
alpar@1
|
741 |
$(x_B)_k$ is $j$-th structural variable ($1\leq j\leq n$), the routine
|
alpar@1
|
742 |
returns $m+j$. Here $m$ is the number of rows and $n$ is the number of
|
alpar@1
|
743 |
columns in the problem object.
|
alpar@1
|
744 |
|
alpar@1
|
745 |
\subsubsection*{Comments}
|
alpar@1
|
746 |
|
alpar@1
|
747 |
Sometimes the application program may need to know which original
|
alpar@1
|
748 |
(auxiliary and structural) variable correspond to a given basic
|
alpar@1
|
749 |
variable, or, that is the same, which column of the augmented constraint
|
alpar@1
|
750 |
matrix $(I\ |-\!A)$ correspond to a given column of the basis matrix
|
alpar@1
|
751 |
$B$.
|
alpar@1
|
752 |
|
alpar@1
|
753 |
\def\arraystretch{1}
|
alpar@1
|
754 |
|
alpar@1
|
755 |
The correspondence is defined as follows:\footnote{For more details see
|
alpar@1
|
756 |
Subsection \ref{basbgd}, page \pageref{basbgd}.}
|
alpar@1
|
757 |
$$\left(\begin{array}{@{}c@{}}x_B\\x_N\\\end{array}\right)=
|
alpar@1
|
758 |
\Pi\left(\begin{array}{@{}c@{}}x_R\\x_S\\\end{array}\right)
|
alpar@1
|
759 |
\ \ \Leftrightarrow
|
alpar@1
|
760 |
\ \ \left(\begin{array}{@{}c@{}}x_R\\x_S\\\end{array}\right)=
|
alpar@1
|
761 |
\Pi^T\left(\begin{array}{@{}c@{}}x_B\\x_N\\\end{array}\right),$$
|
alpar@1
|
762 |
where $x_B$ is the vector of basic variables, $x_N$ is the vector of
|
alpar@1
|
763 |
non-basic variables, $x_R$ is the vector of auxiliary variables
|
alpar@1
|
764 |
following in their original order,\footnote{The original order of
|
alpar@1
|
765 |
auxiliary and structural variables is defined by the ordinal numbers
|
alpar@1
|
766 |
of corresponding rows and columns in the problem object.} $x_S$ is the
|
alpar@1
|
767 |
vector of structural variables following in their original order, $\Pi$
|
alpar@1
|
768 |
is a permutation matrix (which is a component of the basis
|
alpar@1
|
769 |
factorization).
|
alpar@1
|
770 |
|
alpar@1
|
771 |
Thus, if $(x_B)_k=(x_R)_i$ is $i$-th auxiliary variable, the routine
|
alpar@1
|
772 |
returns $i$, and if $(x_B)_k=(x_S)_j$ is $j$-th structural variable,
|
alpar@1
|
773 |
the routine returns $m+j$, where $m$ is the number of rows in the
|
alpar@1
|
774 |
problem object.
|
alpar@1
|
775 |
|
alpar@1
|
776 |
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
|
alpar@1
|
777 |
|
alpar@1
|
778 |
\newpage
|
alpar@1
|
779 |
|
alpar@1
|
780 |
\subsection{glp\_get\_row\_bind---retrieve row index in the basis\\
|
alpar@1
|
781 |
header}
|
alpar@1
|
782 |
|
alpar@1
|
783 |
\subsubsection*{Synopsis}
|
alpar@1
|
784 |
|
alpar@1
|
785 |
\begin{verbatim}
|
alpar@1
|
786 |
int glp_get_row_bind(glp_prob *lp, int i);
|
alpar@1
|
787 |
\end{verbatim}
|
alpar@1
|
788 |
|
alpar@1
|
789 |
\subsubsection*{Returns}
|
alpar@1
|
790 |
|
alpar@1
|
791 |
The routine \verb|glp_get_row_bind| returns the index $k$ of basic
|
alpar@1
|
792 |
variable $(x_B)_k$, $1\leq k\leq m$, which is $i$-th auxiliary variable
|
alpar@1
|
793 |
(that is, the auxiliary variable corresponding to $i$-th row),
|
alpar@1
|
794 |
$1\leq i\leq m$, in the current basis associated with the specified
|
alpar@1
|
795 |
problem object, where $m$ is the number of rows. However, if $i$-th
|
alpar@1
|
796 |
auxiliary variable is non-basic, the routine returns zero.
|
alpar@1
|
797 |
|
alpar@1
|
798 |
\subsubsection*{Comments}
|
alpar@1
|
799 |
|
alpar@1
|
800 |
The routine \verb|glp_get_row_bind| is an inverse to the routine
|
alpar@1
|
801 |
\verb|glp_get_bhead|: if \verb|glp_get_bhead|$(lp,k)$ returns $i$,
|
alpar@1
|
802 |
\verb|glp_get_row_bind|$(lp,i)$ returns $k$, and vice versa.
|
alpar@1
|
803 |
|
alpar@1
|
804 |
\subsection{glp\_get\_col\_bind---retrieve column index in the basis
|
alpar@1
|
805 |
header}
|
alpar@1
|
806 |
|
alpar@1
|
807 |
\subsubsection*{Synopsis}
|
alpar@1
|
808 |
|
alpar@1
|
809 |
\begin{verbatim}
|
alpar@1
|
810 |
int glp_get_col_bind(glp_prob *lp, int j);
|
alpar@1
|
811 |
\end{verbatim}
|
alpar@1
|
812 |
|
alpar@1
|
813 |
\subsubsection*{Returns}
|
alpar@1
|
814 |
|
alpar@1
|
815 |
The routine \verb|glp_get_col_bind| returns the index $k$ of basic
|
alpar@1
|
816 |
variable $(x_B)_k$, $1\leq k\leq m$, which is $j$-th structural
|
alpar@1
|
817 |
variable (that is, the structural variable corresponding to $j$-th
|
alpar@1
|
818 |
column), $1\leq j\leq n$, in the current basis associated with the
|
alpar@1
|
819 |
specified problem object, where $m$ is the number of rows, $n$ is the
|
alpar@1
|
820 |
number of columns. However, if $j$-th structural variable is non-basic,
|
alpar@1
|
821 |
the routine returns zero.
|
alpar@1
|
822 |
|
alpar@1
|
823 |
\subsubsection*{Comments}
|
alpar@1
|
824 |
|
alpar@1
|
825 |
The routine \verb|glp_get_col_bind| is an inverse to the routine
|
alpar@1
|
826 |
\verb|glp_get_bhead|: if \verb|glp_get_bhead|$(lp,k)$ returns $m+j$,
|
alpar@1
|
827 |
\verb|glp_get_col_bind|$(lp,j)$ returns $k$, and vice versa.
|
alpar@1
|
828 |
|
alpar@1
|
829 |
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
|
alpar@1
|
830 |
|
alpar@1
|
831 |
\newpage
|
alpar@1
|
832 |
|
alpar@1
|
833 |
\subsection{glp\_ftran---perform forward transformation}
|
alpar@1
|
834 |
|
alpar@1
|
835 |
\subsubsection*{Synopsis}
|
alpar@1
|
836 |
|
alpar@1
|
837 |
\begin{verbatim}
|
alpar@1
|
838 |
void glp_ftran(glp_prob *lp, double x[]);
|
alpar@1
|
839 |
\end{verbatim}
|
alpar@1
|
840 |
|
alpar@1
|
841 |
\subsubsection*{Description}
|
alpar@1
|
842 |
|
alpar@1
|
843 |
The routine \verb|glp_ftran| performs forward transformation (FTRAN),
|
alpar@1
|
844 |
i.e. it solves the system $Bx=b$, where $B$ is the basis matrix
|
alpar@1
|
845 |
associated with the specified problem object, $x$ is the vector of
|
alpar@1
|
846 |
unknowns to be computed, $b$ is the vector of right-hand sides.
|
alpar@1
|
847 |
|
alpar@1
|
848 |
On entry to the routine elements of the vector $b$ should be stored in
|
alpar@1
|
849 |
locations \verb|x[1]|, \dots, \verb|x[m]|, where $m$ is the number of
|
alpar@1
|
850 |
rows. On exit the routine stores elements of the vector $x$ in the same
|
alpar@1
|
851 |
locations.
|
alpar@1
|
852 |
|
alpar@1
|
853 |
\subsection{glp\_btran---perform backward transformation}
|
alpar@1
|
854 |
|
alpar@1
|
855 |
\subsubsection*{Synopsis}
|
alpar@1
|
856 |
|
alpar@1
|
857 |
\begin{verbatim}
|
alpar@1
|
858 |
void glp_btran(glp_prob *lp, double x[]);
|
alpar@1
|
859 |
\end{verbatim}
|
alpar@1
|
860 |
|
alpar@1
|
861 |
\subsubsection*{Description}
|
alpar@1
|
862 |
|
alpar@1
|
863 |
The routine \verb|glp_btran| performs backward transformation (BTRAN),
|
alpar@1
|
864 |
i.e. it solves the system $B^Tx=b$, where $B^T$ is a matrix transposed
|
alpar@1
|
865 |
to the basis matrix $B$ associated with the specified problem object,
|
alpar@1
|
866 |
$x$ is the vector of unknowns to be computed, $b$ is the vector of
|
alpar@1
|
867 |
right-hand sides.
|
alpar@1
|
868 |
|
alpar@1
|
869 |
On entry to the routine elements of the vector $b$ should be stored in
|
alpar@1
|
870 |
locations \verb|x[1]|, \dots, \verb|x[m]|, where $m$ is the number of
|
alpar@1
|
871 |
rows. On exit the routine stores elements of the vector $x$ in the same
|
alpar@1
|
872 |
locations.
|
alpar@1
|
873 |
|
alpar@1
|
874 |
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
|
alpar@1
|
875 |
|
alpar@1
|
876 |
\newpage
|
alpar@1
|
877 |
|
alpar@1
|
878 |
\subsection{glp\_warm\_up---``warm up'' LP basis}
|
alpar@1
|
879 |
|
alpar@1
|
880 |
\subsubsection*{Synopsis}
|
alpar@1
|
881 |
|
alpar@1
|
882 |
\begin{verbatim}
|
alpar@1
|
883 |
int glp_warm_up(glp_prob *P);
|
alpar@1
|
884 |
\end{verbatim}
|
alpar@1
|
885 |
|
alpar@1
|
886 |
\subsubsection*{Description}
|
alpar@1
|
887 |
|
alpar@1
|
888 |
The routine \verb|glp_warm_up| ``warms up'' the LP basis for the
|
alpar@1
|
889 |
specified problem object using current statuses assigned to rows and
|
alpar@1
|
890 |
columns (that is, to auxiliary and structural variables).
|
alpar@1
|
891 |
|
alpar@1
|
892 |
This operation includes computing factorization of the basis matrix
|
alpar@1
|
893 |
(if it does not exist), computing primal and dual components of basic
|
alpar@1
|
894 |
solution, and determining the solution status.
|
alpar@1
|
895 |
|
alpar@1
|
896 |
\subsubsection*{Returns}
|
alpar@1
|
897 |
|
alpar@1
|
898 |
\begin{tabular}{@{}p{25mm}p{97.3mm}@{}}
|
alpar@1
|
899 |
0 & The operation has been successfully performed.\\
|
alpar@1
|
900 |
\verb|GLP_EBADB| & The basis matrix is invalid, because the number of
|
alpar@1
|
901 |
basic (auxiliary and structural) variables is not the same as the number
|
alpar@1
|
902 |
of rows in the problem object.\\
|
alpar@1
|
903 |
\verb|GLP_ESING| & The basis matrix is singular within the working
|
alpar@1
|
904 |
precision.\\
|
alpar@1
|
905 |
\verb|GLP_ECOND| & The basis matrix is ill-conditioned, i.e. its
|
alpar@1
|
906 |
condition number is too large.\\
|
alpar@1
|
907 |
\end{tabular}
|
alpar@1
|
908 |
|
alpar@1
|
909 |
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
|
alpar@1
|
910 |
|
alpar@1
|
911 |
\newpage
|
alpar@1
|
912 |
|
alpar@1
|
913 |
\section{Simplex tableau routines}
|
alpar@1
|
914 |
|
alpar@1
|
915 |
\subsection{glp\_eval\_tab\_row---compute row of the tableau}
|
alpar@1
|
916 |
|
alpar@1
|
917 |
\subsubsection*{Synopsis}
|
alpar@1
|
918 |
|
alpar@1
|
919 |
\begin{verbatim}
|
alpar@1
|
920 |
int glp_eval_tab_row(glp_prob *lp, int k, int ind[],
|
alpar@1
|
921 |
double val[]);
|
alpar@1
|
922 |
\end{verbatim}
|
alpar@1
|
923 |
|
alpar@1
|
924 |
\subsubsection*{Description}
|
alpar@1
|
925 |
|
alpar@1
|
926 |
The routine \verb|glp_eval_tab_row| computes a row of the current
|
alpar@1
|
927 |
simplex tableau (see Subsection 3.1.1, formula (3.12)), which (row)
|
alpar@1
|
928 |
corresponds to some basic variable specified by the parameter $k$ as
|
alpar@1
|
929 |
follows: if $1\leq k\leq m$, the basic variable is $k$-th auxiliary
|
alpar@1
|
930 |
variable, and if $m+1\leq k\leq m+n$, the basic variable is $(k-m)$-th
|
alpar@1
|
931 |
structural variable, where $m$ is the number of rows and $n$ is the
|
alpar@1
|
932 |
number of columns in the specified problem object. The basis
|
alpar@1
|
933 |
factorization must exist.
|
alpar@1
|
934 |
|
alpar@1
|
935 |
The computed row shows how the specified basic variable depends on
|
alpar@1
|
936 |
non-basic variables:
|
alpar@1
|
937 |
$$x_k=(x_B)_i=\xi_{i1}(x_N)_1+\xi_{i2}(x_N)_2+\dots+\xi_{in}(x_N)_n,$$
|
alpar@1
|
938 |
where $\xi_{i1}$, $\xi_{i2}$, \dots, $\xi_{in}$ are elements of the
|
alpar@1
|
939 |
simplex table row, $(x_N)_1$, $(x_N)_2$, \dots, $(x_N)_n$ are non-basic
|
alpar@1
|
940 |
(auxiliary and structural) variables.
|
alpar@1
|
941 |
|
alpar@1
|
942 |
The routine stores column indices and corresponding numeric values of
|
alpar@1
|
943 |
non-zero elements of the computed row in unordered sparse format in
|
alpar@1
|
944 |
locations \verb|ind[1]|, \dots, \verb|ind[len]| and \verb|val[1]|,
|
alpar@1
|
945 |
\dots, \verb|val[len]|, respectively, where $0\leq{\tt len}\leq n$ is
|
alpar@1
|
946 |
the number of non-zero elements in the row returned on exit.
|
alpar@1
|
947 |
|
alpar@1
|
948 |
Element indices stored in the array \verb|ind| have the same sense as
|
alpar@1
|
949 |
index $k$, i.e. indices 1 to $m$ denote auxiliary variables while
|
alpar@1
|
950 |
indices $m+1$ to $m+n$ denote structural variables (all these variables
|
alpar@1
|
951 |
are obviously non-basic by definition).
|
alpar@1
|
952 |
|
alpar@1
|
953 |
\subsubsection*{Returns}
|
alpar@1
|
954 |
|
alpar@1
|
955 |
The routine \verb|glp_eval_tab_row| returns \verb|len|, which is the
|
alpar@1
|
956 |
number of non-zero elements in the simplex table row stored in the
|
alpar@1
|
957 |
arrays \verb|ind| and \verb|val|.
|
alpar@1
|
958 |
|
alpar@1
|
959 |
\subsubsection*{Comments}
|
alpar@1
|
960 |
|
alpar@1
|
961 |
A row of the simplex table is computed as follows. At first, the
|
alpar@1
|
962 |
routine checks that the specified variable $x_k$ is basic and uses the
|
alpar@1
|
963 |
permutation matrix $\Pi$ (3.7) to determine index $i$ of basic variable
|
alpar@1
|
964 |
$(x_B)_i$, which corresponds to $x_k$.
|
alpar@1
|
965 |
|
alpar@1
|
966 |
The row to be computed is $i$-th row of the matrix $\Xi$ (3.12),
|
alpar@1
|
967 |
therefore:
|
alpar@1
|
968 |
$$\xi_i=e_i^T\Xi=-e_i^TB^{-1}N=-(B^{-T}e_i)^TN,$$
|
alpar@1
|
969 |
where $e_i$ is $i$-th unity vector. So the routine performs BTRAN to
|
alpar@1
|
970 |
obtain $i$-th row of the inverse $B^{-1}$:
|
alpar@1
|
971 |
$$\varrho_i=B^{-T}e_i,$$
|
alpar@1
|
972 |
and then computes elements of the simplex table row as inner products:
|
alpar@1
|
973 |
$$\xi_{ij}=-\varrho_i^TN_j,\ \ j=1,2,\dots,n,$$
|
alpar@1
|
974 |
where $N_j$ is $j$-th column of matrix $N$ (3.9), which (column)
|
alpar@1
|
975 |
corresponds to non-basic variable $(x_N)_j$. The permutation matrix
|
alpar@1
|
976 |
$\Pi$ is used again to convert indices $j$ of non-basic columns to
|
alpar@1
|
977 |
original ordinal numbers of auxiliary and structural variables.
|
alpar@1
|
978 |
|
alpar@1
|
979 |
\subsection{glp\_eval\_tab\_col---compute column of the tableau}
|
alpar@1
|
980 |
|
alpar@1
|
981 |
\subsubsection*{Synopsis}
|
alpar@1
|
982 |
|
alpar@1
|
983 |
\begin{verbatim}
|
alpar@1
|
984 |
int glp_eval_tab_col(glp_prob *lp, int k, int ind[],
|
alpar@1
|
985 |
double val[]);
|
alpar@1
|
986 |
\end{verbatim}
|
alpar@1
|
987 |
|
alpar@1
|
988 |
\subsubsection*{Description}
|
alpar@1
|
989 |
|
alpar@1
|
990 |
The routine \verb|glp_eval_tab_col| computes a column of the current
|
alpar@1
|
991 |
simplex tableau (see Subsection 3.1.1, formula (3.12)), which (column)
|
alpar@1
|
992 |
corresponds to some non-basic variable specified by the parameter $k$:
|
alpar@1
|
993 |
if $1\leq k\leq m$, the non-basic variable is $k$-th auxiliary variable,
|
alpar@1
|
994 |
and if $m+1\leq k\leq m+n$, the non-basic variable is $(k-m)$-th
|
alpar@1
|
995 |
structural variable, where $m$ is the number of rows and $n$ is the
|
alpar@1
|
996 |
number of columns in the specified problem object. The basis
|
alpar@1
|
997 |
factorization must exist.
|
alpar@1
|
998 |
|
alpar@1
|
999 |
The computed column shows how basic variables depends on the specified
|
alpar@1
|
1000 |
non-basic variable $x_k=(x_N)_j$:
|
alpar@1
|
1001 |
$$
|
alpar@1
|
1002 |
\begin{array}{r@{\ }c@{\ }l@{\ }l}
|
alpar@1
|
1003 |
(x_B)_1&=&\dots+\xi_{1j}(x_N)_j&+\dots\\
|
alpar@1
|
1004 |
(x_B)_2&=&\dots+\xi_{2j}(x_N)_j&+\dots\\
|
alpar@1
|
1005 |
.\ \ .&.&.\ \ .\ \ .\ \ .\ \ .\ \ .\ \ .\\
|
alpar@1
|
1006 |
(x_B)_m&=&\dots+\xi_{mj}(x_N)_j&+\dots\\
|
alpar@1
|
1007 |
\end{array}
|
alpar@1
|
1008 |
$$
|
alpar@1
|
1009 |
where $\xi_{1j}$, $\xi_{2j}$, \dots, $\xi_{mj}$ are elements of the
|
alpar@1
|
1010 |
simplex table column, $(x_B)_1$, $(x_B)_2$, \dots, $(x_B)_m$ are basic
|
alpar@1
|
1011 |
(auxiliary and structural) variables.
|
alpar@1
|
1012 |
|
alpar@1
|
1013 |
The routine stores row indices and corresponding numeric values of
|
alpar@1
|
1014 |
non-zero elements of the computed column in unordered sparse format in
|
alpar@1
|
1015 |
locations \verb|ind[1]|, \dots, \verb|ind[len]| and \verb|val[1]|,
|
alpar@1
|
1016 |
\dots, \verb|val[len]|, respectively, where $0\leq{\tt len}\leq m$ is
|
alpar@1
|
1017 |
the number of non-zero elements in the column returned on exit.
|
alpar@1
|
1018 |
|
alpar@1
|
1019 |
Element indices stored in the array \verb|ind| have the same sense as
|
alpar@1
|
1020 |
index $k$, i.e. indices 1 to $m$ denote auxiliary variables while
|
alpar@1
|
1021 |
indices $m+1$ to $m+n$ denote structural variables (all these variables
|
alpar@1
|
1022 |
are obviously basic by definition).
|
alpar@1
|
1023 |
|
alpar@1
|
1024 |
\subsubsection*{Returns}
|
alpar@1
|
1025 |
|
alpar@1
|
1026 |
The routine \verb|glp_eval_tab_col| returns \verb|len|, which is the
|
alpar@1
|
1027 |
number of non-zero elements in the simplex table column stored in the
|
alpar@1
|
1028 |
arrays \verb|ind| and \verb|val|.
|
alpar@1
|
1029 |
|
alpar@1
|
1030 |
\subsubsection*{Comments}
|
alpar@1
|
1031 |
|
alpar@1
|
1032 |
A column of the simplex table is computed as follows. At first, the
|
alpar@1
|
1033 |
routine checks that the specified variable $x_k$ is non-basic and uses
|
alpar@1
|
1034 |
the permutation matrix $\Pi$ (3.7) to determine index $j$ of non-basic
|
alpar@1
|
1035 |
variable $(x_N)_j$, which corresponds to $x_k$.
|
alpar@1
|
1036 |
|
alpar@1
|
1037 |
The column to be computed is $j$-th column of the matrix $\Xi$ (3.12),
|
alpar@1
|
1038 |
therefore:
|
alpar@1
|
1039 |
$$\Xi_j=\Xi e_j=-B^{-1}Ne_j=-B^{-1}N_j,$$
|
alpar@1
|
1040 |
where $e_j$ is $j$-th unity vector, $N_j$ is $j$-th column of matrix
|
alpar@1
|
1041 |
$N$ (3.9). So the routine performs FTRAN to transform $N_j$ to the
|
alpar@1
|
1042 |
simplex table column $\Xi_j=(\xi_{ij})$ and uses the permutation matrix
|
alpar@1
|
1043 |
$\Pi$ to convert row indices $i$ to original ordinal numbers of
|
alpar@1
|
1044 |
auxiliary and structural variables.
|
alpar@1
|
1045 |
|
alpar@1
|
1046 |
\newpage
|
alpar@1
|
1047 |
|
alpar@1
|
1048 |
\subsection{glp\_transform\_row---transform explicitly specified\\
|
alpar@1
|
1049 |
row}
|
alpar@1
|
1050 |
|
alpar@1
|
1051 |
\subsubsection*{Synopsis}
|
alpar@1
|
1052 |
|
alpar@1
|
1053 |
\begin{verbatim}
|
alpar@1
|
1054 |
int glp_transform_row(glp_prob *P, int len, int ind[],
|
alpar@1
|
1055 |
double val[]);
|
alpar@1
|
1056 |
\end{verbatim}
|
alpar@1
|
1057 |
|
alpar@1
|
1058 |
\subsubsection*{Description}
|
alpar@1
|
1059 |
|
alpar@1
|
1060 |
The routine \verb|glp_transform_row| performs the same operation as the
|
alpar@1
|
1061 |
routine \verb|glp_eval_tab_row| with exception that the row to be
|
alpar@1
|
1062 |
transformed is specified explicitly as a sparse vector.
|
alpar@1
|
1063 |
|
alpar@1
|
1064 |
The explicitly specified row may be thought as a linear form:
|
alpar@1
|
1065 |
$$x=a_1x_{m+1}+a_2x_{m+2}+\dots+a_nx_{m+n},$$
|
alpar@1
|
1066 |
where $x$ is an auxiliary variable for this row, $a_j$ are coefficients
|
alpar@1
|
1067 |
of the linear form, $x_{m+j}$ are structural variables.
|
alpar@1
|
1068 |
|
alpar@1
|
1069 |
On entry column indices and numerical values of non-zero coefficients
|
alpar@1
|
1070 |
$a_j$ of the specified row should be placed in locations \verb|ind[1]|,
|
alpar@1
|
1071 |
\dots, \verb|ind[len]| and \verb|val[1]|, \dots, \verb|val[len]|, where
|
alpar@1
|
1072 |
\verb|len| is number of non-zero coefficients.
|
alpar@1
|
1073 |
|
alpar@1
|
1074 |
This routine uses the system of equality constraints and the current
|
alpar@1
|
1075 |
basis in order to express the auxiliary variable $x$ through the current
|
alpar@1
|
1076 |
non-basic variables (as if the transformed row were added to the problem
|
alpar@1
|
1077 |
object and the auxiliary variable $x$ were basic), i.e. the resultant
|
alpar@1
|
1078 |
row has the form:
|
alpar@1
|
1079 |
$$x=\xi_1(x_N)_1+\xi_2(x_N)_2+\dots+\xi_n(x_N)_n,$$
|
alpar@1
|
1080 |
where $\xi_j$ are influence coefficients, $(x_N)_j$ are non-basic
|
alpar@1
|
1081 |
(auxiliary and structural) variables, $n$ is the number of columns in
|
alpar@1
|
1082 |
the problem object.
|
alpar@1
|
1083 |
|
alpar@1
|
1084 |
On exit the routine stores indices and numerical values of non-zero
|
alpar@1
|
1085 |
coefficients $\xi_j$ of the resultant row in locations \verb|ind[1]|,
|
alpar@1
|
1086 |
\dots, \verb|ind[len']| and \verb|val[1]|, \dots, \verb|val[len']|,
|
alpar@1
|
1087 |
where $0\leq{\tt len'}\leq n$ is the number of non-zero coefficients in
|
alpar@1
|
1088 |
the resultant row returned by the routine. Note that indices of
|
alpar@1
|
1089 |
non-basic variables stored in the array \verb|ind| correspond to
|
alpar@1
|
1090 |
original ordinal numbers of variables: indices 1 to $m$ mean auxiliary
|
alpar@1
|
1091 |
variables and indices $m+1$ to $m+n$ mean structural ones.
|
alpar@1
|
1092 |
|
alpar@1
|
1093 |
\subsubsection*{Returns}
|
alpar@1
|
1094 |
|
alpar@1
|
1095 |
The routine \verb|glp_transform_row| returns \verb|len'|, the number of
|
alpar@1
|
1096 |
non-zero coefficients in the resultant row stored in the arrays
|
alpar@1
|
1097 |
\verb|ind| and \verb|val|.
|
alpar@1
|
1098 |
|
alpar@1
|
1099 |
\subsection{glp\_transform\_col---transform explicitly specified\\
|
alpar@1
|
1100 |
column}
|
alpar@1
|
1101 |
|
alpar@1
|
1102 |
\subsubsection*{Synopsis}
|
alpar@1
|
1103 |
|
alpar@1
|
1104 |
\begin{verbatim}
|
alpar@1
|
1105 |
int glp_transform_col(glp_prob *P, int len, int ind[],
|
alpar@1
|
1106 |
double val[]);
|
alpar@1
|
1107 |
\end{verbatim}
|
alpar@1
|
1108 |
|
alpar@1
|
1109 |
\subsubsection*{Description}
|
alpar@1
|
1110 |
|
alpar@1
|
1111 |
The routine \verb|glp_transform_col| performs the same operation as the
|
alpar@1
|
1112 |
routine \verb|glp_eval_tab_col| with exception that the column to be
|
alpar@1
|
1113 |
transformed is specified explicitly as a sparse vector.
|
alpar@1
|
1114 |
|
alpar@1
|
1115 |
The explicitly specified column may be thought as it were added to
|
alpar@1
|
1116 |
the original system of equality constraints:
|
alpar@1
|
1117 |
$$
|
alpar@1
|
1118 |
\begin{array}{l@{\ }c@{\ }r@{\ }c@{\ }r@{\ }c@{\ }r}
|
alpar@1
|
1119 |
x_1&=&a_{11}x_{m+1}&+\dots+&a_{1n}x_{m+n}&+&a_1x \\
|
alpar@1
|
1120 |
x_2&=&a_{21}x_{m+1}&+\dots+&a_{2n}x_{m+n}&+&a_2x \\
|
alpar@1
|
1121 |
\multicolumn{7}{c}
|
alpar@1
|
1122 |
{.\ \ .\ \ .\ \ .\ \ .\ \ .\ \ .\ \ .\ \ .\ \ .\ \ .\ \ .\ \ .\ \ .}\\
|
alpar@1
|
1123 |
x_m&=&a_{m1}x_{m+1}&+\dots+&a_{mn}x_{m+n}&+&a_mx \\
|
alpar@1
|
1124 |
\end{array}
|
alpar@1
|
1125 |
$$
|
alpar@1
|
1126 |
where $x_i$ are auxiliary variables, $x_{m+j}$ are structural variables
|
alpar@1
|
1127 |
(presented in the problem object), $x$ is a structural variable for the
|
alpar@1
|
1128 |
explicitly specified column, $a_i$ are constraint coefficients at $x$.
|
alpar@1
|
1129 |
|
alpar@1
|
1130 |
On entry row indices and numerical values of non-zero coefficients
|
alpar@1
|
1131 |
$a_i$ of the specified column should be placed in locations
|
alpar@1
|
1132 |
\verb|ind[1]|, \dots, \verb|ind[len]| and \verb|val[1]|, \dots,
|
alpar@1
|
1133 |
\verb|val[len]|, where \verb|len| is number of non-zero coefficients.
|
alpar@1
|
1134 |
|
alpar@1
|
1135 |
This routine uses the system of equality constraints and the current
|
alpar@1
|
1136 |
basis in order to express the current basic variables through the
|
alpar@1
|
1137 |
structural variable $x$ (as if the transformed column were added to the
|
alpar@1
|
1138 |
problem object and the variable $x$ were non-basic):
|
alpar@1
|
1139 |
$$
|
alpar@1
|
1140 |
\begin{array}{l@{\ }c@{\ }r}
|
alpar@1
|
1141 |
(x_B)_1&=\dots+&\xi_{1}x\\
|
alpar@1
|
1142 |
(x_B)_2&=\dots+&\xi_{2}x\\
|
alpar@1
|
1143 |
\multicolumn{3}{c}{.\ \ .\ \ .\ \ .\ \ .\ \ .}\\
|
alpar@1
|
1144 |
(x_B)_m&=\dots+&\xi_{m}x\\
|
alpar@1
|
1145 |
\end{array}
|
alpar@1
|
1146 |
$$
|
alpar@1
|
1147 |
where $\xi_i$ are influence coefficients, $x_B$ are basic (auxiliary
|
alpar@1
|
1148 |
and structural) variables, $m$ is the number of rows in the problem
|
alpar@1
|
1149 |
object.
|
alpar@1
|
1150 |
|
alpar@1
|
1151 |
On exit the routine stores indices and numerical values of non-zero
|
alpar@1
|
1152 |
coefficients $\xi_i$ of the resultant column in locations \verb|ind[1]|,
|
alpar@1
|
1153 |
\dots, \verb|ind[len']| and \verb|val[1]|, \dots, \verb|val[len']|,
|
alpar@1
|
1154 |
where $0\leq{\tt len'}\leq m$ is the number of non-zero coefficients in
|
alpar@1
|
1155 |
the resultant column returned by the routine. Note that indices of basic
|
alpar@1
|
1156 |
variables stored in the array \verb|ind| correspond to original ordinal
|
alpar@1
|
1157 |
numbers of variables, i.e. indices 1 to $m$ mean auxiliary variables,
|
alpar@1
|
1158 |
indices $m+1$ to $m+n$ mean structural ones.
|
alpar@1
|
1159 |
|
alpar@1
|
1160 |
\subsubsection*{Returns}
|
alpar@1
|
1161 |
|
alpar@1
|
1162 |
The routine \verb|glp_transform_col| returns \verb|len'|, the number of
|
alpar@1
|
1163 |
non-zero coefficients in the resultant column stored in the arrays
|
alpar@1
|
1164 |
\verb|ind| and \verb|val|.
|
alpar@1
|
1165 |
|
alpar@1
|
1166 |
\subsection{glp\_prim\_rtest---perform primal ratio test}
|
alpar@1
|
1167 |
|
alpar@1
|
1168 |
\subsubsection*{Synopsis}
|
alpar@1
|
1169 |
|
alpar@1
|
1170 |
\begin{verbatim}
|
alpar@1
|
1171 |
int glp_prim_rtest(glp_prob *P, int len, const int ind[],
|
alpar@1
|
1172 |
const double val[], int dir, double eps);
|
alpar@1
|
1173 |
\end{verbatim}
|
alpar@1
|
1174 |
|
alpar@1
|
1175 |
\subsubsection*{Description}
|
alpar@1
|
1176 |
|
alpar@1
|
1177 |
The routine \verb|glp_prim_rtest| performs the primal ratio test using
|
alpar@1
|
1178 |
an explicitly specified column of the simplex table.
|
alpar@1
|
1179 |
|
alpar@1
|
1180 |
The current basic solution associated with the LP problem object must be
|
alpar@1
|
1181 |
primal feasible.
|
alpar@1
|
1182 |
|
alpar@1
|
1183 |
The explicitly specified column of the simplex table shows how the basic
|
alpar@1
|
1184 |
variables $x_B$ depend on some non-basic variable $x$ (which is not
|
alpar@1
|
1185 |
necessarily presented in the problem object):
|
alpar@1
|
1186 |
$$
|
alpar@1
|
1187 |
\begin{array}{l@{\ }c@{\ }r}
|
alpar@1
|
1188 |
(x_B)_1&=\dots+&\xi_{1}x\\
|
alpar@1
|
1189 |
(x_B)_2&=\dots+&\xi_{2}x\\
|
alpar@1
|
1190 |
\multicolumn{3}{c}{.\ \ .\ \ .\ \ .\ \ .\ \ .}\\
|
alpar@1
|
1191 |
(x_B)_m&=\dots+&\xi_{m}x\\
|
alpar@1
|
1192 |
\end{array}
|
alpar@1
|
1193 |
$$
|
alpar@1
|
1194 |
|
alpar@1
|
1195 |
The column is specifed on entry to the routine in sparse format. Ordinal
|
alpar@1
|
1196 |
numbers of basic variables $(x_B)_i$ should be placed in locations
|
alpar@1
|
1197 |
\verb|ind[1]|, \dots, \verb|ind[len]|, where ordinal number 1 to $m$
|
alpar@1
|
1198 |
denote auxiliary variables, and ordinal numbers $m+1$ to $m+n$ denote
|
alpar@1
|
1199 |
structural variables. The corresponding non-zero coefficients $\xi_i$
|
alpar@1
|
1200 |
should be placed in locations \verb|val[1]|, \dots, \verb|val[len]|. The
|
alpar@1
|
1201 |
arrays \verb|ind| and \verb|val| are not changed by the routine.
|
alpar@1
|
1202 |
|
alpar@1
|
1203 |
The parameter \verb|dir| specifies direction in which the variable $x$
|
alpar@1
|
1204 |
changes on entering the basis: $+1$ means increasing, $-1$ means
|
alpar@1
|
1205 |
decreasing.
|
alpar@1
|
1206 |
|
alpar@1
|
1207 |
The parameter \verb|eps| is an absolute tolerance (small positive
|
alpar@1
|
1208 |
number, say, $10^{-9}$) used by the routine to skip $\xi_i$'s whose
|
alpar@1
|
1209 |
magnitude is less than \verb|eps|.
|
alpar@1
|
1210 |
|
alpar@1
|
1211 |
The routine determines which basic variable (among those specified in
|
alpar@1
|
1212 |
\verb|ind[1]|, \dots, \verb|ind[len]|) reaches its (lower or upper)
|
alpar@1
|
1213 |
bound first before any other basic variables do, and which, therefore,
|
alpar@1
|
1214 |
should leave the basis in order to keep primal feasibility.
|
alpar@1
|
1215 |
|
alpar@1
|
1216 |
\subsubsection*{Returns}
|
alpar@1
|
1217 |
|
alpar@1
|
1218 |
The routine \verb|glp_prim_rtest| returns the index, \verb|piv|, in the
|
alpar@1
|
1219 |
arrays \verb|ind| and \verb|val| corresponding to the pivot element
|
alpar@1
|
1220 |
chosen, $1\leq$ \verb|piv| $\leq$ \verb|len|. If the adjacent basic
|
alpar@1
|
1221 |
solution is primal unbounded, and therefore the choice cannot be made,
|
alpar@1
|
1222 |
the routine returns zero.
|
alpar@1
|
1223 |
|
alpar@1
|
1224 |
\subsubsection*{Comments}
|
alpar@1
|
1225 |
|
alpar@1
|
1226 |
If the non-basic variable $x$ is presented in the LP problem object, the
|
alpar@1
|
1227 |
input column can be computed with the routine \verb|glp_eval_tab_col|;
|
alpar@1
|
1228 |
otherwise, it can be computed with the routine \verb|glp_transform_col|.
|
alpar@1
|
1229 |
|
alpar@1
|
1230 |
\subsection{glp\_dual\_rtest---perform dual ratio test}
|
alpar@1
|
1231 |
|
alpar@1
|
1232 |
\subsubsection*{Synopsis}
|
alpar@1
|
1233 |
|
alpar@1
|
1234 |
\begin{verbatim}
|
alpar@1
|
1235 |
int glp_dual_rtest(glp_prob *P, int len, const int ind[],
|
alpar@1
|
1236 |
const double val[], int dir, double eps);
|
alpar@1
|
1237 |
\end{verbatim}
|
alpar@1
|
1238 |
|
alpar@1
|
1239 |
\subsubsection*{Description}
|
alpar@1
|
1240 |
|
alpar@1
|
1241 |
The routine \verb|glp_dual_rtest| performs the dual ratio test using
|
alpar@1
|
1242 |
an explicitly specified row of the simplex table.
|
alpar@1
|
1243 |
|
alpar@1
|
1244 |
The current basic solution associated with the LP problem object must be
|
alpar@1
|
1245 |
dual feasible.
|
alpar@1
|
1246 |
|
alpar@1
|
1247 |
The explicitly specified row of the simplex table is a linear form
|
alpar@1
|
1248 |
that shows how some basic variable $x$ (which is not necessarily
|
alpar@1
|
1249 |
presented in the problem object) depends on non-basic variables $x_N$:
|
alpar@1
|
1250 |
$$x=\xi_1(x_N)_1+\xi_2(x_N)_2+\dots+\xi_n(x_N)_n.$$
|
alpar@1
|
1251 |
|
alpar@1
|
1252 |
The row is specified on entry to the routine in sparse format. Ordinal
|
alpar@1
|
1253 |
numbers of non-basic variables $(x_N)_j$ should be placed in locations
|
alpar@1
|
1254 |
\verb|ind[1]|, \dots, \verb|ind[len]|, where ordinal numbers 1 to $m$
|
alpar@1
|
1255 |
denote auxiliary variables, and ordinal numbers $m+1$ to $m+n$ denote
|
alpar@1
|
1256 |
structural variables. The corresponding non-zero coefficients $\xi_j$
|
alpar@1
|
1257 |
should be placed in locations \verb|val[1]|, \dots, \verb|val[len]|.
|
alpar@1
|
1258 |
The arrays \verb|ind| and \verb|val| are not changed by the routine.
|
alpar@1
|
1259 |
|
alpar@1
|
1260 |
The parameter \verb|dir| specifies direction in which the variable $x$
|
alpar@1
|
1261 |
changes on leaving the basis: $+1$ means that $x$ goes on its lower
|
alpar@1
|
1262 |
bound, so its reduced cost (dual variable) is increasing (minimization)
|
alpar@1
|
1263 |
or decreasing (maximization); $-1$ means that $x$ goes on its upper
|
alpar@1
|
1264 |
bound, so its reduced cost is decreasing (minimization) or increasing
|
alpar@1
|
1265 |
(maximization).
|
alpar@1
|
1266 |
|
alpar@1
|
1267 |
The parameter \verb|eps| is an absolute tolerance (small positive
|
alpar@1
|
1268 |
number, say, $10^{-9}$) used by the routine to skip $\xi_j$'s whose
|
alpar@1
|
1269 |
magnitude is less than \verb|eps|.
|
alpar@1
|
1270 |
|
alpar@1
|
1271 |
The routine determines which non-basic variable (among those specified
|
alpar@1
|
1272 |
in \verb|ind[1]|, \dots, \verb|ind[len]|) should enter the basis in
|
alpar@1
|
1273 |
order to keep dual feasibility, because its reduced cost reaches the
|
alpar@1
|
1274 |
(zero) bound first before this occurs for any other non-basic variables.
|
alpar@1
|
1275 |
|
alpar@1
|
1276 |
\subsubsection*{Returns}
|
alpar@1
|
1277 |
|
alpar@1
|
1278 |
The routine \verb|glp_dual_rtest| returns the index, \verb|piv|, in the
|
alpar@1
|
1279 |
arrays \verb|ind| and \verb|val| corresponding to the pivot element
|
alpar@1
|
1280 |
chosen, $1\leq$ \verb|piv| $\leq$ \verb|len|. If the adjacent basic
|
alpar@1
|
1281 |
solution is dual unbounded, and therefore the choice cannot be made,
|
alpar@1
|
1282 |
the routine returns zero.
|
alpar@1
|
1283 |
|
alpar@1
|
1284 |
\subsubsection*{Comments}
|
alpar@1
|
1285 |
|
alpar@1
|
1286 |
If the basic variable $x$ is presented in the LP problem object, the
|
alpar@1
|
1287 |
input row can be computed with the routine \verb|glp_eval_tab_row|;
|
alpar@1
|
1288 |
otherwise, it can be computed with the routine \verb|glp_transform_row|.
|
alpar@1
|
1289 |
|
alpar@1
|
1290 |
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
|
alpar@1
|
1291 |
|
alpar@1
|
1292 |
\newpage
|
alpar@1
|
1293 |
|
alpar@1
|
1294 |
\section{Post-optimal analysis routines}
|
alpar@1
|
1295 |
|
alpar@1
|
1296 |
\subsection{glp\_analyze\_bound---analyze active bound of non-basic
|
alpar@1
|
1297 |
variable}
|
alpar@1
|
1298 |
|
alpar@1
|
1299 |
\subsubsection*{Synopsis}
|
alpar@1
|
1300 |
|
alpar@1
|
1301 |
\begin{verbatim}
|
alpar@1
|
1302 |
void glp_analyze_bound(glp_prob *P, int k, double *limit1,
|
alpar@1
|
1303 |
int *var1, double *limit2, int *var2);
|
alpar@1
|
1304 |
\end{verbatim}
|
alpar@1
|
1305 |
|
alpar@1
|
1306 |
\subsubsection*{Description}
|
alpar@1
|
1307 |
|
alpar@1
|
1308 |
The routine \verb|glp_analyze_bound| analyzes the effect of varying the
|
alpar@1
|
1309 |
active bound of specified non-basic variable.
|
alpar@1
|
1310 |
|
alpar@1
|
1311 |
The non-basic variable is specified by the parameter $k$, where
|
alpar@1
|
1312 |
$1\leq k\leq m$ means auxiliary variable of corresponding row, and
|
alpar@1
|
1313 |
$m+1\leq k\leq m+n$ means structural variable (column).
|
alpar@1
|
1314 |
|
alpar@1
|
1315 |
Note that the current basic solution must be optimal, and the basis
|
alpar@1
|
1316 |
factorization must exist.
|
alpar@1
|
1317 |
|
alpar@1
|
1318 |
Results of the analysis have the following meaning.
|
alpar@1
|
1319 |
|
alpar@1
|
1320 |
\verb|value1| is the minimal value of the active bound, at which the
|
alpar@1
|
1321 |
basis still remains primal feasible and thus optimal. \verb|-DBL_MAX|
|
alpar@1
|
1322 |
means that the active bound has no lower limit.
|
alpar@1
|
1323 |
|
alpar@1
|
1324 |
\verb|var1| is the ordinal number of an auxiliary (1 to $m$) or
|
alpar@1
|
1325 |
structural ($m+1$ to $m+n$) basic variable, which reaches its bound
|
alpar@1
|
1326 |
first and thereby limits further decreasing the active bound being
|
alpar@1
|
1327 |
analyzed. if \verb|value1| = \verb|-DBL_MAX|, \verb|var1| is set to 0.
|
alpar@1
|
1328 |
|
alpar@1
|
1329 |
\verb|value2| is the maximal value of the active bound, at which the
|
alpar@1
|
1330 |
basis still remains primal feasible and thus optimal. \verb|+DBL_MAX|
|
alpar@1
|
1331 |
means that the active bound has no upper limit.
|
alpar@1
|
1332 |
|
alpar@1
|
1333 |
\verb|var2| is the ordinal number of an auxiliary (1 to $m$) or
|
alpar@1
|
1334 |
structural ($m+1$ to $m+n$) basic variable, which reaches its bound
|
alpar@1
|
1335 |
first and thereby limits further increasing the active bound being
|
alpar@1
|
1336 |
analyzed. if \verb|value2| = \verb|+DBL_MAX|, \verb|var2| is set to 0.
|
alpar@1
|
1337 |
|
alpar@1
|
1338 |
The parameters \verb|value1|, \verb|var1|, \verb|value2|, \verb|var2|
|
alpar@1
|
1339 |
can be specified as \verb|NULL|, in which case corresponding information
|
alpar@1
|
1340 |
is not stored.
|
alpar@1
|
1341 |
|
alpar@1
|
1342 |
\newpage
|
alpar@1
|
1343 |
|
alpar@1
|
1344 |
\subsection{glp\_analyze\_coef---analyze objective coefficient at basic
|
alpar@1
|
1345 |
variable}
|
alpar@1
|
1346 |
|
alpar@1
|
1347 |
\subsubsection*{Synopsis}
|
alpar@1
|
1348 |
|
alpar@1
|
1349 |
\begin{verbatim}
|
alpar@1
|
1350 |
void glp_analyze_coef(glp_prob *P, int k, double *coef1,
|
alpar@1
|
1351 |
int *var1, double *value1, double *coef2, int *var2,
|
alpar@1
|
1352 |
double *value2);
|
alpar@1
|
1353 |
\end{verbatim}
|
alpar@1
|
1354 |
|
alpar@1
|
1355 |
\subsubsection*{Description}
|
alpar@1
|
1356 |
|
alpar@1
|
1357 |
The routine \verb|glp_analyze_coef| analyzes the effect of varying the
|
alpar@1
|
1358 |
objective coefficient at specified basic variable.
|
alpar@1
|
1359 |
|
alpar@1
|
1360 |
The basic variable is specified by the parameter $k$, where
|
alpar@1
|
1361 |
$1\leq k\leq m$ means auxiliary variable of corresponding row, and
|
alpar@1
|
1362 |
$m+1\leq k\leq m+n$ means structural variable (column).
|
alpar@1
|
1363 |
|
alpar@1
|
1364 |
Note that the current basic solution must be optimal, and the basis
|
alpar@1
|
1365 |
factorization must exist.
|
alpar@1
|
1366 |
|
alpar@1
|
1367 |
Results of the analysis have the following meaning.
|
alpar@1
|
1368 |
|
alpar@1
|
1369 |
\verb|coef1| is the minimal value of the objective coefficient, at
|
alpar@1
|
1370 |
which the basis still remains dual feasible and thus optimal.
|
alpar@1
|
1371 |
\verb|-DBL_MAX| means that the objective coefficient has no lower limit.
|
alpar@1
|
1372 |
|
alpar@1
|
1373 |
\verb|var1| is the ordinal number of an auxiliary (1 to $m$) or
|
alpar@1
|
1374 |
structural ($m+1$ to $m+n$) non-basic variable, whose reduced cost
|
alpar@1
|
1375 |
reaches its zero bound first and thereby limits further decreasing the
|
alpar@1
|
1376 |
objective coefficient being analyzed. If \verb|coef1| = \verb|-DBL_MAX|,
|
alpar@1
|
1377 |
\verb|var1| is set to 0.
|
alpar@1
|
1378 |
|
alpar@1
|
1379 |
\verb|value1| is value of the basic variable being analyzed in an
|
alpar@1
|
1380 |
adjacent basis, which is defined as follows. Let the objective
|
alpar@1
|
1381 |
coefficient reaches its minimal value (\verb|coef1|) and continues
|
alpar@1
|
1382 |
decreasing. Then the reduced cost of the limiting non-basic variable
|
alpar@1
|
1383 |
(\verb|var1|) becomes dual infeasible and the current basis becomes
|
alpar@1
|
1384 |
non-optimal that forces the limiting non-basic variable to enter the
|
alpar@1
|
1385 |
basis replacing there some basic variable that leaves the basis to keep
|
alpar@1
|
1386 |
primal feasibility. Should note that on determining the adjacent basis
|
alpar@1
|
1387 |
current bounds of the basic variable being analyzed are ignored as if
|
alpar@1
|
1388 |
it were free (unbounded) variable, so it cannot leave the basis. It may
|
alpar@1
|
1389 |
happen that no dual feasible adjacent basis exists, in which case
|
alpar@1
|
1390 |
\verb|value1| is set to \verb|-DBL_MAX| or \verb|+DBL_MAX|.
|
alpar@1
|
1391 |
|
alpar@1
|
1392 |
\verb|coef2| is the maximal value of the objective coefficient, at
|
alpar@1
|
1393 |
which the basis still remains dual feasible and thus optimal.
|
alpar@1
|
1394 |
\verb|+DBL_MAX| means that the objective coefficient has no upper limit.
|
alpar@1
|
1395 |
|
alpar@1
|
1396 |
\verb|var2| is the ordinal number of an auxiliary (1 to $m$) or
|
alpar@1
|
1397 |
structural ($m+1$ to $m+n$) non-basic variable, whose reduced cost
|
alpar@1
|
1398 |
reaches its zero bound first and thereby limits further increasing the
|
alpar@1
|
1399 |
objective coefficient being analyzed. If \verb|coef2| = \verb|+DBL_MAX|,
|
alpar@1
|
1400 |
\verb|var2| is set to 0.
|
alpar@1
|
1401 |
|
alpar@1
|
1402 |
\verb|value2| is value of the basic variable being analyzed in an
|
alpar@1
|
1403 |
adjacent basis, which is defined exactly in the same way as
|
alpar@1
|
1404 |
\verb|value1| above with exception that now the objective coefficient
|
alpar@1
|
1405 |
is increasing.
|
alpar@1
|
1406 |
|
alpar@1
|
1407 |
The parameters \verb|coef1|, \verb|var1|, \verb|value1|, \verb|coef2|,
|
alpar@1
|
1408 |
\verb|var2|, \verb|value2| can be specified as \verb|NULL|, in which
|
alpar@1
|
1409 |
case corresponding information is not stored.
|
alpar@1
|
1410 |
|
alpar@1
|
1411 |
%* eof *%
|