/home/coin/SVN-release/OS-2.0.1/OS/src/OSSolverInterfaces/OSIpoptSolver.cpp

Go to the documentation of this file.
00001 /* $Id: OSIpoptSolver.cpp 2831 2009-07-21 07:34:44Z kmartin $ */
00020 //#define DEBUG
00021 
00022 #include "OSIpoptSolver.h"
00023 #include "OSDataStructures.h"
00024 #include "OSParameters.h" 
00025 #include "OSMathUtil.h"
00026 #include "CoinFinite.hpp"
00027 
00028 
00029 using std::cout; 
00030 using std::endl; 
00031 using std::ostringstream;
00032 //using namespace Ipopt;
00033 
00034 
00035 IpoptSolver::IpoptSolver() {
00036         osrlwriter = new OSrLWriter();
00037         osresult = new OSResult();
00038         m_osilreader = NULL;
00039         m_osolreader = NULL;
00040         ipoptErrorMsg = new std::string("");
00041 } 
00042 
00043 IpoptSolver::~IpoptSolver() {
00044         #ifdef DEBUG
00045         cout << "inside IpoptSolver destructor" << endl;
00046         #endif
00047         if(m_osilreader != NULL) delete m_osilreader;
00048         m_osilreader = NULL;
00049         if(m_osolreader != NULL) delete m_osolreader;
00050         m_osolreader = NULL;
00051         delete osresult;
00052         osresult = NULL;
00053         delete osrlwriter;
00054         osrlwriter = NULL;
00055         //delete osinstance;
00056         //osinstance = NULL;
00057         delete ipoptErrorMsg;
00058         #ifdef DEBUG
00059         cout << "leaving IpoptSolver destructor" << endl;
00060         #endif
00061 }
00062 
00063 // returns the size of the problem
00064 bool IpoptProblem::get_nlp_info(Index& n, Index& m, Index& nnz_jac_g,
00065                              Index& nnz_h_lag, IndexStyleEnum& index_style)
00066 {
00067         try{
00068                 if(osinstance->getObjectiveNumber() <= 0) throw ErrorClass("Ipopt NEEDS AN OBJECTIVE FUNCTION");     
00069                 if( (osinstance->getNumberOfIntegerVariables() + osinstance->getNumberOfBinaryVariables()) > 0 )  
00070                         throw ErrorClass("Ipopt does not solve integer programs -- please try Bonmin or Couenne");
00071                 // number of variables
00072                 n = osinstance->getVariableNumber();
00073                 // number of constraints
00074                 m = osinstance->getConstraintNumber();
00075         #ifdef DEBUG
00076                 cout << "number variables  !!!!!!!!!!!!!!!!!!!!!!!!!!!" << n << endl;
00077                 cout << "number constraints  !!!!!!!!!!!!!!!!!!!!!!!!!!!" << m << endl;
00078         #endif
00079                 try{
00080                         osinstance->initForAlgDiff( );
00081                 }
00082                 catch(const ErrorClass& eclass){
00083         #ifdef DEBUG
00084                         cout << "error in OSIpoptSolver, line 78:\n" << eclass.errormsg << endl;
00085         #endif
00086                         *ipoptErrorMsg = eclass.errormsg;
00087                         throw;  
00088                 }       
00089                 // use the OS Expression tree for function evaluations instead of CppAD
00090                 osinstance->bUseExpTreeForFunEval = true;
00091                 //std::cout << "Call sparse jacobian" << std::endl;
00092                 SparseJacobianMatrix *sparseJacobian = NULL;
00093                 try{
00094                         sparseJacobian = osinstance->getJacobianSparsityPattern();
00095                 }
00096                 catch(const ErrorClass& eclass){
00097         #ifdef DEBUG
00098                         cout << "error in OSIpoptSolver, line 91:\n" << eclass.errormsg << endl;
00099         #endif
00100                         *ipoptErrorMsg = eclass.errormsg;
00101                         throw;  
00102                 }
00103                 //std::cout << "Done calling sparse jacobian" << std::endl;
00104                 nnz_jac_g = sparseJacobian->valueSize;
00105         #ifdef DEBUG
00106                 cout << "nnz_jac_g  !!!!!!!!!!!!!!!!!!!!!!!!!!!" << nnz_jac_g << endl;  
00107         #endif
00108                 // nonzeros in upper hessian
00109                 
00110                 if( (osinstance->getNumberOfNonlinearExpressions() == 0) && (osinstance->getNumberOfQuadraticTerms() == 0) ) {
00111         #ifdef DEBUG
00112                         cout << "This is a linear program"  << endl;
00113         #endif
00114                         nnz_h_lag = 0;
00115                 }
00116                 else{
00117                         //std::cout << "Get Lagrangain Hessian Sparsity Pattern " << std::endl;
00118                         SparseHessianMatrix *sparseHessian = osinstance->getLagrangianHessianSparsityPattern();
00119                         //std::cout << "Done Getting Lagrangain Hessian Sparsity Pattern " << std::endl;
00120                         nnz_h_lag = sparseHessian->hessDimension;
00121                 }
00122         #ifdef DEBUG
00123                 cout << "print nnz_h_lag (OSIpoptSolver.cpp)" << endl;  
00124                 cout << "nnz_h_lag  !!!!!!!!!!!!!!!!!!!!!!!!!!!" << nnz_h_lag << endl;  
00125                 cout << "set index_style (OSIpoptSolver.cpp)" << endl;  
00126         #endif
00127                 // use the C style indexing (0-based)
00128                 index_style = TNLP::C_STYLE;
00129         #ifdef DEBUG
00130                 cout << "return from get_nlp_info (OSIpoptSolver.cpp)" << nnz_h_lag << endl;    
00131         #endif
00132           
00134 
00135           return true;
00136         }
00137         catch(const ErrorClass& eclass){
00138 
00139                 *ipoptErrorMsg = eclass.errormsg;
00140                 throw;  
00141         }
00142         
00143 }//get_nlp_info
00144 
00145 
00146 bool  IpoptProblem::get_bounds_info(Index n, Number* x_l, Number* x_u,
00147                                 Index m, Number* g_l, Number* g_u){
00148         int i; 
00149         double * mdVarLB = osinstance->getVariableLowerBounds();
00150         //std::cout << "GET BOUNDS INFORMATION FOR IPOPT !!!!!!!!!!!!!!!!!" << std::endl;
00151         // variables upper bounds
00152         double * mdVarUB = osinstance->getVariableUpperBounds();
00153 
00154         for(i = 0; i < n; i++){
00155                 x_l[ i] = mdVarLB[ i];
00156                 x_u[ i] = mdVarUB[ i];
00157                 //cout << "x_l !!!!!!!!!!!!!!!!!!!!!!!!!!!" << x_l[i] << endl;
00158                 //cout << "x_u !!!!!!!!!!!!!!!!!!!!!!!!!!!" << x_u[i] << endl;
00159         }
00160         // Ipopt interprets any number greater than nlp_upper_bound_inf as
00161         // infinity. The default value of nlp_upper_bound_inf and nlp_lower_bound_inf
00162         // is 1e19 and can be changed through ipopt options.
00163         // e.g. g_u[0] = 2e19;
00164 
00165         //constraint lower bounds
00166         double * mdConLB = osinstance->getConstraintLowerBounds();
00167         //constraint upper bounds
00168         double * mdConUB = osinstance->getConstraintUpperBounds();
00169         
00170         for(int i = 0; i < m; i++){
00171                 g_l[ i] = mdConLB[ i];
00172                 g_u[ i] = mdConUB[ i];
00173                 //cout << "lower !!!!!!!!!!!!!!!!!!!!!!!!!!!" << g_l[i] << endl;
00174                 //cout << "upper !!!!!!!!!!!!!!!!!!!!!!!!!!!" << g_u[i] << endl;
00175         }  
00176         return true;
00177 }//get_bounds_info
00178 
00179 
00180 // returns the initial point for the problem
00181 bool IpoptProblem::get_starting_point(Index n, bool init_x, Number* x,
00182      bool init_z, Number* z_L, Number* z_U, Index m, bool init_lambda,
00183      Number* lambda) {
00184         // Here, we assume we only have starting values for x, if you code
00185         // your own NLP, you can provide starting values for the dual variables
00186         // if you wish
00187         assert(init_x == true);
00188         assert(init_z == false);
00189         assert(init_lambda == false);
00190         int i, m1, n1;
00191 
00192 #ifdef DEBUG
00193         cout << "get initial values !!!!!!!!!!!!!!!!!!!!!!!!!! " << endl;
00194 #endif
00195         
00196         //now set initial values
00197 #ifdef DEBUG
00198         cout << "get number of initial values !!!!!!!!!!!!!!!!!!!!!!!!!! " << endl;
00199 #endif
00200         int k;
00201         if (osoption != NULL)
00202                 m1 = osoption->getNumberOfInitVarValues();
00203         else
00204                 m1 = 0;
00205 #ifdef DEBUG
00206         cout << "number of variables initialed: " << m1 << endl;
00207 #endif
00208 
00209         n1 = osinstance->getVariableNumber();
00210         bool* initialed;
00211         initialed = new bool[n1];
00212 #ifdef DEBUG
00213         cout << "number of variables in total: " << n1 << endl;
00214 #endif
00215 
00216         for(k = 0; k < n1; k++)
00217                 initialed[k] = false;
00218 
00219         if (m1 > 0)
00220         {
00221 #ifdef DEBUG
00222                         cout << "get initial values " << endl;
00223 #endif
00224 
00225                 InitVarValue**  initVarVector = osoption->getInitVarValuesSparse();
00226 #ifdef DEBUG
00227                 cout << "done " << endl;
00228 #endif
00229                 try
00230                 {
00231                         double initval;
00232                         for(k = 0; k < m1; k++)
00233                         {       i = initVarVector[k]->idx;
00234                                 if (initVarVector[k]->idx > n1)
00235                                         throw ErrorClass ("Illegal index value in variable initialization");
00236 
00237                                         initval = initVarVector[k]->value;
00238                                         if (osinstance->instanceData->variables->var[k]->ub == OSDBL_MAX)
00239                                         {       if (osinstance->instanceData->variables->var[k]->lb > initval)
00240                                                         throw ErrorClass ("Initial value outside of bounds");
00241                                         }
00242                                         else
00243                                                 if (osinstance->instanceData->variables->var[k]->lb == -OSDBL_MAX)
00244                                                 {       if (osinstance->instanceData->variables->var[k]->ub < initval)
00245                                                                 throw ErrorClass ("Initial value outside of bounds");
00246                                                 }
00247                                                 else
00248                                                 {       if ((osinstance->instanceData->variables->var[k]->lb > initval) ||
00249                                                                 (osinstance->instanceData->variables->var[k]->ub < initval))
00250                                                                 throw ErrorClass ("Initial value outside of bounds");
00251                                                 }
00252 
00253                                 x[initVarVector[k]->idx] = initval;
00254                                 initialed[initVarVector[k]->idx] = true;
00255                         }
00256                 }
00257                 catch(const ErrorClass& eclass)
00258                 {       cout << "Error in IpoptProblem::get_starting_point (OSIpoptSolver.cpp, line 247)";
00259                         cout << endl << endl << endl;
00260                 }       
00261         }  //  end if (m1 > 0)          
00262 
00263         double default_initval;
00264         default_initval = 1.7171;
00265         //default_initval = 0;
00266 
00267         for(k = 0; k < n1; k++)
00268         {       if (!initialed[k])
00269                         if (osinstance->instanceData->variables->var[k]->ub == OSDBL_MAX)
00270                                 if (osinstance->instanceData->variables->var[k]->lb <= default_initval)
00271                                         x[k] = default_initval;
00272                                 else
00273                                         x[k] = osinstance->instanceData->variables->var[k]->lb;
00274                         else
00275                                 if (osinstance->instanceData->variables->var[k]->lb == -OSDBL_MAX)
00276                                         if (osinstance->instanceData->variables->var[k]->ub >= default_initval)
00277                                                 x[k] = default_initval;
00278                                         else
00279                                                 x[k] = osinstance->instanceData->variables->var[k]->ub;
00280                                 else
00281                                         if ((osinstance->instanceData->variables->var[k]->lb <= default_initval) && 
00282                                                 (osinstance->instanceData->variables->var[k]->ub >= default_initval))
00283                                                 x[k] = default_initval;
00284                                         else
00285                                                 if (osinstance->instanceData->variables->var[k]->lb > default_initval)
00286                                                         x[k] = osinstance->instanceData->variables->var[k]->lb;
00287                                                 else
00288                                                         x[k] = osinstance->instanceData->variables->var[k]->ub;
00289         }
00290 
00291 #ifdef DEBUG
00292         for(i = 0; i < n1; i++){
00293                 std::cout << "INITIAL VALUE !!!!!!!!!!!!!!!!!!!!  " << x[ i] << std::endl;
00294         }
00295 #endif
00296         //make sure objvalue is initialized
00297         osinstance->calculateAllObjectiveFunctionValues( x, true);
00298         delete[] initialed;
00299 
00300         return true;
00301 }//get_starting_point
00302 
00303 // returns the value of the objective function
00304 bool IpoptProblem::eval_f(Index n, const Number* x, bool new_x, Number& obj_value){
00305         try{
00306         
00307                 //the following is a kludge for ipopt, new_x does not seem to get initilized if there are no constraints.
00308                 //if(osinstance->getConstraintNumber() <= 0) new_x = true;
00309                 if(new_x == false) obj_value  = osinstance->calculateAllObjectiveFunctionValues( const_cast<double*>(x), false)[ 0];
00310                         else obj_value = osinstance->calculateAllObjectiveFunctionValues( const_cast<double*>(x), NULL, NULL, true, 0 )[ 0];
00311                 //if( CoinIsnan( (double)obj_value) ) return false;
00312                 if( CoinIsnan( obj_value ) )return false;
00313 
00314         }
00315         catch(const ErrorClass& eclass){
00316 
00317                 *ipoptErrorMsg = eclass.errormsg;
00318                 throw;  
00319         }
00320         
00321         return true;
00322 }
00323 
00324 bool IpoptProblem::eval_grad_f(Index n, const Number* x, bool new_x, Number* grad_f){
00325         int i;
00326         double *objGrad;
00327         try{
00328                 //objGrad = osinstance->calculateAllObjectiveFunctionGradients( const_cast<double*>(x), NULL, NULL,  new_x, 1)[ 0];
00329                 objGrad = osinstance->calculateObjectiveFunctionGradient( const_cast<double*>(x), NULL, NULL, -1,  new_x, 1);
00330         }
00331         catch(const ErrorClass& eclass){
00332 #ifdef DEBUG
00333                 cout << "error in OSIpoptSolver, line 314:\n" << eclass.errormsg << endl;
00334 #endif
00335                 *ipoptErrorMsg = eclass.errormsg;
00336                 throw;  
00337         }
00338         for(i = 0; i < n; i++){
00339                 grad_f[ i]  = objGrad[ i];
00340         }
00341         return true;
00342 }//eval_grad_f
00343 
00344 // return the value of the constraints: g(x)
00345 bool IpoptProblem::eval_g(Index n, const Number* x, bool new_x, Index m, Number* g) {
00346         try{
00347                 double *conVals = osinstance->calculateAllConstraintFunctionValues( const_cast<double*>(x), NULL, NULL, new_x, 0 );
00348                 int i;
00349                 for(i = 0; i < m; i++){
00350                         if( CoinIsnan( (double)conVals[ i] ) ) return false;
00351                         g[i] = conVals[ i]  ;   
00352                 } 
00353                 return true;
00354         }
00355         catch(const ErrorClass& eclass){
00356 #ifdef DEBUG
00357                 cout << "error in OSIpoptSolver, line 338:\n" << eclass.errormsg << endl;
00358 #endif
00359                 *ipoptErrorMsg = eclass.errormsg;
00360                 throw;  
00361         }
00362 }//eval_g
00363 
00364 
00365 // return the structure or values of the jacobian
00366 bool IpoptProblem::eval_jac_g(Index n, const Number* x, bool new_x,
00367                            Index m, Index nele_jac, Index* iRow, Index *jCol,
00368                            Number* values){
00369         SparseJacobianMatrix *sparseJacobian;
00370         if (values == NULL) {
00371                 // return the values of the jacobian of the constraints
00372                 //cout << "n: " << n << endl;
00373                 //cout << "m: " << m << endl;
00374                 //cout << "nele_jac: " <<  nele_jac << endl;
00375                 // return the structure of the jacobian
00376                 try{
00377                         sparseJacobian = osinstance->getJacobianSparsityPattern();
00378                 }
00379                 catch(const ErrorClass& eclass){
00380 #ifdef DEBUG
00381                 cout << "error in OSIpoptSolver, line 362:\n" << eclass.errormsg << endl;
00382 #endif
00383                         *ipoptErrorMsg =  eclass.errormsg; 
00384                         throw; 
00385                 }
00386                 int i = 0;
00387                 int k, idx;
00388                 for(idx = 0; idx < m; idx++){
00389                         for(k = *(sparseJacobian->starts + idx); k < *(sparseJacobian->starts + idx + 1); k++){
00390                                 iRow[i] = idx;
00391                                 jCol[i] = *(sparseJacobian->indexes + k);
00392                                 //cout << "ROW IDX  !!!!!!!!!!!!!!!!!!!!!!!!!!!"  << iRow[i] << endl;
00393                                 //cout << "COL IDX  !!!!!!!!!!!!!!!!!!!!!!!!!!!"  << jCol[i] << endl;
00394                                 i++;
00395                         }
00396                 }
00397         }
00398         else {
00399                 //std::cout << "EVALUATING JACOBIAN" << std::endl; 
00400                 try{
00401                         sparseJacobian = osinstance->calculateAllConstraintFunctionGradients( const_cast<double*>(x), NULL, NULL,  new_x, 1);
00402                 }
00403                 catch(const ErrorClass& eclass){
00404 #ifdef DEBUG
00405                 cout << "error in OSIpoptSolver, line 386:\n" << eclass.errormsg << endl;
00406 #endif
00407                         *ipoptErrorMsg = eclass.errormsg;
00408                         throw;  
00409                 }
00410                 //osinstance->getIterateResults( (double*)x, 0.0, NULL, -1, new_x,  1);
00411                 for(int i = 0; i < nele_jac; i++){
00412                         values[ i] = sparseJacobian->values[i];
00413                         //values[ i] = osinstance->m_mdJacValue[ i];
00414                         //cout << "values[i]:!!!!!!!!!!!!  " <<  values[ i] << endl;
00415                         //cout << "m_mdJacValue[ i]:!!!!!!!!!!!!  " <<  osinstance->m_mdJacValue[ i] << endl;
00416                 }
00417         }
00418   return true;
00419 }//eval_jac_g
00420 
00421 //return the structure or values of the hessian
00422 bool IpoptProblem::eval_h(Index n, const Number* x, bool new_x,
00423                        Number obj_factor, Index m, const Number* lambda,
00424                        bool new_lambda, Index nele_hess, Index* iRow,
00425                        Index* jCol, Number* values){
00426 
00428         SparseHessianMatrix *sparseHessian;
00429         
00430         int i;
00431         if (values == NULL) {
00432                 // return the structure. This is a symmetric matrix, fill the lower left triangle only.
00433                 //cout << "get structure of HESSIAN !!!!!!!!!!!!!!!!!!!!!!!!!! "  << endl;
00434                 try{
00435                         sparseHessian = osinstance->getLagrangianHessianSparsityPattern( );
00436                 }
00437                 catch(const ErrorClass& eclass){
00438 
00439                         *ipoptErrorMsg = eclass.errormsg;
00440                         throw;  
00441                 }
00442 //              cout << "got structure of HESSIAN !!!!!!!!!!!!!!!!!!!!!!!!!! "  << endl;
00443                 for(i = 0; i < nele_hess; i++){
00444                         iRow[i] = *(sparseHessian->hessColIdx + i);
00445                         jCol[i] = *(sparseHessian->hessRowIdx + i);
00446                                 //cout << "ROW HESS IDX  !!!!!!!!!!!!!!!!!!!!!!!!!!!"  << iRow[i] << endl;
00447                                 //cout << "COL HESS IDX  !!!!!!!!!!!!!!!!!!!!!!!!!!!"  << jCol[i] << endl;
00448                 }
00449         }
00450         else {
00451                 //std::cout << "EVALUATING HESSIAN" << std::endl; 
00452                 // return the values. This is a symmetric matrix, fill the lower left triangle only
00453                 double* objMultipliers = new double[1];
00454                 objMultipliers[0] = obj_factor;
00455                 try{
00456                         sparseHessian = osinstance->calculateLagrangianHessian( const_cast<double*>(x), objMultipliers, const_cast<double*>(lambda) ,  new_x, 2);
00457                 delete[]  objMultipliers;
00458                 }
00459                 catch(const ErrorClass& eclass){
00460 #ifdef DEBUG
00461                 cout << "error in OSIpoptSolver, line 444:\n" << eclass.errormsg << endl;
00462 #endif
00463                         *ipoptErrorMsg = eclass.errormsg;
00464                         delete[]  objMultipliers;
00465                         throw;  
00466                 }
00467                 for(i = 0; i < nele_hess; i++){
00468                         values[ i]  = *(sparseHessian->hessValues + i);
00469                 }
00470         }
00472         return true;
00473 }//eval_h
00474 
00475 bool IpoptProblem::get_scaling_parameters(Number& obj_scaling,
00476                         bool& use_x_scaling, Index n,
00477                         Number* x_scaling,
00478                     bool& use_g_scaling, Index m,
00479                     Number* g_scaling){
00480         if( osinstance->instanceData->objectives->obj[ 0]->maxOrMin.compare("min") != 0){
00481                 obj_scaling = -1;
00482         }
00483     else obj_scaling = 1;
00484     use_x_scaling = false;
00485     use_g_scaling = false;
00486         return true;
00487 }//get_scaling_parameters
00488 
00489 void IpoptProblem::finalize_solution(SolverReturn status,
00490                                Index n, const Number* x, const Number* z_L, const Number* z_U,
00491                                   Index m, const Number* g, const Number* lambda,
00492                                   Number obj_value,
00493                                    const IpoptData* ip_data,
00494                                    IpoptCalculatedQuantities* ip_cq)
00495 {
00496           // here is where we would store the solution to variables, or write to a file, etc
00497           // so we could use the solution.
00498           // For this example, we write the solution to the console
00499         OSrLWriter *osrlwriter ;
00500         osrlwriter = new OSrLWriter();
00501 #ifdef DEBUG
00502           printf("\n\nSolution of the primal variables, x\n");
00503           for (Index i=0; i<n; i++) {
00504             printf("x[%d] = %e\n", i, x[i]);
00505           }
00506         
00507           printf("\n\nSolution of the bound multipliers, z_L and z_U\n");
00508           for (Index i=0; i<n; i++) {
00509             printf("z_L[%d] = %e\n", i, z_L[i]);
00510           }
00511           for (Index i=0; i<n; i++) {
00512             printf("z_U[%d] = %e\n", i, z_U[i]);
00513           }
00514 #endif
00515           printf("\nObjective value f(x*) = %e\n", obj_value);
00516           
00517         int solIdx = 0;
00518         int numberOfOtherVariableResults;
00519         int otherIdx;
00520         ostringstream outStr;
00521         
00522         std::string *rcost = NULL;
00523         int* idx = NULL;
00524         double* mdObjValues = new double[1];
00525         std::string message = "Ipopt solver finishes to the end.";
00526         std::string solutionDescription = "";   
00527 
00528         try{
00529                 // resultHeader infomration
00530                 if(osresult->setServiceName( "Ipopt solver service") != true)
00531                         throw ErrorClass("OSResult error: setServiceName");
00532                 if(osresult->setInstanceName(  osinstance->getInstanceName()) != true)
00533                         throw ErrorClass("OSResult error: setInstanceName");
00534 
00535                 //if(osresult->setJobID( osoption->jobID) != true)
00536                 //      throw ErrorClass("OSResult error: setJobID");
00537 
00538                 // set basic problem parameters
00539                 if(osresult->setVariableNumber( osinstance->getVariableNumber()) != true)
00540                         throw ErrorClass("OSResult error: setVariableNumer");
00541                 if(osresult->setObjectiveNumber( 1) != true)
00542                         throw ErrorClass("OSResult error: setObjectiveNumber");
00543                 if(osresult->setConstraintNumber( osinstance->getConstraintNumber()) != true)
00544                         throw ErrorClass("OSResult error: setConstraintNumber");
00545                 if(osresult->setSolutionNumber(  1) != true)
00546                         throw ErrorClass("OSResult error: setSolutionNumer");   
00547 
00548 
00549                 if(osresult->setGeneralMessage( message) != true)
00550                         throw ErrorClass("OSResult error: setGeneralMessage");
00551 
00552                 switch( status){
00553                         case SUCCESS:
00554                                 solutionDescription = "SUCCESS[IPOPT]: Algorithm terminated normally at a locally optimal point, satisfying the convergence tolerances.";
00555                                 osresult->setSolutionStatus(solIdx,  "locallyOptimal", solutionDescription);
00556                                 osresult->setPrimalVariableValuesDense(solIdx, const_cast<double*>(x)); 
00557                                 osresult->setDualVariableValuesDense(solIdx, const_cast<double*>( lambda)); 
00558                                 mdObjValues[0] = obj_value;
00559                                 osresult->setObjectiveValuesDense(solIdx, mdObjValues); 
00560                                 
00561                                 
00562                                 // set other
00563                                 
00564                                 numberOfOtherVariableResults = 2;
00565                                 osresult->setNumberOfOtherVariableResults(solIdx, numberOfOtherVariableResults);
00566                                 
00567                         
00568                                 rcost = new std::string[ osinstance->getVariableNumber()];
00569                                 idx = new int[ osinstance->getVariableNumber()];
00570                                 
00571 
00572                                 for (Index i = 0; i < n; i++) {
00573                                     rcost[ i] =  os_dtoa_format( z_L[i]); 
00574                                         idx[ i] = i;
00575                                  }
00576                                 otherIdx = 0;
00577                                 osresult->setAnOtherVariableResultSparse(solIdx, otherIdx, "varL", "", "Lagrange Multiplier on the Variable Lower Bound",  idx,  rcost,  osinstance->getVariableNumber() );
00578                                 
00579                                 
00580                                 for (Index i = 0; i < n; i++) {
00581                                     rcost[ i] =  os_dtoa_format( z_U[i]); 
00582                                  }
00583                                 otherIdx = 1;
00584                                 osresult->setAnOtherVariableResultSparse(solIdx, otherIdx, "varU", "", "Lagrange Multiplier on the Variable Upper Bound", idx,  rcost, osinstance->getVariableNumber() );
00585                                 
00586                                 delete[] rcost;
00587                                 delete[] idx;
00588                                 
00589                                 // done with other
00590                                 
00591                                 
00592                         break;
00593                         case MAXITER_EXCEEDED:
00594                                 solutionDescription = "MAXITER_EXCEEDED[IPOPT]: Maximum number of iterations exceeded.";
00595                                 osresult->setSolutionStatus(solIdx,  "stoppedByLimit", solutionDescription);
00596                                 osresult->setPrimalVariableValuesDense(solIdx, const_cast<double*>(x)); 
00597                                 osresult->setDualVariableValuesDense(solIdx, const_cast<double*>( lambda)); 
00598                                 mdObjValues[0] = obj_value;
00599                                 osresult->setObjectiveValuesDense(solIdx, mdObjValues); 
00600                         break;
00601                         case STOP_AT_TINY_STEP:
00602                                 solutionDescription = "STOP_AT_TINY_STEP[IPOPT]: Algorithm proceeds with very little progress.";
00603                                 osresult->setSolutionStatus(solIdx,  "stoppedByLimit", solutionDescription);
00604                                 osresult->setPrimalVariableValuesDense(solIdx, const_cast<double*>( x)); 
00605                                 osresult->setDualVariableValuesDense(solIdx, const_cast<double*>( lambda)); 
00606                                 mdObjValues[0] = obj_value;
00607                                 osresult->setObjectiveValuesDense(solIdx, mdObjValues); 
00608                         break;
00609                         case STOP_AT_ACCEPTABLE_POINT:
00610                                 solutionDescription = "STOP_AT_ACCEPTABLE_POINT[IPOPT]: Algorithm stopped at a point that was converged, not to _desired_ tolerances, but to _acceptable_ tolerances";
00611                                 osresult->setSolutionStatus(solIdx,  "IpoptAccetable", solutionDescription);
00612                                 osresult->setPrimalVariableValuesDense(solIdx, const_cast<double*>(x)); 
00613                                 osresult->setDualVariableValuesDense(solIdx, const_cast<double*>( lambda)); 
00614                                 mdObjValues[0] = obj_value;
00615                                 osresult->setObjectiveValuesDense(solIdx, mdObjValues); 
00616                         break;
00617                         case LOCAL_INFEASIBILITY:
00618                                 solutionDescription = "LOCAL_INFEASIBILITY[IPOPT]: Algorithm converged to a point of local infeasibility. Problem may be infeasible.";
00619                                 osresult->setSolutionStatus(solIdx,  "infeasible", solutionDescription);
00620                         break;
00621                         case USER_REQUESTED_STOP:
00622                                 solutionDescription = "USER_REQUESTED_STOP[IPOPT]: The user call-back function  intermediate_callback returned false, i.e., the user code requested a premature termination of the optimization.";
00623                                 osresult->setSolutionStatus(solIdx,  "error", solutionDescription);
00624                         break;
00625                         case DIVERGING_ITERATES:
00626                                 solutionDescription = "DIVERGING_ITERATES[IPOPT]: It seems that the iterates diverge.";
00627                                 osresult->setSolutionStatus(solIdx,  "unbounded", solutionDescription);
00628                         break;
00629                         case RESTORATION_FAILURE:
00630                                 solutionDescription = "RESTORATION_FAILURE[IPOPT]: Restoration phase failed, algorithm doesn't know how to proceed.";
00631                                 osresult->setSolutionStatus(solIdx,  "error", solutionDescription);
00632                         break;
00633                         case ERROR_IN_STEP_COMPUTATION:
00634                                 solutionDescription = "ERROR_IN_STEP_COMPUTATION[IPOPT]: An unrecoverable error occurred while IPOPT tried to compute the search direction.";
00635                                 osresult->setSolutionStatus(solIdx,  "error", solutionDescription);
00636                         break;
00637                         case INVALID_NUMBER_DETECTED:
00638                                 solutionDescription = "INVALID_NUMBER_DETECTED[IPOPT]: Algorithm received an invalid number (such as NaN or Inf) from the NLP; see also option check_derivatives_for_naninf.";
00639                                 osresult->setSolutionStatus(solIdx,  "error", solutionDescription);
00640                         break;
00641                         case INTERNAL_ERROR:
00642                                 solutionDescription = "INTERNAL_ERROR[IPOPT]: An unknown internal error occurred. Please contact the IPOPT authors through the mailing list.";
00643                                 osresult->setSolutionStatus(solIdx,  "error", solutionDescription);
00644                         break;
00645                         default:
00646                                 solutionDescription = "OTHER[IPOPT]: other unknown solution status from Ipopt solver";
00647                                 osresult->setSolutionStatus(solIdx,  "other", solutionDescription);
00648                 }
00649                 osresult->setGeneralStatusType("normal");
00650                 delete osrlwriter;
00651                 delete[] mdObjValues;
00652                 osrlwriter = NULL;
00653 
00654         }
00655         catch(const ErrorClass& eclass){
00656 #ifdef DEBUG
00657                 cout << "error in OSIpoptSolver, line 636:\n" << eclass.errormsg << endl;
00658 #endif
00659                 osresult->setGeneralMessage( eclass.errormsg);
00660                 osresult->setGeneralStatusType( "error");
00661                 std::string osrl = osrlwriter->writeOSrL( osresult);
00662                 delete osrlwriter;
00663                 osrlwriter = NULL;
00664                 throw ErrorClass(  osrl) ;
00665                 delete[] mdObjValues;
00666                 mdObjValues = NULL;
00667         }
00669 }
00670 
00671 
00672 void IpoptSolver::setSolverOptions() throw (ErrorClass) {
00673         try{
00674                 if(osinstance->getObjectiveNumber() <= 0) throw ErrorClass("Ipopt NEEDS AN OBJECTIVE FUNCTION");
00675                 this->bSetSolverOptions = true;
00676                 /* set the default options */   
00677                 //app->Options()->SetNumericValue("tol", 1e-9);
00678                 app->Options()->SetIntegerValue("print_level", 0);
00679                 app->Options()->SetIntegerValue("max_iter", 20000);
00680                 app->Options()->SetNumericValue("bound_relax_factor", 0, true, true);
00681                 app->Options()->SetStringValue("mu_strategy", "adaptive", true, true);
00682                 //app->Options()->SetStringValue("output_file", "ipopt.out");
00683                 app->Options()->SetStringValue("check_derivatives_for_naninf", "yes");
00684                 // hessian constant for an LP
00685                 if( (osinstance-> getNumberOfNonlinearExpressions() <= 0) && (osinstance->getNumberOfQuadraticTerms() <= 0) ){
00686                         app->Options()->SetStringValue("hessian_constant", "yes", true, true);
00687                 }
00688                 if( osinstance->instanceData->objectives->obj[ 0]->maxOrMin.compare("min") != 0){
00689                         app->Options()->SetStringValue("nlp_scaling_method", "user-scaling");
00690                 }
00691                 /* end of the default options, now get options from OSoL */
00692         
00693                 
00694                 if(osoption == NULL && osol.length() > 0)
00695                 {
00696                         m_osolreader = new OSoLReader();
00697                         osoption = m_osolreader->readOSoL( osol);
00698                 }
00699 
00700                 if( osoption != NULL  &&  osoption->getNumberOfSolverOptions() > 0 ){
00701 //                      std::cout << "number of solver options "  <<  osoption->getNumberOfSolverOptions() << std::endl;
00702                         std::vector<SolverOption*> optionsVector;
00703                         optionsVector = osoption->getSolverOptions( "ipopt");
00704                         char *pEnd;
00705                         int i;
00706                         int num_ipopt_options = optionsVector.size();
00707                         for(i = 0; i < num_ipopt_options; i++){
00708 //                              std::cout << "ipopt solver option  "  << optionsVector[ i]->name << std::endl;
00709                                 if(optionsVector[ i]->type == "numeric" ){
00710 //                                      std::cout << "FOUND A NUMERIC OPTION  "  <<  os_strtod( optionsVector[ i]->value.c_str(), &pEnd ) << std::endl;
00711                                         app->Options()->SetNumericValue(optionsVector[ i]->name, os_strtod( optionsVector[ i]->value.c_str(), &pEnd ) );        
00712                                 }
00713                                 else if(optionsVector[ i]->type == "integer" ){
00714 //                                      std::cout << "FOUND AN INTEGER OPTION  "  << atoi( optionsVector[ i]->value.c_str() ) << std::endl;
00715                                         app->Options()->SetIntegerValue(optionsVector[ i]->name, atoi( optionsVector[ i]->value.c_str() ) );
00716                                 }
00717                                 else if(optionsVector[ i]->type == "string" ){
00718 //                                      std::cout << "FOUND A STRING OPTION  "  <<  optionsVector[ i]->value.c_str() << std::endl;
00719                                         app->Options()->SetStringValue(optionsVector[ i]->name, optionsVector[ i]->value);
00720                                 }
00721                         }       
00722                 }
00723         }
00724         catch(const ErrorClass& eclass){
00725 #ifdef DEBUG
00726                 cout << "error in OSIpoptSolver, line 695:\n" << eclass.errormsg << endl;
00727 #endif
00728                 std::cout << "THERE IS AN ERROR" << std::endl;
00729                 osresult->setGeneralMessage( eclass.errormsg);
00730                 osresult->setGeneralStatusType( "error");
00731                 osrl = osrlwriter->writeOSrL( osresult);
00732                 throw ErrorClass( osrl) ;
00733         }                               
00734 }//end setSolverOptions() 
00735 
00736 
00737 
00738 void IpoptSolver::buildSolverInstance() throw (ErrorClass) {
00739         try{
00740                 
00741                 if(osil.length() == 0 && osinstance == NULL) throw ErrorClass("there is no instance");
00742                 if(osinstance == NULL){
00743                         m_osilreader = new OSiLReader();
00744                         osinstance = m_osilreader->readOSiL( osil);
00745                 }
00746                 // Create a new instance of your nlp 
00747                 nlp = new IpoptProblem( osinstance, osoption, osresult, ipoptErrorMsg);
00748                 app = new IpoptApplication();
00749                 this->bCallbuildSolverInstance = true;
00750         }
00751         catch(const ErrorClass& eclass){
00752 #ifdef DEBUG
00753                 cout << "error in OSIpoptSolver, line 722:\n" << eclass.errormsg << endl;
00754 #endif
00755                 std::cout << "THERE IS AN ERROR" << std::endl;
00756                 osresult->setGeneralMessage( eclass.errormsg);
00757                 osresult->setGeneralStatusType( "error");
00758                 osrl = osrlwriter->writeOSrL( osresult);
00759                 throw ErrorClass( osrl) ;
00760         }                               
00761 }//end buildSolverInstance() 
00762 
00763 
00764 
00765 void IpoptSolver::solve() throw (ErrorClass) {
00766         if( this->bCallbuildSolverInstance == false) buildSolverInstance();
00767         if( this->bSetSolverOptions == false) setSolverOptions();
00768         try{
00769                 clock_t start, finish;
00770                 double duration;
00771                 start = clock();
00772                 //OSiLWriter osilwriter;
00773                 //cout << osilwriter.writeOSiL( osinstance) << endl;
00774                 if(osinstance->getVariableNumber() <= 0)throw ErrorClass("Ipopt requires decision variables");
00775                 finish = clock();
00776                 duration = (double) (finish - start) / CLOCKS_PER_SEC;
00777                 //cout << "Parsing took (seconds): " << duration << endl; 
00778                 //dataEchoCheck();
00779                 /***************now the ipopt invokation*********************/
00780                 // see if we have a linear program
00781                 if(osinstance->getObjectiveNumber() <= 0) throw ErrorClass("Ipopt NEEDS AN OBJECTIVE FUNCTION");
00782                 // Intialize the IpoptApplication and process the options
00783 //              std::cout << "Call Ipopt Initialize" << std::endl;
00784                 app->Initialize();
00785 //              std::cout << "Finished Ipopt Initialize" << std::endl;
00786                 //nlp->osinstance = this->osinstance;
00787                 // Ask Ipopt to solve the problem
00788 //              std::cout << "Call Ipopt Optimize" << std::endl;
00789                 ApplicationReturnStatus status = app->OptimizeTNLP( nlp);
00790 //              std::cout << "Finish Ipopt Optimize" << std::endl;
00791                 osrl = osrlwriter->writeOSrL( osresult);
00792 //              std::cout << "Finish writing the osrl" << std::endl;
00793                 //if (status != Solve_Succeeded) {
00794                 if (status < -2) {
00795                         throw ErrorClass("Ipopt FAILED TO SOLVE THE PROBLEM: " + *ipoptErrorMsg);
00796                 }       
00797         }
00798         catch(const ErrorClass& eclass){
00799 #ifdef DEBUG
00800                 cout << "error in OSIpoptSolver, line 775:\n" << eclass.errormsg << endl;
00801 #endif
00802                 osresult->setGeneralMessage( eclass.errormsg);
00803                 osresult->setGeneralStatusType( "error");
00804                 osrl = osrlwriter->writeOSrL( osresult);
00805                 throw ErrorClass( osrl) ;
00806         }
00807 }//solve
00808 
00809 
00810 void IpoptSolver::dataEchoCheck(){
00811         
00812         int i;
00813         
00814         // print out problem parameters
00815         cout << "This is problem:  " << osinstance->getInstanceName() << endl;
00816         cout << "The problem source is:  " << osinstance->getInstanceSource() << endl;
00817         cout << "The problem description is:  " << osinstance->getInstanceDescription() << endl;
00818         cout << "number of variables = " << osinstance->getVariableNumber() << endl;
00819         cout << "number of Rows = " << osinstance->getConstraintNumber() << endl;
00820 
00821         // print out the variable information
00822         if(osinstance->getVariableNumber() > 0){
00823                 for(i = 0; i < osinstance->getVariableNumber(); i++){
00824                         if(osinstance->getVariableNames() != NULL) cout << "variable Names  " << osinstance->getVariableNames()[ i]  << endl;
00825                         if(osinstance->getVariableTypes() != NULL) cout << "variable Types  " << osinstance->getVariableTypes()[ i]  << endl;
00826                         if(osinstance->getVariableLowerBounds() != NULL) cout << "variable Lower Bounds  " << osinstance->getVariableLowerBounds()[ i]  << endl;
00827                         if(osinstance->getVariableUpperBounds() != NULL) cout << "variable Upper Bounds  " <<  osinstance->getVariableUpperBounds()[i] << endl;
00828                 }
00829         }
00830         
00831         // print out objective function information
00832         if(osinstance->getVariableNumber() > 0 || osinstance->instanceData->objectives->obj != NULL || osinstance->instanceData->objectives->numberOfObjectives > 0){
00833                 if( osinstance->getObjectiveMaxOrMins()[0] == "min")  cout <<  "problem is a minimization" << endl;
00834                 else cout <<  "problem is a maximization" << endl;
00835                 for(i = 0; i < osinstance->getVariableNumber(); i++){
00836                         cout << "OBJ COEFFICIENT =  " <<  osinstance->getDenseObjectiveCoefficients()[0][i] << endl;
00837                 }
00838         }
00839         // print out constraint information
00840         if(osinstance->getConstraintNumber() > 0){
00841                 for(i = 0; i < osinstance->getConstraintNumber(); i++){
00842                         if(osinstance->getConstraintNames() != NULL) cout << "row name = " << osinstance->getConstraintNames()[i] <<  endl;
00843                         if(osinstance->getConstraintLowerBounds() != NULL) cout << "row lower bound = " << osinstance->getConstraintLowerBounds()[i] <<  endl;
00844                         if(osinstance->getConstraintUpperBounds() != NULL) cout << "row upper bound = " << osinstance->getConstraintUpperBounds()[i] <<  endl; 
00845                 }
00846         }
00847         
00848         // print out linear constraint data
00849         cout << endl;
00850         cout << "number of nonzeros =  " << osinstance->getLinearConstraintCoefficientNumber() << endl;
00851         for(i = 0; i <= osinstance->getVariableNumber(); i++){
00852                 cout << "Start Value =  " << osinstance->getLinearConstraintCoefficientsInColumnMajor()->starts[ i] << endl;
00853         }
00854         cout << endl;
00855         for(i = 0; i < osinstance->getLinearConstraintCoefficientNumber(); i++){
00856                 cout << "Index Value =  " << osinstance->getLinearConstraintCoefficientsInColumnMajor()->indexes[i] << endl;
00857                 cout << "Nonzero Value =  " << osinstance->getLinearConstraintCoefficientsInColumnMajor()->values[i] << endl;
00858         }
00859 
00860         // print out quadratic data
00861         cout << "number of qterms =  " <<  osinstance->getNumberOfQuadraticTerms() << endl;
00862         for(int i = 0; i <  osinstance->getNumberOfQuadraticTerms(); i++){
00863                 cout << "Row Index = " <<  osinstance->getQuadraticTerms()->rowIndexes[i] << endl;
00864                 cout << "Var Index 1 = " << osinstance->getQuadraticTerms()->varOneIndexes[ i] << endl;
00865                 cout << "Var Index 2 = " << osinstance->getQuadraticTerms()->varTwoIndexes[ i] << endl;
00866                 cout << "Coefficient = " << osinstance->getQuadraticTerms()->coefficients[ i] << endl;
00867         }
00868 } // end dataEchoCheck
00869 
00870 
00871 IpoptProblem::IpoptProblem(OSInstance *osinstance_,  OSOption *osoption_, OSResult *osresult_, std::string* ipoptErrorMsg_) {
00872         osinstance = osinstance_;
00873         osoption = osoption_;
00874         osresult = osresult_;
00875         ipoptErrorMsg = ipoptErrorMsg_;
00876 }
00877 
00878 IpoptProblem::~IpoptProblem() {
00879 
00880 }
00881 
00882 
00883 

Generated on Thu Oct 8 03:03:01 2009 by  doxygen 1.4.7