1 |
/* ASCEND modelling environment |
2 |
Copyright (C) 2007 Carnegie Mellon University |
3 |
|
4 |
This program is free software; you can redistribute it and/or modify |
5 |
it under the terms of the GNU General Public License as published by |
6 |
the Free Software Foundation; either version 2, or (at your option) |
7 |
any later version. |
8 |
|
9 |
This program is distributed in the hope that it will be useful, |
10 |
but WITHOUT ANY WARRANTY; without even the implied warranty of |
11 |
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the |
12 |
GNU General Public License for more details. |
13 |
|
14 |
You should have received a copy of the GNU General Public License |
15 |
along with this program; if not, write to the Free Software |
16 |
Foundation, Inc., 59 Temple Place - Suite 330, |
17 |
Boston, MA 02111-1307, USA. |
18 |
*//** |
19 |
@file |
20 |
Connection of the IPOPT optimisation solver into ASCEND. |
21 |
|
22 |
THIS IS STILL VERY MUCH UNDER DEVELOPMENT AND INCOMPLETE. I'VE ACTUALLY |
23 |
ONLY JUST STARTED WRITING IT by starting with asc_tron.c and modifying. |
24 |
|
25 |
The IPOPT solver is documented at http://projects.coin-or.org/Ipopt/ |
26 |
*//* |
27 |
ASCEND wrapper for IPOPT originally by John Pye, Jun 2007. |
28 |
*/ |
29 |
|
30 |
#include <utilities/config.h> |
31 |
|
32 |
#ifndef ASC_WITH_IPOPT |
33 |
# error "ASC_WITH_IPOPT must be defined in order to build this." |
34 |
#endif |
35 |
|
36 |
#include <solver/solver.h> |
37 |
|
38 |
#include <system/calc.h> |
39 |
#include <system/relman.h> |
40 |
#include <system/slv_stdcalls.h> |
41 |
#include <system/block.h> |
42 |
|
43 |
#include <utilities/ascConfig.h> |
44 |
#include <utilities/ascPanic.h> |
45 |
#include <utilities/ascMalloc.h> |
46 |
#include <utilities/ascDynaLoad.h> |
47 |
#include <utilities/mem.h> |
48 |
#include <utilities/ascEnvVar.h> |
49 |
#include <general/tm_time.h> |
50 |
#include <general/env.h> |
51 |
|
52 |
#include <ipopt/IpStdCInterface.h> |
53 |
|
54 |
ASC_DLLSPEC SolverRegisterFn ipopt_register; |
55 |
|
56 |
/*------------------------------------------------------------------------------ |
57 |
DATA STRUCTURES AND FORWARD DEFS |
58 |
*/ |
59 |
|
60 |
/** |
61 |
Documentation of solver options for IPOPT is at |
62 |
http://www.coin-or.org/Ipopt/documentation/node1.html |
63 |
*/ |
64 |
enum{ |
65 |
IPOPT_PARAM_TOL |
66 |
,IPOPT_PARAM_MAX_ITER |
67 |
,IPOPT_PARAM_SAFECALC |
68 |
,IPOPT_PARAM_MU_STRATEGY |
69 |
,IPOPT_PARAMS |
70 |
}; |
71 |
|
72 |
#define SYS(s) ((IpoptSystem *)(s)) |
73 |
|
74 |
struct IpoptSystemStruct{ |
75 |
|
76 |
/* |
77 |
Problem definition |
78 |
*/ |
79 |
slv_system_t slv; /* slv_system_t back-link */ |
80 |
struct rel_relation *obj; /* Objective function: NULL = none */ |
81 |
struct rel_relation *old_obj;/* Objective function: NULL = none */ |
82 |
struct var_variable **vlist; /* Variable list (NULL terminated) */ |
83 |
struct rel_relation **rlist; /* Relation list (NULL terminated) */ |
84 |
|
85 |
/* |
86 |
Solver information |
87 |
*/ |
88 |
int32 presolved; /* ? Has the system been presolved */ |
89 |
int32 resolve; /* ? Has the system been resolved */ |
90 |
slv_parameters_t p; /* Parameters */ |
91 |
slv_status_t s; /* Status (as of iteration end) */ |
92 |
|
93 |
int32 cap; /* Order of matrix/vectors */ |
94 |
int32 rank; /* Symbolic rank of problem */ |
95 |
int32 vused; /* Free and incident variables */ |
96 |
int32 vtot; /* length of varlist */ |
97 |
int32 rused; /* Included relations */ |
98 |
int32 rtot; /* length of rellist */ |
99 |
double clock; /* CPU time */ |
100 |
|
101 |
int32 calc_ok; |
102 |
|
103 |
void *parm_array[IPOPT_PARAMS]; |
104 |
struct slv_parameter pa[IPOPT_PARAMS]; |
105 |
|
106 |
/* |
107 |
IPOPT DATA |
108 |
*/ |
109 |
Index n; /* number of variables */ |
110 |
Index m; /* number of constraints */ |
111 |
|
112 |
int nnzJ; /* number of non zeros in the jacobian of the constraints */ |
113 |
|
114 |
Number* x_L; /* lower bounds on x */ |
115 |
Number* x_U; /* upper bounds on x */ |
116 |
Number* g_L; /* lower bounds on g */ |
117 |
Number* g_U; /* upper bounds on g */ |
118 |
|
119 |
IpoptProblem nlp; /* IpoptProblem */ |
120 |
|
121 |
enum ApplicationReturnStatus status; /* Solve return code */ |
122 |
Number* x; /* starting point and solution vector */ |
123 |
Number* mult_x_L; /* lower bound multipliers at the solution */ |
124 |
Number* mult_x_U; /* upper bound multipliers at the solution */ |
125 |
Index i; /* generic counter */ |
126 |
}; |
127 |
|
128 |
typedef struct IpoptSystemStruct IpoptSystem; |
129 |
|
130 |
static int ipopt_get_default_parameters(slv_system_t server, SlvClientToken asys, slv_parameters_t *parameters); |
131 |
|
132 |
/*------------------------------------------------------------------------------ |
133 |
SYSTEM SETUP/DESTROY AND SOLVER ELIGIBILITY |
134 |
*/ |
135 |
|
136 |
static SlvClientToken ipopt_create(slv_system_t server, int32*statusindex){ |
137 |
IpoptSystem *sys; |
138 |
|
139 |
sys = ASC_NEW_CLEAR(IpoptSystem); |
140 |
if(sys==NULL){ |
141 |
*statusindex = 1; |
142 |
return sys; |
143 |
} |
144 |
|
145 |
sys->p.parms = sys->pa; |
146 |
sys->p.dynamic_parms = 0; |
147 |
ipopt_get_default_parameters(server,(SlvClientToken)sys,&(sys->p)); |
148 |
sys->p.whose = (*statusindex); |
149 |
|
150 |
sys->presolved = 0; |
151 |
sys->resolve = 0; |
152 |
|
153 |
sys->n = -1; |
154 |
sys->m = -1; |
155 |
|
156 |
sys->s.ok = TRUE; |
157 |
sys->s.calc_ok = TRUE; |
158 |
sys->s.costsize = 0; |
159 |
sys->s.cost = NULL; /*redundant, but sanity preserving */ |
160 |
|
161 |
sys->vlist = slv_get_solvers_var_list(server); |
162 |
sys->rlist = slv_get_solvers_rel_list(server); |
163 |
sys->obj = slv_get_obj_relation(server); |
164 |
if(sys->vlist == NULL) { |
165 |
ASC_FREE(sys); |
166 |
ERROR_REPORTER_HERE(ASC_PROG_ERR,"IPOPT called with no variables."); |
167 |
*statusindex = -2; |
168 |
return NULL; |
169 |
} |
170 |
if(sys->rlist == NULL && sys->obj == NULL) { |
171 |
ASC_FREE(sys); |
172 |
ERROR_REPORTER_HERE(ASC_PROG_ERR,"IPOPT called with no relations or objective."); |
173 |
*statusindex = -1; |
174 |
return NULL; |
175 |
} |
176 |
|
177 |
/* do nothing with the objective list, pars, bounds, extrels, etc etc */ |
178 |
|
179 |
slv_check_var_initialization(server); |
180 |
*statusindex = 0; |
181 |
return((SlvClientToken)sys); |
182 |
} |
183 |
|
184 |
static int32 ipopt_destroy(slv_system_t server, SlvClientToken asys){ |
185 |
UNUSED_PARAMETER(server); |
186 |
ERROR_REPORTER_HERE(ASC_PROG_ERR,"Not implemented"); |
187 |
return 1; |
188 |
} |
189 |
|
190 |
static int32 ipopt_eligible_solver(slv_system_t server){ |
191 |
UNUSED_PARAMETER(server); |
192 |
ERROR_REPORTER_HERE(ASC_PROG_ERR,"Not implemented"); |
193 |
return 1; |
194 |
} |
195 |
|
196 |
/*------------------------------------------------------------------------------ |
197 |
SOLVER PARAMETERS |
198 |
*/ |
199 |
|
200 |
static |
201 |
int32 ipopt_get_default_parameters(slv_system_t server, SlvClientToken asys |
202 |
,slv_parameters_t *parameters |
203 |
){ |
204 |
IpoptSystem *sys = NULL; |
205 |
struct slv_parameter *new_parms = NULL; |
206 |
int32 make_macros = 0; |
207 |
|
208 |
if(server != NULL && asys != NULL) { |
209 |
sys = SYS(asys); |
210 |
make_macros = 1; |
211 |
} |
212 |
|
213 |
if(parameters->parms == NULL) { |
214 |
new_parms = ASC_NEW_ARRAY_OR_NULL(struct slv_parameter,IPOPT_PARAMS); |
215 |
if(new_parms == NULL) { |
216 |
return -1; |
217 |
} |
218 |
parameters->parms = new_parms; |
219 |
parameters->dynamic_parms = 1; |
220 |
} |
221 |
|
222 |
parameters->num_parms = 0; |
223 |
|
224 |
slv_param_int(parameters,IPOPT_PARAM_MAX_ITER |
225 |
,(SlvParameterInitInt){{"max_iter" |
226 |
,"Maximum number of iterations",2 |
227 |
,"The algorithm terminates with an error message if the number of iterations exceeded this number." |
228 |
}, 3000, 0, 100000000} |
229 |
); |
230 |
|
231 |
slv_param_real(parameters,IPOPT_PARAM_TOL |
232 |
,(SlvParameterInitReal){{"tol" |
233 |
,"Desired convergence tolerance (relative)",2 |
234 |
,"Determines the convergence tolerance for the algorithm. The algorithm" |
235 |
" terminates successfully, if the (scaled) NLP error becomes smaller" |
236 |
" than this value, and if the (absolute) criteria according to " |
237 |
" 'dual_inf_tol', 'primal_inf_tol', and 'cmpl_inf_tol' are met. (This" |
238 |
" is epsilon_tol in Eqn. (6) in implementation paper). See also " |
239 |
" 'acceptable_tol' as a second termination criterion. Note, some other" |
240 |
" algorithmic features also use this quantity to determine thresholds" |
241 |
" etc." |
242 |
}, 1e-8, 0, 1e20} |
243 |
); |
244 |
|
245 |
slv_param_char(parameters,IPOPT_PARAM_MU_STRATEGY |
246 |
,(SlvParameterInitChar){{"mu_strategy" |
247 |
,"Update strategy for barrier parameter",6 |
248 |
,"Determines which barrier parameter update strategy is to be used." |
249 |
" 'MONOTONE' is the monotone (Fiacco-McCormick) strategy;" |
250 |
" 'ADAPTIVE' is the adaptive update strategy." |
251 |
}, "MONOTONE"}, (char *[]){ |
252 |
"MONOTONE","ADAPTIVE",NULL |
253 |
} |
254 |
); |
255 |
|
256 |
asc_assert(parameters->num_parms==IPOPT_PARAMS); |
257 |
|
258 |
return 1; |
259 |
} |
260 |
|
261 |
static void ipopt_get_parameters(slv_system_t server, SlvClientToken asys |
262 |
, slv_parameters_t *parameters |
263 |
){ |
264 |
IpoptSystem *sys; |
265 |
UNUSED_PARAMETER(server); |
266 |
|
267 |
sys = SYS(asys); |
268 |
//if(check_system(sys)) return; |
269 |
mem_copy_cast(&(sys->p),parameters,sizeof(slv_parameters_t)); |
270 |
} |
271 |
|
272 |
|
273 |
static void ipopt_set_parameters(slv_system_t server, SlvClientToken asys |
274 |
,slv_parameters_t *parameters |
275 |
){ |
276 |
IpoptSystem *sys; |
277 |
UNUSED_PARAMETER(server); |
278 |
|
279 |
sys = SYS(asys); |
280 |
//if (check_system(sys)) return; |
281 |
mem_copy_cast(parameters,&(sys->p),sizeof(slv_parameters_t)); |
282 |
} |
283 |
|
284 |
/*------------------------------------------------------------------------------ |
285 |
EVALUATION FUNCTIONS |
286 |
*/ |
287 |
|
288 |
/** |
289 |
update the model with new 'x' vector. |
290 |
@return 0 on success. |
291 |
*/ |
292 |
int ipopt_update_model(IpoptSystem *sys, const double *x){ |
293 |
ERROR_REPORTER_HERE(ASC_PROG_ERR,"Not implemented"); |
294 |
return 1; /* error */ |
295 |
} |
296 |
|
297 |
/** Function to evaluate the objective function f(x). |
298 |
@return 1 on success, 0 on failure |
299 |
|
300 |
@param n (in), the number of variables in the problem (dimension of 'x'). |
301 |
@param x (in), the values for the primal variables, 'x' , at which 'f(x)' is to be evaluated. |
302 |
@param new_x (in), false if any evaluation method was previously called with the same values in 'x', true otherwise. |
303 |
@param obj_value (out) the value of the objective function ('f(x)'). |
304 |
*/ |
305 |
Bool ipopt_eval_f(Index n, Number *x, Bool new_x, Number *obj_value, void *user_data){ |
306 |
IpoptSystem *sys; |
307 |
sys = SYS(user_data); |
308 |
int res; |
309 |
|
310 |
asc_assert(n==sys->n); |
311 |
asc_assert(sys->obj!=NULL); |
312 |
|
313 |
if(new_x){ |
314 |
res = ipopt_update_model(sys,x); |
315 |
if(res)return 0; /* fail model update */ |
316 |
} |
317 |
|
318 |
sys->calc_ok = TRUE; |
319 |
|
320 |
*obj_value = relman_eval(sys->obj,&(sys->calc_ok),SLV_PARAM_BOOL(&(sys->p),IPOPT_PARAM_SAFECALC)); |
321 |
|
322 |
return sys->calc_ok; |
323 |
} |
324 |
|
325 |
/** |
326 |
@return 1 on success |
327 |
*/ |
328 |
Bool ipopt_eval_grad_f(Index n, Number* x, Bool new_x, Number* grad_f, void *user_data){ |
329 |
IpoptSystem *sys; |
330 |
sys = SYS(user_data); |
331 |
int j, res, len; |
332 |
int count; |
333 |
double *derivatives; |
334 |
int *variables; |
335 |
static var_filter_t vfilter = { |
336 |
VAR_ACTIVE | VAR_INCIDENT | VAR_SVAR |
337 |
,VAR_ACTIVE | VAR_INCIDENT | VAR_SVAR | VAR_FIXED |
338 |
}; |
339 |
|
340 |
asc_assert(n==sys->n); |
341 |
asc_assert(sys->obj); |
342 |
|
343 |
if(new_x){ |
344 |
res = ipopt_update_model(sys,x); |
345 |
if(res)return 0; /* fail model update */ |
346 |
} |
347 |
|
348 |
|
349 |
/* evaluate grad_f(x) somehow */ |
350 |
for(j=0; j<n; ++j){ |
351 |
grad_f[j] = 0; |
352 |
} |
353 |
|
354 |
len = rel_n_incidences(sys->obj); |
355 |
variables = ASC_NEW_ARRAY(int,len); |
356 |
derivatives = ASC_NEW_ARRAY(double,len); |
357 |
|
358 |
relman_diff2( |
359 |
sys->obj,&vfilter,derivatives,variables |
360 |
, &count,SLV_PARAM_BOOL(&(sys->p),IPOPT_PARAM_SAFECALC) |
361 |
); |
362 |
|
363 |
for(j=0; j<len; ++j){ |
364 |
grad_f[variables[j]] = derivatives[j]; |
365 |
} |
366 |
|
367 |
ASC_FREE(variables); |
368 |
ASC_FREE(derivatives); |
369 |
|
370 |
return 1; /* success, presumably */ |
371 |
} |
372 |
|
373 |
Bool ipopt_eval_g(Index n, Number* x, Bool new_x, Index m, Number *g, void *user_data){ |
374 |
IpoptSystem *sys; |
375 |
sys = SYS(user_data); |
376 |
int i, res; |
377 |
|
378 |
asc_assert(n==sys->n); |
379 |
asc_assert(m==sys->m); |
380 |
|
381 |
if(new_x){ |
382 |
res = ipopt_update_model(sys,x); |
383 |
if(res)return 0; /* fail model update */ |
384 |
} |
385 |
|
386 |
for(i=0; i<m; ++i){ |
387 |
g[i] = 0; |
388 |
} |
389 |
|
390 |
return 0; /* fail: not yet implemented */ |
391 |
} |
392 |
|
393 |
Bool ipopt_eval_jac_g(Index n, Number* x, Bool new_x, Index m |
394 |
, Index nele_jac, Index* iRow, Index *jCol, Number* values |
395 |
, void *user_data |
396 |
){ |
397 |
IpoptSystem *sys; |
398 |
sys = SYS(user_data); |
399 |
int i,j,res; |
400 |
|
401 |
asc_assert(n==sys->n); |
402 |
asc_assert(nele_jac==sys->nnzJ); |
403 |
asc_assert(m==sys->m); |
404 |
|
405 |
if(new_x){ |
406 |
res = ipopt_update_model(sys,x); |
407 |
if(res)return 0; /* fail model update */ |
408 |
} |
409 |
|
410 |
/* loop through the constraints */ |
411 |
for(i=0; i<m; ++i){ |
412 |
/* get derivatives for constraint i */ |
413 |
/* insert the derivatives into the matrix in row i, columns j */ |
414 |
} |
415 |
|
416 |
return 0; /* fail: not yet implemented */ |
417 |
} |
418 |
|
419 |
Bool ipopt_eval_h(Index n, Number* x, Bool new_x |
420 |
, Number obj_factor, Index m, Number* lambda |
421 |
, Bool new_lambda, Index nele_hess, Index* iRow |
422 |
, Index* jCol, Number* values |
423 |
){ |
424 |
/* not sure about this one yet: do all the 2nd deriv things work for this? */ |
425 |
return 0; /* fail: not yet implemented */ |
426 |
} |
427 |
|
428 |
/*------------------------------------------------------------------------------ |
429 |
SOLVE ROUTINES |
430 |
*/ |
431 |
|
432 |
static int ipopt_solve(slv_system_t server, SlvClientToken asys){ |
433 |
IpoptSystem *sys; |
434 |
UNUSED_PARAMETER(server); |
435 |
int i; |
436 |
|
437 |
sys = SYS(asys); |
438 |
|
439 |
/* set the number of variables and allocate space for the bounds */ |
440 |
sys->x_L = ASC_NEW_ARRAY(Number,sys->n); |
441 |
sys->x_U = ASC_NEW_ARRAY(Number,sys->n); |
442 |
|
443 |
/* set the values for the variable bounds */ |
444 |
|
445 |
#if 0 |
446 |
// sample code for the C interface |
447 |
|
448 |
int main(){ |
449 |
|
450 |
/* set the number of constraints and allocate space for the bounds */ |
451 |
m=2; |
452 |
g_L = (Number*)malloc(sizeof(Number)*m); |
453 |
g_U = (Number*)malloc(sizeof(Number)*m); |
454 |
/* set the values of the constraint bounds */ |
455 |
g_L[0] = 25; g_U[0] = 2e19; |
456 |
g_L[1] = 40; g_U[1] = 40; |
457 |
|
458 |
/* create the IpoptProblem */ |
459 |
nlp = CreateIpoptProblem(n, x_L, x_U, m, g_L, g_U, 8, 10, 0, |
460 |
&eval_f, &eval_g, &eval_grad_f, |
461 |
&eval_jac_g, &eval_h); |
462 |
|
463 |
/* We can free the memory now - the values for the bounds have been |
464 |
copied internally in CreateIpoptProblem */ |
465 |
free(x_L); |
466 |
free(x_U); |
467 |
free(g_L); |
468 |
free(g_U); |
469 |
|
470 |
/* set some options */ |
471 |
AddIpoptNumOption(nlp, "tol", 1e-9); |
472 |
AddIpoptStrOption(nlp, "mu_strategy", "adaptive"); |
473 |
|
474 |
/* allocate space for the initial point and set the values */ |
475 |
x = (Number*)malloc(sizeof(Number)*n); |
476 |
x[0] = 1.0; |
477 |
x[1] = 5.0; |
478 |
x[2] = 5.0; |
479 |
x[3] = 1.0; |
480 |
|
481 |
/* allocate space to store the bound multipliers at the solution */ |
482 |
mult_x_L = (Number*)malloc(sizeof(Number)*n); |
483 |
mult_x_U = (Number*)malloc(sizeof(Number)*n); |
484 |
|
485 |
/* solve the problem */ |
486 |
status = IpoptSolve(nlp, x, NULL, &obj, NULL, mult_x_L, mult_x_U, NULL); |
487 |
|
488 |
if (status == Solve_Succeeded) { |
489 |
printf("\n\nSolution of the primal variables, x\n"); |
490 |
for (i=0; i<n; i++) { |
491 |
printf("x[%d] = %e\n", i, x[i]); |
492 |
} |
493 |
|
494 |
printf("\n\nSolution of the bound multipliers, z_L and z_U\n"); |
495 |
for (i=0; i<n; i++) { |
496 |
printf("z_L[%d] = %e\n", i, mult_x_L[i]); |
497 |
} |
498 |
for (i=0; i<n; i++) { |
499 |
printf("z_U[%d] = %e\n", i, mult_x_U[i]); |
500 |
} |
501 |
|
502 |
printf("\n\nObjective value\n"); |
503 |
printf("f(x*) = %e\n", obj); |
504 |
} |
505 |
|
506 |
/* free allocated memory */ |
507 |
FreeIpoptProblem(nlp); |
508 |
free(x); |
509 |
free(mult_x_L); |
510 |
free(mult_x_U); |
511 |
|
512 |
return 0; |
513 |
} |
514 |
#endif |
515 |
|
516 |
|
517 |
ERROR_REPORTER_HERE(ASC_PROG_ERR,"Not implemented"); |
518 |
return 1; |
519 |
} |
520 |
|
521 |
static int ipopt_presolve(slv_system_t server, SlvClientToken asys){ |
522 |
IpoptSystem *sys; |
523 |
struct var_variable **vp; |
524 |
struct rel_relation **rp; |
525 |
int cap, ind; |
526 |
int matrix_creation_needed = 1; |
527 |
int *cntvect, temp; |
528 |
|
529 |
CONSOLE_DEBUG("PRESOLVE"); |
530 |
|
531 |
sys = SYS(asys); |
532 |
//iteration_begins(sys); |
533 |
//check_system(sys); |
534 |
|
535 |
asc_assert(sys->vlist && sys->rlist); |
536 |
|
537 |
sys->obj = slv_get_obj_relation(server); /*may have changed objective*/ |
538 |
if(!sys->obj){ |
539 |
ERROR_REPORTER_HERE(ASC_PROG_ERR,"No objective function was specified"); |
540 |
return -3; |
541 |
} |
542 |
|
543 |
#if 0 |
544 |
if(sys->presolved > 0) { /* system has been presolved before */ |
545 |
if(!conopt_dof_changed(sys) /*no changes in fixed or included flags*/ |
546 |
&& sys->p.partition == sys->J.old_partition |
547 |
&& sys->obj == sys->old_obj |
548 |
){ |
549 |
matrix_creation_needed = 0; |
550 |
CONOPT_CONSOLE_DEBUG("YOU JUST AVOIDED MATRIX DESTRUCTION/CREATION"); |
551 |
} |
552 |
} |
553 |
#endif |
554 |
|
555 |
#if 0 |
556 |
// check all this... |
557 |
|
558 |
/* set all relations as being 'unsatisfied' to start with... */ |
559 |
rp=sys->rlist; |
560 |
for( ind = 0; ind < sys->rtot; ++ind ) { |
561 |
rel_set_satisfied(rp[ind],FALSE); |
562 |
} |
563 |
|
564 |
if( matrix_creation_needed ) { |
565 |
|
566 |
cap = slv_get_num_solvers_rels(server); |
567 |
sys->cap = slv_get_num_solvers_vars(server); |
568 |
sys->cap = MAX(sys->cap,cap); |
569 |
vp=sys->vlist; |
570 |
for( ind = 0; ind < sys->vtot; ++ind ) { |
571 |
var_set_in_block(vp[ind],FALSE); |
572 |
} |
573 |
rp=sys->rlist; |
574 |
for( ind = 0; ind < sys->rtot; ++ind ) { |
575 |
rel_set_in_block(rp[ind],FALSE); |
576 |
/* rel_set_satisfied(rp[ind],FALSE); */ |
577 |
} |
578 |
|
579 |
sys->presolved = 1; /* full presolve recognized here */ |
580 |
sys->resolve = 0; /* initialize resolve flag */ |
581 |
|
582 |
/* provide nnz for jacobian matrix of constraints */ |
583 |
|
584 |
/* provide sparsity structure for jacobian */ |
585 |
|
586 |
/* provide nnz for hessian matrix */ |
587 |
|
588 |
/* provide sparsity structure for hessian matrix */ |
589 |
|
590 |
sys->J.old_partition = sys->p.partition; |
591 |
sys->old_obj = sys->obj; |
592 |
|
593 |
slv_sort_rels_and_vars(server,&(sys->con.m),&(sys->con.n)); |
594 |
CONOPT_CONSOLE_DEBUG("FOUND %d CONSTRAINTS AND %d VARS",sys->con.m,sys->con.n); |
595 |
if (sys->obj != NULL) { |
596 |
CONOPT_CONSOLE_DEBUG("ADDING OBJECT AS A ROW"); |
597 |
sys->con.m++; /* treat objective as a row */ |
598 |
} |
599 |
|
600 |
cntvect = ASC_NEW_ARRAY(int,COIDEF_Size()); |
601 |
COIDEF_Ini(cntvect); |
602 |
sys->con.cntvect = cntvect; |
603 |
CONOPT_CONSOLE_DEBUG("NUMBER OF CONSTRAINTS = %d",sys->con.m); |
604 |
COIDEF_NumVar(cntvect, &(sys->con.n)); |
605 |
COIDEF_NumCon(cntvect, &(sys->con.m)); |
606 |
sys->con.nz = num_jacobian_nonzeros(sys, &(sys->con.maxrow)); |
607 |
COIDEF_NumNZ(cntvect, &(sys->con.nz)); |
608 |
COIDEF_NumNlNz(cntvect, &(sys->con.nz)); |
609 |
|
610 |
sys->con.base = 0; |
611 |
COIDEF_Base(cntvect,&(sys->con.base)); |
612 |
COIDEF_ErrLim(cntvect, &(DOMLIM)); |
613 |
COIDEF_ItLim(cntvect, &(ITER_LIMIT)); |
614 |
|
615 |
if(sys->obj!=NULL){ |
616 |
sys->con.optdir = relman_obj_direction(sys->obj); |
617 |
sys->con.objcon = sys->con.m - 1; /* objective will be last row */ |
618 |
CONOPT_CONSOLE_DEBUG("SETTING OBJECTIVE CONSTRAINT TO BE %d",sys->con.objcon); |
619 |
}else{ |
620 |
sys->con.optdir = 0; |
621 |
sys->con.objcon = 0; |
622 |
} |
623 |
COIDEF_OptDir(cntvect, &(sys->con.optdir)); |
624 |
COIDEF_ObjCon(cntvect, &(sys->con.objcon)); |
625 |
|
626 |
temp = 0; |
627 |
COIDEF_StdOut(cntvect, &temp); |
628 |
|
629 |
COIDEF_ReadMatrix(cntvect, &conopt_readmatrix); |
630 |
COIDEF_FDEval(cntvect, &conopt_fdeval); |
631 |
COIDEF_Option(cntvect, &conopt_option); |
632 |
COIDEF_Solution(cntvect, &conopt_solution); |
633 |
COIDEF_Status(cntvect, &conopt_status); |
634 |
COIDEF_Message(cntvect, &asc_conopt_message); |
635 |
COIDEF_ErrMsg(cntvect, &conopt_errmsg); |
636 |
COIDEF_Progress(cntvect, &asc_conopt_progress); |
637 |
|
638 |
int debugfv = 1; |
639 |
COIDEF_DebugFV(cntvect, &debugfv); |
640 |
|
641 |
destroy_vectors(sys); |
642 |
destroy_matrices(sys); |
643 |
create_matrices(server,sys); |
644 |
create_vectors(sys); |
645 |
|
646 |
sys->s.block.current_reordered_block = -2; |
647 |
} |
648 |
|
649 |
/* Reset status */ |
650 |
sys->con.optimized = 0; |
651 |
sys->s.iteration = 0; |
652 |
sys->s.cpu_elapsed = 0.0; |
653 |
sys->s.converged = sys->s.diverged = sys->s.inconsistent = FALSE; |
654 |
sys->s.block.previous_total_size = 0; |
655 |
sys->s.costsize = 1+sys->s.block.number_of; |
656 |
|
657 |
if( matrix_creation_needed ) { |
658 |
destroy_array(sys->s.cost); |
659 |
sys->s.cost = create_zero_array(sys->s.costsize,struct slv_block_cost); |
660 |
for( ind = 0; ind < sys->s.costsize; ++ind ) { |
661 |
sys->s.cost[ind].reorder_method = -1; |
662 |
} |
663 |
} else { |
664 |
reset_cost(sys->s.cost,sys->s.costsize); |
665 |
} |
666 |
|
667 |
/* set to go to first unconverged block */ |
668 |
sys->s.block.current_block = -1; |
669 |
sys->s.block.current_size = 0; |
670 |
sys->s.calc_ok = TRUE; |
671 |
sys->s.block.iteration = 0; |
672 |
sys->objective = MAXDOUBLE/2000.0; |
673 |
|
674 |
update_status(sys); |
675 |
iteration_ends(sys); |
676 |
sys->s.cost[sys->s.block.number_of].time=sys->s.cpu_elapsed; |
677 |
#endif |
678 |
|
679 |
return 0; |
680 |
} |
681 |
|
682 |
static int ipopt_iterate(slv_system_t server, SlvClientToken asys){ |
683 |
UNUSED_PARAMETER(server); |
684 |
ERROR_REPORTER_HERE(ASC_PROG_ERR,"Not implemented"); |
685 |
return 1; |
686 |
} |
687 |
|
688 |
static int ipopt_resolve(slv_system_t server, SlvClientToken asys){ |
689 |
UNUSED_PARAMETER(server); |
690 |
|
691 |
/* if implementing this, use the 'warm start' thing in IPOPT */ |
692 |
|
693 |
/* provide initial values of the 'multipliers' */ |
694 |
|
695 |
ERROR_REPORTER_HERE(ASC_PROG_ERR,"Not implemented"); |
696 |
return 1; |
697 |
} |
698 |
|
699 |
static const SlvFunctionsT ipopt_internals = { |
700 |
67 |
701 |
,"IPOPT" |
702 |
,ipopt_create |
703 |
,ipopt_destroy |
704 |
,ipopt_eligible_solver |
705 |
,ipopt_get_default_parameters |
706 |
,ipopt_get_parameters |
707 |
,ipopt_set_parameters |
708 |
,NULL |
709 |
,ipopt_solve |
710 |
,ipopt_presolve |
711 |
,ipopt_iterate |
712 |
,ipopt_resolve |
713 |
,NULL |
714 |
,NULL |
715 |
,NULL |
716 |
}; |
717 |
|
718 |
int ipopt_register(void){ |
719 |
return solver_register(&ipopt_internals); |
720 |
} |