/home/coin/SVN-release/OS-2.0.1/OS/v2.0/OSIpoptSolver.cpp

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

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