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