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

Go to the documentation of this file.
00001 /* $Id: OSIpoptSolver.cpp 3479 2010-04-26 07:14:52Z 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                 if(osinstance->getObjectiveNumber() > 0){
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         }
00316         catch(const ErrorClass& eclass){
00317 
00318                 *ipoptErrorMsg = eclass.errormsg;
00319                 throw;  
00320         }
00321         
00322         return true;
00323 }
00324 
00325 bool IpoptProblem::eval_grad_f(Index n, const Number* x, bool new_x, Number* grad_f){
00326         int i;
00327         double *objGrad = NULL;
00328         if(osinstance->getObjectiveNumber() > 0){
00329                 try{
00330                         //objGrad = osinstance->calculateAllObjectiveFunctionGradients( const_cast<double*>(x), NULL, NULL,  new_x, 1)[ 0];
00331                         objGrad = osinstance->calculateObjectiveFunctionGradient( const_cast<double*>(x), NULL, NULL, -1,  new_x, 1);
00332                 }
00333                 catch(const ErrorClass& eclass){
00334         #ifdef DEBUG
00335                         cout << "error in OSIpoptSolver, line 314:\n" << eclass.errormsg << endl;
00336         #endif
00337                         *ipoptErrorMsg = eclass.errormsg;
00338                         throw;  
00339                 }
00340                 for(i = 0; i < n; i++){
00341                         grad_f[ i]  = objGrad[ i];
00342                 }
00343         }
00344         return true;
00345 }//eval_grad_f
00346 
00347 // return the value of the constraints: g(x)
00348 bool IpoptProblem::eval_g(Index n, const Number* x, bool new_x, Index m, Number* g) {
00349         try{
00350                 double *conVals = osinstance->calculateAllConstraintFunctionValues( const_cast<double*>(x), NULL, NULL, new_x, 0 );
00351                 int i;
00352                 for(i = 0; i < m; i++){
00353                         if( CoinIsnan( (double)conVals[ i] ) ) return false;
00354                         g[i] = conVals[ i]  ;   
00355                 } 
00356                 return true;
00357         }
00358         catch(const ErrorClass& eclass){
00359 #ifdef DEBUG
00360                 cout << "error in OSIpoptSolver, line 338:\n" << eclass.errormsg << endl;
00361 #endif
00362                 *ipoptErrorMsg = eclass.errormsg;
00363                 throw;  
00364         }
00365 }//eval_g
00366 
00367 
00368 // return the structure or values of the jacobian
00369 bool IpoptProblem::eval_jac_g(Index n, const Number* x, bool new_x,
00370                            Index m, Index nele_jac, Index* iRow, Index *jCol,
00371                            Number* values){
00372         SparseJacobianMatrix *sparseJacobian;
00373         if (values == NULL) {
00374                 // return the values of the jacobian of the constraints
00375                 //cout << "n: " << n << endl;
00376                 //cout << "m: " << m << endl;
00377                 //cout << "nele_jac: " <<  nele_jac << endl;
00378                 // return the structure of the jacobian
00379                 try{
00380                         sparseJacobian = osinstance->getJacobianSparsityPattern();
00381                 }
00382                 catch(const ErrorClass& eclass){
00383 #ifdef DEBUG
00384                 cout << "error in OSIpoptSolver, line 362:\n" << eclass.errormsg << endl;
00385 #endif
00386                         *ipoptErrorMsg =  eclass.errormsg; 
00387                         throw; 
00388                 }
00389                 int i = 0;
00390                 int k, idx;
00391                 for(idx = 0; idx < m; idx++){
00392                         for(k = *(sparseJacobian->starts + idx); k < *(sparseJacobian->starts + idx + 1); k++){
00393                                 iRow[i] = idx;
00394                                 jCol[i] = *(sparseJacobian->indexes + k);
00395                                 //cout << "ROW IDX  !!!!!!!!!!!!!!!!!!!!!!!!!!!"  << iRow[i] << endl;
00396                                 //cout << "COL IDX  !!!!!!!!!!!!!!!!!!!!!!!!!!!"  << jCol[i] << endl;
00397                                 i++;
00398                         }
00399                 }
00400         }
00401         else {
00402                 //std::cout << "EVALUATING JACOBIAN" << std::endl; 
00403                 try{
00404                         sparseJacobian = osinstance->calculateAllConstraintFunctionGradients( const_cast<double*>(x), NULL, NULL,  new_x, 1);
00405                 }
00406                 catch(const ErrorClass& eclass){
00407 #ifdef DEBUG
00408                 cout << "error in OSIpoptSolver, line 386:\n" << eclass.errormsg << endl;
00409 #endif
00410                         *ipoptErrorMsg = eclass.errormsg;
00411                         throw;  
00412                 }
00413                 //osinstance->getIterateResults( (double*)x, 0.0, NULL, -1, new_x,  1);
00414                 for(int i = 0; i < nele_jac; i++){
00415                         values[ i] = sparseJacobian->values[i];
00416                         //values[ i] = osinstance->m_mdJacValue[ i];
00417                         //cout << "values[i]:!!!!!!!!!!!!  " <<  values[ i] << endl;
00418                         //cout << "m_mdJacValue[ i]:!!!!!!!!!!!!  " <<  osinstance->m_mdJacValue[ i] << endl;
00419                 }
00420         }
00421   return true;
00422 }//eval_jac_g
00423 
00424 //return the structure or values of the hessian
00425 bool IpoptProblem::eval_h(Index n, const Number* x, bool new_x,
00426                        Number obj_factor, Index m, const Number* lambda,
00427                        bool new_lambda, Index nele_hess, Index* iRow,
00428                        Index* jCol, Number* values){
00429 
00431         SparseHessianMatrix *sparseHessian;
00432         
00433         int i;
00434         if (values == NULL) {
00435                 // return the structure. This is a symmetric matrix, fill the lower left triangle only.
00436                 //cout << "get structure of HESSIAN !!!!!!!!!!!!!!!!!!!!!!!!!! "  << endl;
00437                 try{
00438                         sparseHessian = osinstance->getLagrangianHessianSparsityPattern( );
00439                 }
00440                 catch(const ErrorClass& eclass){
00441 
00442                         *ipoptErrorMsg = eclass.errormsg;
00443                         throw;  
00444                 }
00445 //              cout << "got structure of HESSIAN !!!!!!!!!!!!!!!!!!!!!!!!!! "  << endl;
00446                 for(i = 0; i < nele_hess; i++){
00447                         iRow[i] = *(sparseHessian->hessColIdx + i);
00448                         jCol[i] = *(sparseHessian->hessRowIdx + i);
00449                                 //cout << "ROW HESS IDX  !!!!!!!!!!!!!!!!!!!!!!!!!!!"  << iRow[i] << endl;
00450                                 //cout << "COL HESS IDX  !!!!!!!!!!!!!!!!!!!!!!!!!!!"  << jCol[i] << endl;
00451                 }
00452         }
00453         else {
00454                 //std::cout << "EVALUATING HESSIAN" << std::endl; 
00455                 // return the values. This is a symmetric matrix, fill the lower left triangle only
00456                 double* objMultipliers = new double[1];
00457                 objMultipliers[0] = obj_factor;
00458                 try{
00459                         sparseHessian = osinstance->calculateLagrangianHessian( const_cast<double*>(x), objMultipliers, const_cast<double*>(lambda) ,  new_x, 2);
00460                 delete[]  objMultipliers;
00461                 }
00462                 catch(const ErrorClass& eclass){
00463 #ifdef DEBUG
00464                 cout << "error in OSIpoptSolver, line 444:\n" << eclass.errormsg << endl;
00465 #endif
00466                         *ipoptErrorMsg = eclass.errormsg;
00467                         delete[]  objMultipliers;
00468                         throw;  
00469                 }
00470                 for(i = 0; i < nele_hess; i++){
00471                         values[ i]  = *(sparseHessian->hessValues + i);
00472                 }
00473         }
00475         return true;
00476 }//eval_h
00477 
00478 bool IpoptProblem::get_scaling_parameters(Number& obj_scaling,
00479                         bool& use_x_scaling, Index n,
00480                         Number* x_scaling,
00481                     bool& use_g_scaling, Index m,
00482                     Number* g_scaling){
00483         if( osinstance->instanceData->objectives->obj[ 0]->maxOrMin.compare("min") != 0){
00484                 obj_scaling = -1;
00485         }
00486     else obj_scaling = 1;
00487     use_x_scaling = false;
00488     use_g_scaling = false;
00489         return true;
00490 }//get_scaling_parameters
00491 
00492 void IpoptProblem::finalize_solution(SolverReturn status,
00493         Index n, const Number* x, const Number* z_L, const Number* z_U,
00494         Index m, const Number* g, const Number* lambda, Number obj_value,
00495         const IpoptData* ip_data, IpoptCalculatedQuantities* ip_cq){
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         if(osinstance->getObjectiveNumber() > 0)
00500         obj_value = osinstance->calculateAllObjectiveFunctionValues( const_cast<double*>(x), true)[ 0];
00501         OSrLWriter *osrlwriter ;
00502         osrlwriter = new OSrLWriter();
00503 #ifdef DEBUG
00504           printf("\n\nSolution of the primal variables, x\n");
00505           for (Index i=0; i<n; i++) {
00506             printf("x[%d] = %e\n", i, x[i]);
00507           }
00508         
00509           printf("\n\nSolution of the bound multipliers, z_L and z_U\n");
00510           for (Index i=0; i<n; i++) {
00511             printf("z_L[%d] = %e\n", i, z_L[i]);
00512           }
00513           for (Index i=0; i<n; i++) {
00514             printf("z_U[%d] = %e\n", i, z_U[i]);
00515           }
00516 #endif
00517         if(osinstance->getObjectiveNumber() > 0) printf("\nObjective value f(x*) = %e\n", obj_value  );
00518           
00519         int solIdx = 0;
00520         int numberOfOtherVariableResults;
00521         int otherIdx;
00522         int numCon = osinstance->getConstraintNumber();
00523         ostringstream outStr;
00524         
00525         
00526         //make sure the sign on the dual values is correct
00527         double *dualValue = NULL;
00528         std::string *rcost = NULL;
00529         int* idx = NULL;
00530         double*  mdObjValues = NULL;
00531         if(osinstance->getObjectiveNumber() > 0){
00532                 mdObjValues = new double[1];
00533         }
00534 
00535         std::string message = "Ipopt solver finishes to the end.";
00536         std::string solutionDescription = "";   
00537 
00538         try{
00539                 // resultHeader infomration
00540                 if(osresult->setSolverInvoked( "COIN-OR Ipopt") != true)
00541                         throw ErrorClass("OSResult error: setSolverInvoked");
00542         if(osresult->setServiceName( getVersionInfo()) != true)
00543                         throw ErrorClass("OSResult error: setServiceName");
00544                 if(osresult->setInstanceName(  osinstance->getInstanceName()) != true)
00545                         throw ErrorClass("OSResult error: setInstanceName");
00546 
00547                 //if(osresult->setJobID( osoption->jobID) != true)
00548                 //      throw ErrorClass("OSResult error: setJobID");
00549 
00550                 // set basic problem parameters
00551                 if(osresult->setVariableNumber( osinstance->getVariableNumber()) != true)
00552                         throw ErrorClass("OSResult error: setVariableNumer");
00553                 if(osresult->setObjectiveNumber( 1) != true)
00554                         throw ErrorClass("OSResult error: setObjectiveNumber");
00555                 if(osresult->setConstraintNumber( osinstance->getConstraintNumber()) != true)
00556                         throw ErrorClass("OSResult error: setConstraintNumber");
00557                 if(osresult->setSolutionNumber(  1) != true)
00558                         throw ErrorClass("OSResult error: setSolutionNumer");   
00559 
00560 
00561                 if(osresult->setGeneralMessage( message) != true)
00562                         throw ErrorClass("OSResult error: setGeneralMessage");
00563 
00564                 switch( status){
00565                         case SUCCESS:
00566                                 solutionDescription = "SUCCESS[IPOPT]: Algorithm terminated normally at a locally optimal point, satisfying the convergence tolerances.";
00567                                 osresult->setSolutionStatus(solIdx,  "locallyOptimal", solutionDescription);
00568 
00569 
00570                                 if(osinstance->getVariableNumber() > 0) osresult->setPrimalVariableValuesDense(solIdx, const_cast<double*>(x)); 
00571 
00572                                 if(numCon > 0){
00573                                         dualValue = new double[ numCon];
00574                                         for (Index i=0; i < numCon; i++) {
00575                                                 dualValue[ i] = -lambda[ i];
00576                                         }
00577                                         //osresult->setDualVariableValuesDense(solIdx, const_cast<double*>(lambda)); 
00578                                         osresult->setDualVariableValuesDense(solIdx, dualValue); 
00579                                 }
00580 
00581                                 if(osinstance->getObjectiveNumber() > 0){
00582                                         mdObjValues[0] = obj_value ;
00583                                         osresult->setObjectiveValuesDense(solIdx, mdObjValues); 
00584                                 }
00585                                 // set other
00586 
00587                                 if(osinstance->getVariableNumber() > 0){
00588                                         numberOfOtherVariableResults = 1;
00589                                         osresult->setNumberOfOtherVariableResults(solIdx, numberOfOtherVariableResults);
00590                                         rcost = new std::string[ osinstance->getVariableNumber()];
00591                                         idx = new int[ osinstance->getVariableNumber()];
00592                                         for (Index i = 0; i < n; i++) {
00593                                         rcost[ i] =  os_dtoa_format( z_L[i] - z_U[i]); 
00594                                                 idx[ i] = i;
00595                                         }
00596                                         otherIdx = 0;
00597                                         osresult->setAnOtherVariableResultSparse(solIdx, otherIdx, "reduced costs", "", "the variable reduced costs", 
00598                                         idx, rcost, osinstance->getVariableNumber());   
00599                                 }
00600                                 //set dual values on variable upper and lower bounds
00601                                 /*
00602                                 for (Index i = 0; i < n; i++) {
00603                                     rcost[ i] =  os_dtoa_format( z_L[i]); 
00604                                         idx[ i] = i;
00605                                  }
00606                                 otherIdx = 0;
00607                                 osresult->setAnOtherVariableResultSparse(solIdx, otherIdx, "varL", "", "Lagrange Multiplier on the Variable Lower Bound",  idx,  rcost,  osinstance->getVariableNumber() );
00608                                 
00609                                 
00610                                 for (Index i = 0; i < n; i++) {
00611                                     rcost[ i] =  os_dtoa_format( z_U[i]); 
00612                                  }
00613                                 otherIdx = 1;
00614                                 osresult->setAnOtherVariableResultSparse(solIdx, otherIdx, "varU", "", "Lagrange Multiplier on the Variable Upper Bound", idx,  rcost, osinstance->getVariableNumber() );
00615                                 */
00616                                 if(osinstance->getVariableNumber() > 0) delete[] rcost;
00617                                 if(osinstance->getVariableNumber() > 0) delete[] idx;
00618                                 if(osinstance->getConstraintNumber() > 0) delete[] dualValue;
00619                                 
00620                                 // done with other
00621                                 
00622                                 
00623                         break;
00624                         case MAXITER_EXCEEDED:
00625                                 solutionDescription = "MAXITER_EXCEEDED[IPOPT]: Maximum number of iterations exceeded.";
00626                                 osresult->setSolutionStatus(solIdx,  "stoppedByLimit", solutionDescription);
00627                                 if(x != NULL) osresult->setPrimalVariableValuesDense(solIdx, const_cast<double*>(x)); 
00628                                 if(lambda != NULL) osresult->setDualVariableValuesDense(solIdx, const_cast<double*>( lambda)); 
00629                                 if(osinstance->getObjectiveNumber() > 0){
00630                                         mdObjValues[0] = obj_value ;
00631                                         osresult->setObjectiveValuesDense(solIdx, mdObjValues); 
00632                                 }
00633                         break;
00634                         case STOP_AT_TINY_STEP:
00635                                 solutionDescription = "STOP_AT_TINY_STEP[IPOPT]: Algorithm proceeds with very little progress.";
00636                                 osresult->setSolutionStatus(solIdx,  "stoppedByLimit", solutionDescription);
00637                                 if(x != NULL) osresult->setPrimalVariableValuesDense(solIdx, const_cast<double*>( x)); 
00638                                 if(lambda != NULL) osresult->setDualVariableValuesDense(solIdx, const_cast<double*>( lambda));
00639                                 if(osinstance->getObjectiveNumber() > 0){ 
00640                                         mdObjValues[0] = obj_value ;
00641                                         osresult->setObjectiveValuesDense(solIdx, mdObjValues); 
00642                                 }
00643                         break;
00644                         case STOP_AT_ACCEPTABLE_POINT:
00645                                 solutionDescription = "STOP_AT_ACCEPTABLE_POINT[IPOPT]: Algorithm stopped at a point that was converged, not to _desired_ tolerances, but to _acceptable_ tolerances";
00646                                 osresult->setSolutionStatus(solIdx,  "IpoptAccetable", solutionDescription);
00647                                 if(lambda != NULL) osresult->setDualVariableValuesDense(solIdx, const_cast<double*>( lambda));
00648                                 if(x != NULL)osresult->setPrimalVariableValuesDense(solIdx, const_cast<double*>(x)); 
00649                                 if(osinstance->getObjectiveNumber() > 0){ 
00650                                         mdObjValues[0] = obj_value ;
00651                                         osresult->setObjectiveValuesDense(solIdx, mdObjValues); 
00652                                 }
00653                         break;
00654                         case LOCAL_INFEASIBILITY:
00655                                 solutionDescription = "LOCAL_INFEASIBILITY[IPOPT]: Algorithm converged to a point of local infeasibility. Problem may be infeasible.";
00656                                 osresult->setSolutionStatus(solIdx,  "infeasible", solutionDescription);
00657                                 if( osinstance->getVariableNumber() == 0) osresult->setSolutionMessage(solIdx, "Warning: this problem has zero decision variables!");
00658                         break;
00659                         case USER_REQUESTED_STOP:
00660                                 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.";
00661                                 osresult->setSolutionStatus(solIdx,  "error", solutionDescription);
00662                                 if( osinstance->getVariableNumber() == 0) osresult->setSolutionMessage(solIdx, "Warning: this problem has zero decision variables!");
00663                         break;
00664                         case DIVERGING_ITERATES:
00665                                 solutionDescription = "DIVERGING_ITERATES[IPOPT]: It seems that the iterates diverge.";
00666                                 osresult->setSolutionStatus(solIdx,  "unbounded", solutionDescription);
00667                                 if( osinstance->getVariableNumber() == 0) osresult->setSolutionMessage(solIdx, "Warning: this problem has zero decision variables!");
00668                         break;
00669                         case RESTORATION_FAILURE:
00670                                 solutionDescription = "RESTORATION_FAILURE[IPOPT]: Restoration phase failed, algorithm doesn't know how to proceed.";
00671                                 osresult->setSolutionStatus(solIdx,  "error", solutionDescription);
00672                                 if( osinstance->getVariableNumber() == 0) osresult->setSolutionMessage(solIdx, "Warning: this problem has zero decision variables!");
00673                         break;
00674                         case ERROR_IN_STEP_COMPUTATION:
00675                                 solutionDescription = "ERROR_IN_STEP_COMPUTATION[IPOPT]: An unrecoverable error occurred while IPOPT tried to compute the search direction.";
00676                                 osresult->setSolutionStatus(solIdx,  "error", solutionDescription);
00677                                 if( osinstance->getVariableNumber() == 0) osresult->setSolutionMessage(solIdx, "Warning: this problem has zero decision variables!");
00678                         break;
00679                         case INVALID_NUMBER_DETECTED:
00680                                 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.";
00681                                 osresult->setSolutionStatus(solIdx,  "error", solutionDescription);
00682                                 if( osinstance->getVariableNumber() == 0) osresult->setSolutionMessage(solIdx, "Warning: this problem has zero decision variables!");
00683                         break;
00684                         case INTERNAL_ERROR:
00685                                 solutionDescription = "INTERNAL_ERROR[IPOPT]: An unknown internal error occurred. Please contact the IPOPT authors through the mailing list.";
00686                                 osresult->setSolutionStatus(solIdx,  "error", solutionDescription);
00687                                 if( osinstance->getVariableNumber() == 0) osresult->setSolutionMessage(solIdx, "Warning: this problem has zero decision variables!");
00688                         break;
00689                         default:
00690                                 solutionDescription = "OTHER[IPOPT]: other unknown solution status from Ipopt solver";
00691                                 osresult->setSolutionStatus(solIdx,  "other", solutionDescription);
00692                                 if( osinstance->getVariableNumber() == 0) osresult->setSolutionMessage(solIdx, "Warning: this problem has zero decision variables!");
00693                 }
00694                 osresult->setGeneralStatusType("normal");
00695                 delete osrlwriter;
00696                 if(osinstance->getObjectiveNumber() > 0){
00697                         delete[] mdObjValues;
00698                 }
00699                 osrlwriter = NULL;
00700 
00701         }
00702         catch(const ErrorClass& eclass){
00703 #ifdef DEBUG
00704                 cout << "error in OSIpoptSolver, line 636:\n" << eclass.errormsg << endl;
00705 #endif
00706                 osresult->setGeneralMessage( eclass.errormsg);
00707                 osresult->setGeneralStatusType( "error");
00708                 std::string osrl = osrlwriter->writeOSrL( osresult);
00709                 delete osrlwriter;
00710                 osrlwriter = NULL;
00711                 throw ErrorClass(  osrl) ;
00712                 if(osinstance->getObjectiveNumber() > 0){
00713                         delete[] mdObjValues;
00714                 }
00715                 mdObjValues = NULL;
00716         }
00718 }
00719 
00720 
00721 void IpoptSolver::setSolverOptions() throw (ErrorClass) {
00722         try{
00723                 //if(osinstance->getObjectiveNumber() <= 0) throw ErrorClass("Ipopt NEEDS AN OBJECTIVE FUNCTION");
00724                 this->bSetSolverOptions = true;
00725                 /* set the default options */   
00726                 //app->Options()->SetNumericValue("tol", 1e-9);
00727                 app->Options()->SetIntegerValue("print_level", 0);
00728                 app->Options()->SetIntegerValue("max_iter", 20000);
00729                 app->Options()->SetNumericValue("bound_relax_factor", 0, true, true);
00730                 app->Options()->SetStringValue("mu_strategy", "adaptive", true, true);
00731                 //app->Options()->SetStringValue("output_file", "ipopt.out");
00732                 app->Options()->SetStringValue("check_derivatives_for_naninf", "yes");
00733                 // hessian constant for an LP
00734                 if( (osinstance-> getNumberOfNonlinearExpressions() <= 0) && (osinstance->getNumberOfQuadraticTerms() <= 0) ){
00735                         app->Options()->SetStringValue("hessian_constant", "yes", true, true);
00736                 }
00737                 if(osinstance->getObjectiveNumber() > 0){
00738                         if( osinstance->instanceData->objectives->obj[ 0]->maxOrMin.compare("min") != 0){
00739                                 app->Options()->SetStringValue("nlp_scaling_method", "user-scaling");
00740                         }
00741                 }
00742                 /* end of the default options, now get options from OSoL */
00743         
00744                 
00745                 if(osoption == NULL && osol.length() > 0)
00746                 {
00747                         m_osolreader = new OSoLReader();
00748                         osoption = m_osolreader->readOSoL( osol);
00749                 }
00750 
00751                 if( osoption != NULL  &&  osoption->getNumberOfSolverOptions() > 0 ){
00752 //                      std::cout << "number of solver options "  <<  osoption->getNumberOfSolverOptions() << std::endl;
00753                         std::vector<SolverOption*> optionsVector;
00754                         optionsVector = osoption->getSolverOptions( "ipopt");
00755                         char *pEnd;
00756                         int i;
00757                         int num_ipopt_options = optionsVector.size();
00758                         for(i = 0; i < num_ipopt_options; i++){
00759 //                              std::cout << "ipopt solver option  "  << optionsVector[ i]->name << std::endl;
00760                                 if(optionsVector[ i]->type == "numeric" ){
00761 //                                      std::cout << "FOUND A NUMERIC OPTION  "  <<  os_strtod( optionsVector[ i]->value.c_str(), &pEnd ) << std::endl;
00762                                         app->Options()->SetNumericValue(optionsVector[ i]->name, os_strtod( optionsVector[ i]->value.c_str(), &pEnd ) );        
00763                                 }
00764                                 else if(optionsVector[ i]->type == "integer" ){
00765 //                                      std::cout << "FOUND AN INTEGER OPTION  "  << atoi( optionsVector[ i]->value.c_str() ) << std::endl;
00766                                         app->Options()->SetIntegerValue(optionsVector[ i]->name, atoi( optionsVector[ i]->value.c_str() ) );
00767                                 }
00768                                 else if(optionsVector[ i]->type == "string" ){
00769 //                                      std::cout << "FOUND A STRING OPTION  "  <<  optionsVector[ i]->value.c_str() << std::endl;
00770                                         app->Options()->SetStringValue(optionsVector[ i]->name, optionsVector[ i]->value);
00771                                 }
00772                         }       
00773                 }
00774         }
00775         catch(const ErrorClass& eclass){
00776 #ifdef DEBUG
00777                 cout << "error in OSIpoptSolver, line 695:\n" << eclass.errormsg << endl;
00778 #endif
00779                 std::cout << "THERE IS AN ERROR" << std::endl;
00780                 osresult->setGeneralMessage( eclass.errormsg);
00781                 osresult->setGeneralStatusType( "error");
00782                 osrl = osrlwriter->writeOSrL( osresult);
00783                 throw ErrorClass( osrl) ;
00784         }                               
00785 }//end setSolverOptions() 
00786 
00787 
00788 
00789 void IpoptSolver::buildSolverInstance() throw (ErrorClass) {
00790         try{
00791                 
00792                 if(osil.length() == 0 && osinstance == NULL) throw ErrorClass("there is no instance");
00793                 if(osinstance == NULL){
00794                         m_osilreader = new OSiLReader();
00795                         osinstance = m_osilreader->readOSiL( osil);
00796                 }
00797                 // Create a new instance of your nlp 
00798                 nlp = new IpoptProblem( osinstance, osoption, osresult, ipoptErrorMsg);
00799                 app = new IpoptApplication();
00800                 this->bCallbuildSolverInstance = true;
00801         }
00802         catch(const ErrorClass& eclass){
00803 #ifdef DEBUG
00804                 cout << "error in OSIpoptSolver, line 722:\n" << eclass.errormsg << endl;
00805 #endif
00806                 std::cout << "THERE IS AN ERROR" << std::endl;
00807                 osresult->setGeneralMessage( eclass.errormsg);
00808                 osresult->setGeneralStatusType( "error");
00809                 osrl = osrlwriter->writeOSrL( osresult);
00810                 throw ErrorClass( osrl) ;
00811         }                               
00812 }//end buildSolverInstance() 
00813 
00814 
00815 
00816 void IpoptSolver::solve() throw (ErrorClass) {
00817         if( this->bCallbuildSolverInstance == false) buildSolverInstance();
00818         if( this->bSetSolverOptions == false) setSolverOptions();
00819         try{
00820                 clock_t start, finish;
00821                 double duration;
00822                 start = clock();
00823                 //OSiLWriter osilwriter;
00824                 //cout << osilwriter.writeOSiL( osinstance) << endl;
00825                 //if(osinstance->getVariableNumber() <= 0)throw ErrorClass("Ipopt requires decision variables");
00826                 finish = clock();
00827                 duration = (double) (finish - start) / CLOCKS_PER_SEC;
00828                 //cout << "Parsing took (seconds): " << duration << endl; 
00829                 //dataEchoCheck();
00830                 /***************now the ipopt invokation*********************/
00831                 // see if we have a linear program
00832                 //if(osinstance->getObjectiveNumber() <= 0) throw ErrorClass("Ipopt NEEDS AN OBJECTIVE FUNCTION");
00833                 // Intialize the IpoptApplication and process the options
00834 //              std::cout << "Call Ipopt Initialize" << std::endl;
00835                 app->Initialize();
00836 //              std::cout << "Finished Ipopt Initialize" << std::endl;
00837                 //nlp->osinstance = this->osinstance;
00838                 // Ask Ipopt to solve the problem
00839 //              std::cout << "Call Ipopt Optimize" << std::endl;
00840                 ApplicationReturnStatus status = app->OptimizeTNLP( nlp);
00841 //              std::cout << "Finish Ipopt Optimize" << std::endl;
00842                 osrl = osrlwriter->writeOSrL( osresult);
00843 //              std::cout << "Finish writing the osrl" << std::endl;
00844                 //if (status != Solve_Succeeded) {
00845                 if (status < -2) {
00846                         throw ErrorClass("Ipopt FAILED TO SOLVE THE PROBLEM: " + *ipoptErrorMsg);
00847                 }       
00848         }
00849         catch(const ErrorClass& eclass){
00850 #ifdef DEBUG
00851                 cout << "error in OSIpoptSolver, line 775:\n" << eclass.errormsg << endl;
00852 #endif
00853                 osresult->setGeneralMessage( eclass.errormsg);
00854                 osresult->setGeneralStatusType( "error");
00855                 osrl = osrlwriter->writeOSrL( osresult);
00856                 throw ErrorClass( osrl) ;
00857         }
00858 }//solve
00859 
00860 
00861 void IpoptSolver::dataEchoCheck(){
00862         
00863         int i;
00864         
00865         // print out problem parameters
00866         cout << "This is problem:  " << osinstance->getInstanceName() << endl;
00867         cout << "The problem source is:  " << osinstance->getInstanceSource() << endl;
00868         cout << "The problem description is:  " << osinstance->getInstanceDescription() << endl;
00869         cout << "number of variables = " << osinstance->getVariableNumber() << endl;
00870         cout << "number of Rows = " << osinstance->getConstraintNumber() << endl;
00871 
00872         // print out the variable information
00873         if(osinstance->getVariableNumber() > 0){
00874                 for(i = 0; i < osinstance->getVariableNumber(); i++){
00875                         if(osinstance->getVariableNames() != NULL) cout << "variable Names  " << osinstance->getVariableNames()[ i]  << endl;
00876                         if(osinstance->getVariableTypes() != NULL) cout << "variable Types  " << osinstance->getVariableTypes()[ i]  << endl;
00877                         if(osinstance->getVariableLowerBounds() != NULL) cout << "variable Lower Bounds  " << osinstance->getVariableLowerBounds()[ i]  << endl;
00878                         if(osinstance->getVariableUpperBounds() != NULL) cout << "variable Upper Bounds  " <<  osinstance->getVariableUpperBounds()[i] << endl;
00879                 }
00880         }
00881         
00882         // print out objective function information
00883         if(osinstance->getVariableNumber() > 0 || osinstance->instanceData->objectives->obj != NULL || osinstance->instanceData->objectives->numberOfObjectives > 0){
00884                 if( osinstance->getObjectiveMaxOrMins()[0] == "min")  cout <<  "problem is a minimization" << endl;
00885                 else cout <<  "problem is a maximization" << endl;
00886                 for(i = 0; i < osinstance->getVariableNumber(); i++){
00887                         cout << "OBJ COEFFICIENT =  " <<  osinstance->getDenseObjectiveCoefficients()[0][i] << endl;
00888                 }
00889         }
00890         // print out constraint information
00891         if(osinstance->getConstraintNumber() > 0){
00892                 for(i = 0; i < osinstance->getConstraintNumber(); i++){
00893                         if(osinstance->getConstraintNames() != NULL) cout << "row name = " << osinstance->getConstraintNames()[i] <<  endl;
00894                         if(osinstance->getConstraintLowerBounds() != NULL) cout << "row lower bound = " << osinstance->getConstraintLowerBounds()[i] <<  endl;
00895                         if(osinstance->getConstraintUpperBounds() != NULL) cout << "row upper bound = " << osinstance->getConstraintUpperBounds()[i] <<  endl; 
00896                 }
00897         }
00898         
00899         // print out linear constraint data
00900         cout << endl;
00901         cout << "number of nonzeros =  " << osinstance->getLinearConstraintCoefficientNumber() << endl;
00902         for(i = 0; i <= osinstance->getVariableNumber(); i++){
00903                 cout << "Start Value =  " << osinstance->getLinearConstraintCoefficientsInColumnMajor()->starts[ i] << endl;
00904         }
00905         cout << endl;
00906         for(i = 0; i < osinstance->getLinearConstraintCoefficientNumber(); i++){
00907                 cout << "Index Value =  " << osinstance->getLinearConstraintCoefficientsInColumnMajor()->indexes[i] << endl;
00908                 cout << "Nonzero Value =  " << osinstance->getLinearConstraintCoefficientsInColumnMajor()->values[i] << endl;
00909         }
00910 
00911         // print out quadratic data
00912         cout << "number of qterms =  " <<  osinstance->getNumberOfQuadraticTerms() << endl;
00913         for(int i = 0; i <  osinstance->getNumberOfQuadraticTerms(); i++){
00914                 cout << "Row Index = " <<  osinstance->getQuadraticTerms()->rowIndexes[i] << endl;
00915                 cout << "Var Index 1 = " << osinstance->getQuadraticTerms()->varOneIndexes[ i] << endl;
00916                 cout << "Var Index 2 = " << osinstance->getQuadraticTerms()->varTwoIndexes[ i] << endl;
00917                 cout << "Coefficient = " << osinstance->getQuadraticTerms()->coefficients[ i] << endl;
00918         }
00919 } // end dataEchoCheck
00920 
00921 
00922 IpoptProblem::IpoptProblem(OSInstance *osinstance_,  OSOption *osoption_, OSResult *osresult_, std::string* ipoptErrorMsg_) {
00923         osinstance = osinstance_;
00924         osoption = osoption_;
00925         osresult = osresult_;
00926         ipoptErrorMsg = ipoptErrorMsg_;
00927 }
00928 
00929 IpoptProblem::~IpoptProblem() {
00930 
00931 }
00932 
00933 
00934 

Generated on Thu Aug 5 03:03:04 2010 by  doxygen 1.4.7