/[ascend]/trunk/solvers/ipopt/asc_ipopt.c
ViewVC logotype

Contents of /trunk/solvers/ipopt/asc_ipopt.c

Parent Directory Parent Directory | Revision Log Revision Log


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

john.pye@anu.edu.au
ViewVC Help
Powered by ViewVC 1.1.22