/Users/kmartin/Documents/files/code/cpp/OScpp/COIN-OS/OS/src/OSSolverInterfaces/OSIpoptSolver.cpp

Go to the documentation of this file.
00001 
00020 #include "OSIpoptSolver.h"
00021 #include "OSCommonUtil.h"
00022 
00023 
00024 using std::cout; 
00025 using std::endl; 
00026 using std::ostringstream;
00027 using namespace Ipopt;
00028 
00029 
00030 IpoptSolver::IpoptSolver() {
00031         osrlwriter = new OSrLWriter();
00032         osresult = new OSResult();
00033         m_osilreader = NULL;
00034         ipoptErrorMsg = "";
00035 
00036 } 
00037 
00038 IpoptSolver::~IpoptSolver() {
00039         #ifdef DEBUG
00040         cout << "inside IpoptSolver destructor" << endl;
00041         #endif
00042         if(m_osilreader != NULL) delete m_osilreader;
00043         m_osilreader = NULL;
00044         delete osresult;
00045         osresult = NULL;
00046         delete osrlwriter;
00047         osrlwriter = NULL;
00048         //delete osinstance;
00049         //osinstance = NULL;
00050         #ifdef DEBUG
00051         cout << "leaving IpoptSolver destructor" << endl;
00052         #endif
00053 }
00054 
00055 // returns the size of the problem
00056 bool IpoptProblem::get_nlp_info(Index& n, Index& m, Index& nnz_jac_g,
00057                              Index& nnz_h_lag, IndexStyleEnum& index_style)
00058 {
00059         if(osinstance->getObjectiveNumber() <= 0) throw ErrorClass("Ipopt NEEDS AN OBJECTIVE FUNCTION");                        
00060         // number of variables
00061         n = osinstance->getVariableNumber();
00062         // number of constraints
00063         m = osinstance->getConstraintNumber();
00064         cout << "number variables  !!!!!!!!!!!!!!!!!!!!!!!!!!!" << n << endl;
00065         cout << "number constraints  !!!!!!!!!!!!!!!!!!!!!!!!!!!" << m << endl;
00066         try{
00067                 osinstance->initForAlgDiff( );
00068         }
00069         catch(const ErrorClass& eclass){
00070                 ipoptErrorMsg = eclass.errormsg;
00071                 throw;  
00072         }       
00073         // use the OS Expression tree for function evaluations instead of CppAD
00074         osinstance->bUseExpTreeForFunEval = true;
00075         //std::cout << "Call sparse jacobian" << std::endl;
00076         SparseJacobianMatrix *sparseJacobian = NULL;
00077         try{
00078                 sparseJacobian = osinstance->getJacobianSparsityPattern();
00079         }
00080         catch(const ErrorClass& eclass){
00081                 ipoptErrorMsg = eclass.errormsg;
00082                 throw;  
00083         }
00084         //std::cout << "Done calling sparse jacobian" << std::endl;
00085         nnz_jac_g = sparseJacobian->valueSize;
00086         cout << "nnz_jac_g  !!!!!!!!!!!!!!!!!!!!!!!!!!!" << nnz_jac_g << endl;  
00087         // nonzeros in upper hessian
00088         
00089         if( (osinstance->getNumberOfNonlinearExpressions() == 0) && (osinstance->getNumberOfQuadraticTerms() == 0) ) {
00090                 cout << "This is a linear program"  << endl;
00091                 nnz_h_lag = 0;
00092         }
00093         else{
00094                 //std::cout << "Get Lagrangain Hessian Sparsity Pattern " << std::endl;
00095                 SparseHessianMatrix *sparseHessian = osinstance->getLagrangianHessianSparsityPattern();
00096                 //std::cout << "Done Getting Lagrangain Hessian Sparsity Pattern " << std::endl;
00097                 nnz_h_lag = sparseHessian->hessDimension;
00098         }
00099         cout << "nnz_h_lag  !!!!!!!!!!!!!!!!!!!!!!!!!!!" << nnz_h_lag << endl;  
00100         // use the C style indexing (0-based)
00101         index_style = TNLP::C_STYLE;
00102   
00104 
00105   return true;
00106 }//get_nlp_info
00107 
00108 
00109 bool  IpoptProblem::get_bounds_info(Index n, Number* x_l, Number* x_u,
00110                                 Index m, Number* g_l, Number* g_u){
00111         int i; 
00112         double * mdVarLB = osinstance->getVariableLowerBounds();
00113         //std::cout << "GET BOUNDS INFORMATION FOR IPOPT !!!!!!!!!!!!!!!!!" << std::endl;
00114         // variables upper bounds
00115         double * mdVarUB = osinstance->getVariableUpperBounds();
00116 
00117         for(i = 0; i < n; i++){
00118                 x_l[ i] = mdVarLB[ i];
00119                 x_u[ i] = mdVarUB[ i];
00120                 //cout << "x_l !!!!!!!!!!!!!!!!!!!!!!!!!!!" << x_l[i] << endl;
00121                 //cout << "x_u !!!!!!!!!!!!!!!!!!!!!!!!!!!" << x_u[i] << endl;
00122         }
00123         // Ipopt interprets any number greater than nlp_upper_bound_inf as
00124         // infinity. The default value of nlp_upper_bound_inf and nlp_lower_bound_inf
00125         // is 1e19 and can be changed through ipopt options.
00126         // e.g. g_u[0] = 2e19;
00127 
00128         //constraint lower bounds
00129         double * mdConLB = osinstance->getConstraintLowerBounds();
00130         //constraint upper bounds
00131         double * mdConUB = osinstance->getConstraintUpperBounds();
00132         
00133         for(int i = 0; i < m; i++){
00134                 g_l[ i] = mdConLB[ i];
00135                 g_u[ i] = mdConUB[ i];
00136                 //cout << "lower !!!!!!!!!!!!!!!!!!!!!!!!!!!" << g_l[i] << endl;
00137                 //cout << "upper !!!!!!!!!!!!!!!!!!!!!!!!!!!" << g_u[i] << endl;
00138         }  
00139         return true;
00140 }//get_bounds_info
00141 
00142 
00143 // returns the initial point for the problem
00144 bool IpoptProblem::get_starting_point(Index n, bool init_x, Number* x,
00145      bool init_z, Number* z_L, Number* z_U, Index m, bool init_lambda,
00146      Number* lambda) {
00147         // Here, we assume we only have starting values for x, if you code
00148         // your own NLP, you can provide starting values for the dual variables
00149         // if you wish
00150         assert(init_x == true);
00151         assert(init_z == false);
00152         assert(init_lambda == false);
00153         int i;
00154         //cout << "get initial values !!!!!!!!!!!!!!!!!!!!!!!!!! " << endl;
00155         double *mdXInit = osinstance->getVariableInitialValues(); 
00156         if( mdXInit != NULL) {
00157                 for(i = 0; i < n; i++){
00158                         if( CommonUtil::ISOSNAN( mdXInit[ i]) == true){ 
00159                                 x[ i] = 1.7171; 
00160                                 //std::cout << "INITIAL VALUE !!!!!!!!!!!!!!!!!!!!  " << x[ i] << std::endl;
00161                         }
00162                         else x[ i] = mdXInit[ i];
00163                         //std::cout << "INITIAL VALUE !!!!!!!!!!!!!!!!!!!!  " << x[ i] << std::endl;    
00164                 }       
00165         }
00166         else{
00167                 for(i = 0; i < n; i++){
00168                         x[ i] = 1.7171;
00169                 }
00170         }
00171         //cout << "got initial values !!!!!!!!!!!!!!!!!!!!!!!!!! " << endl;
00172         return true;
00173 }//get_starting_point
00174 
00175 // returns the value of the objective function
00176 bool IpoptProblem::eval_f(Index n, const Number* x, bool new_x, Number& obj_value){
00177         try{
00178                 obj_value  = osinstance->calculateAllObjectiveFunctionValues( const_cast<double*>(x), NULL, NULL, new_x, 0 )[ 0];
00179         }
00180         catch(const ErrorClass& eclass){
00181                 ipoptErrorMsg = eclass.errormsg;
00182                 throw;  
00183         }
00184         if( CommonUtil::ISOSNAN( (double)obj_value) ) return false;
00185         return true;
00186 }
00187 
00188 bool IpoptProblem::eval_grad_f(Index n, const Number* x, bool new_x, Number* grad_f){
00189         int i;
00190         double *objGrad;
00191         try{
00192                 //objGrad = osinstance->calculateAllObjectiveFunctionGradients( const_cast<double*>(x), NULL, NULL,  new_x, 1)[ 0];
00193                 objGrad = osinstance->calculateObjectiveFunctionGradient( const_cast<double*>(x), NULL, NULL, -1,  new_x, 1);
00194         }
00195         catch(const ErrorClass& eclass){
00196                 ipoptErrorMsg = eclass.errormsg;
00197                 throw;  
00198         }
00199         for(i = 0; i < n; i++){
00200                 grad_f[ i]  = objGrad[ i];
00201         }
00202         return true;
00203 }//eval_grad_f
00204 
00205 // return the value of the constraints: g(x)
00206 bool IpoptProblem::eval_g(Index n, const Number* x, bool new_x, Index m, Number* g) {
00207         try{
00208                 double *conVals = osinstance->calculateAllConstraintFunctionValues( const_cast<double*>(x), NULL, NULL, new_x, 0 );
00209                 int i;
00210                 for(i = 0; i < m; i++){
00211                         if( CommonUtil::ISOSNAN( (double)conVals[ i] ) ) return false;
00212                         g[i] = conVals[ i]  ;   
00213                 } 
00214                 return true;
00215         }
00216         catch(const ErrorClass& eclass){
00217                 ipoptErrorMsg = eclass.errormsg;
00218                 throw;  
00219         }
00220 }//eval_g
00221 
00222 
00223 // return the structure or values of the jacobian
00224 bool IpoptProblem::eval_jac_g(Index n, const Number* x, bool new_x,
00225                            Index m, Index nele_jac, Index* iRow, Index *jCol,
00226                            Number* values){
00227         SparseJacobianMatrix *sparseJacobian;
00228         if (values == NULL) {
00229                 // return the values of the jacobian of the constraints
00230                 //cout << "n: " << n << endl;
00231                 //cout << "m: " << m << endl;
00232                 //cout << "nele_jac: " <<  nele_jac << endl;
00233                 // return the structure of the jacobian
00234                 try{
00235                         sparseJacobian = osinstance->getJacobianSparsityPattern();
00236                 }
00237                 catch(const ErrorClass& eclass){
00238                         ipoptErrorMsg =  eclass.errormsg; 
00239                         throw; 
00240                 }
00241                 int i = 0;
00242                 int k, idx;
00243                 for(idx = 0; idx < m; idx++){
00244                         for(k = *(sparseJacobian->starts + idx); k < *(sparseJacobian->starts + idx + 1); k++){
00245                                 iRow[i] = idx;
00246                                 jCol[i] = *(sparseJacobian->indexes + k);
00247                                 //cout << "ROW IDX  !!!!!!!!!!!!!!!!!!!!!!!!!!!"  << iRow[i] << endl;
00248                                 //cout << "COL IDX  !!!!!!!!!!!!!!!!!!!!!!!!!!!"  << jCol[i] << endl;
00249                                 i++;
00250                         }
00251                 }
00252         }
00253         else {
00254                 //std::cout << "EVALUATING JACOBIAN" << std::endl; 
00255                 try{
00256                         sparseJacobian = osinstance->calculateAllConstraintFunctionGradients( const_cast<double*>(x), NULL, NULL,  new_x, 1);
00257                 }
00258                 catch(const ErrorClass& eclass){
00259                         ipoptErrorMsg = eclass.errormsg;
00260                         throw;  
00261                 }
00262                 //osinstance->getIterateResults( (double*)x, 0.0, NULL, -1, new_x,  1);
00263                 for(int i = 0; i < nele_jac; i++){
00264                         values[ i] = sparseJacobian->values[i];
00265                         //values[ i] = osinstance->m_mdJacValue[ i];
00266                         //cout << "values[i]:!!!!!!!!!!!!  " <<  values[ i] << endl;
00267                         //cout << "m_mdJacValue[ i]:!!!!!!!!!!!!  " <<  osinstance->m_mdJacValue[ i] << endl;
00268                 }
00269         }
00270   return true;
00271 }//eval_jac_g
00272 
00273 //return the structure or values of the hessian
00274 bool IpoptProblem::eval_h(Index n, const Number* x, bool new_x,
00275                        Number obj_factor, Index m, const Number* lambda,
00276                        bool new_lambda, Index nele_hess, Index* iRow,
00277                        Index* jCol, Number* values){
00278 
00280         SparseHessianMatrix *sparseHessian;
00281         
00282         int i;
00283         if (values == NULL) {
00284                 // return the structure. This is a symmetric matrix, fill the lower left triangle only.
00285                 //cout << "get structure of HESSIAN !!!!!!!!!!!!!!!!!!!!!!!!!! "  << endl;
00286                 try{
00287                         sparseHessian = osinstance->getLagrangianHessianSparsityPattern( );
00288                 }
00289                 catch(const ErrorClass& eclass){
00290                         ipoptErrorMsg = eclass.errormsg;
00291                         throw;  
00292                 }
00293                 //cout << "got structure of HESSIAN !!!!!!!!!!!!!!!!!!!!!!!!!! "  << endl;
00294                 for(i = 0; i < nele_hess; i++){
00295                         iRow[i] = *(sparseHessian->hessColIdx + i);
00296                         jCol[i] = *(sparseHessian->hessRowIdx + i);
00297                                 //cout << "ROW HESS IDX  !!!!!!!!!!!!!!!!!!!!!!!!!!!"  << iRow[i] << endl;
00298                                 //cout << "COL HESS IDX  !!!!!!!!!!!!!!!!!!!!!!!!!!!"  << jCol[i] << endl;
00299                 }
00300         }
00301         else {
00302                 //std::cout << "EVALUATING HESSIAN" << std::endl; 
00303                 // return the values. This is a symmetric matrix, fill the lower left triangle only
00304                 double* objMultipliers = new double[1];
00305                 objMultipliers[0] = obj_factor;
00306                 try{
00307                         sparseHessian = osinstance->calculateLagrangianHessian( const_cast<double*>(x), objMultipliers, (double*)lambda ,  new_x, 2);
00308                 delete[]  objMultipliers;
00309                 }
00310                 catch(const ErrorClass& eclass){
00311                         ipoptErrorMsg = eclass.errormsg;
00312                         delete[]  objMultipliers;
00313                         throw;  
00314                 }
00315                 for(i = 0; i < nele_hess; i++){
00316                         values[ i]  = *(sparseHessian->hessValues + i);
00317                 }
00318         }
00320         return true;
00321 }//eval_h
00322 
00323 bool IpoptProblem::get_scaling_parameters(Number& obj_scaling,
00324                         bool& use_x_scaling, Index n,
00325                         Number* x_scaling,
00326                     bool& use_g_scaling, Index m,
00327                     Number* g_scaling){
00328         if( osinstance->instanceData->objectives->obj[ 0]->maxOrMin.compare("min") != 0){
00329                 obj_scaling = -1;
00330         }
00331     else obj_scaling = 1;
00332     use_x_scaling = false;
00333     use_g_scaling = false;
00334         return true;
00335 }//get_scaling_parameters
00336 
00337 void IpoptProblem::finalize_solution(SolverReturn status,
00338                                Index n, const Number* x, const Number* z_L, const Number* z_U,
00339                                   Index m, const Number* g, const Number* lambda,
00340                                   Number obj_value,
00341                                    const IpoptData* ip_data,
00342                                    IpoptCalculatedQuantities* ip_cq)
00343 {
00344           // here is where we would store the solution to variables, or write to a file, etc
00345           // so we could use the solution.
00346           // For this example, we write the solution to the console
00347         OSrLWriter *osrlwriter ;
00348         osrlwriter = new OSrLWriter();
00349 #ifdef DEBUG
00350           printf("\n\nSolution of the primal variables, x\n");
00351           for (Index i=0; i<n; i++) {
00352             printf("x[%d] = %e\n", i, x[i]);
00353           }
00354         
00355           printf("\n\nSolution of the bound multipliers, z_L and z_U\n");
00356           for (Index i=0; i<n; i++) {
00357             printf("z_L[%d] = %e\n", i, z_L[i]);
00358           }
00359           for (Index i=0; i<n; i++) {
00360             printf("z_U[%d] = %e\n", i, z_U[i]);
00361           }
00362 #endif
00363           printf("\n\nObjective value\n");
00364           printf("f(x*) = %e\n", obj_value);
00365         int solIdx = 0;
00366         ostringstream outStr;
00367         double* mdObjValues = new double[1];
00368         std::string message = "Ipopt solver finishes to the end.";
00369         std::string solutionDescription = "";   
00370 
00371         try{
00372                 // resultHeader infomration
00373                 if(osresult->setServiceName( "Ipopt solver service") != true)
00374                         throw ErrorClass("OSResult error: setServiceName");
00375                 if(osresult->setInstanceName(  osinstance->getInstanceName()) != true)
00376                         throw ErrorClass("OSResult error: setInstanceName");
00377 
00378                 //if(osresult->setJobID( osoption->jobID) != true)
00379                 //      throw ErrorClass("OSResult error: setJobID");
00380 
00381                 // set basic problem parameters
00382                 if(osresult->setVariableNumber( osinstance->getVariableNumber()) != true)
00383                         throw ErrorClass("OSResult error: setVariableNumer");
00384                 if(osresult->setObjectiveNumber( 1) != true)
00385                         throw ErrorClass("OSResult error: setObjectiveNumber");
00386                 if(osresult->setConstraintNumber( osinstance->getConstraintNumber()) != true)
00387                         throw ErrorClass("OSResult error: setConstraintNumber");
00388                 if(osresult->setSolutionNumber(  1) != true)
00389                         throw ErrorClass("OSResult error: setSolutionNumer");   
00390 
00391 
00392                 if(osresult->setGeneralMessage( message) != true)
00393                         throw ErrorClass("OSResult error: setGeneralMessage");
00394 
00395                 switch( status){
00396                         case SUCCESS:
00397                                 solutionDescription = "SUCCESS[IPOPT]: Algorithm terminated successfully at a locally optimal point, satisfying the convergence tolerances.";
00398                                 osresult->setSolutionStatus(solIdx,  "locallyOptimal", solutionDescription);
00399                                 osresult->setPrimalVariableValues(solIdx, const_cast<double*>(x));
00400                                 osresult->setDualVariableValues(solIdx, const_cast<double*>( lambda));
00401                                 mdObjValues[0] = obj_value;
00402                                 osresult->setObjectiveValues(solIdx, mdObjValues);
00403                         break;
00404                         case MAXITER_EXCEEDED:
00405                                 solutionDescription = "MAXITER_EXCEEDED[IPOPT]: Maximum number of iterations exceeded.";
00406                                 osresult->setSolutionStatus(solIdx,  "stoppedByLimit", solutionDescription);
00407                                 osresult->setPrimalVariableValues(solIdx, const_cast<double*>(x));
00408                                 osresult->setDualVariableValues(solIdx, const_cast<double*>( lambda));
00409                                 mdObjValues[0] = obj_value;
00410                                 osresult->setObjectiveValues(solIdx, mdObjValues);
00411                         break;
00412                         case STOP_AT_TINY_STEP:
00413                                 solutionDescription = "STOP_AT_TINY_STEP[IPOPT]: Algorithm proceeds with very little progress.";
00414                                 osresult->setSolutionStatus(solIdx,  "stoppedByLimit", solutionDescription);
00415                                 osresult->setPrimalVariableValues(solIdx, const_cast<double*>( x));
00416                                 osresult->setDualVariableValues(solIdx, const_cast<double*>( lambda));
00417                                 mdObjValues[0] = obj_value;
00418                                 osresult->setObjectiveValues(solIdx, mdObjValues);
00419                         break;
00420                         case STOP_AT_ACCEPTABLE_POINT:
00421                                 solutionDescription = "STOP_AT_ACCEPTABLE_POINT[IPOPT]: Algorithm stopped at a point that was converged, not to _desired_ tolerances, but to _acceptable_ tolerances";
00422                                 osresult->setSolutionStatus(solIdx,  "IpoptAccetable", solutionDescription);
00423                                 osresult->setPrimalVariableValues(solIdx, const_cast<double*>(x));
00424                                 osresult->setDualVariableValues(solIdx, const_cast<double*>( lambda));
00425                                 mdObjValues[0] = obj_value;
00426                                 osresult->setObjectiveValues(solIdx, mdObjValues);
00427                         break;
00428                         case LOCAL_INFEASIBILITY:
00429                                 solutionDescription = "LOCAL_INFEASIBILITY[IPOPT]: Algorithm converged to a point of local infeasibility. Problem may be infeasible.";
00430                                 osresult->setSolutionStatus(solIdx,  "infeasible", solutionDescription);
00431                         break;
00432                         case USER_REQUESTED_STOP:
00433                                 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.";
00434                                 osresult->setSolutionStatus(solIdx,  "error", solutionDescription);
00435                         break;
00436                         case DIVERGING_ITERATES:
00437                                 solutionDescription = "DIVERGING_ITERATES[IPOPT]: It seems that the iterates diverge.";
00438                                 osresult->setSolutionStatus(solIdx,  "unbounded", solutionDescription);
00439                         break;
00440                         case RESTORATION_FAILURE:
00441                                 solutionDescription = "RESTORATION_FAILURE[IPOPT]: Restoration phase failed, algorithm doesn't know how to proceed.";
00442                                 osresult->setSolutionStatus(solIdx,  "error", solutionDescription);
00443                         break;
00444                         case ERROR_IN_STEP_COMPUTATION:
00445                                 solutionDescription = "ERROR_IN_STEP_COMPUTATION[IPOPT]: An unrecoverable error occurred while IPOPT tried to compute the search direction.";
00446                                 osresult->setSolutionStatus(solIdx,  "error", solutionDescription);
00447                         break;
00448                         case INVALID_NUMBER_DETECTED:
00449                                 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.";
00450                                 osresult->setSolutionStatus(solIdx,  "error", solutionDescription);
00451                         break;
00452                         case INTERNAL_ERROR:
00453                                 solutionDescription = "INTERNAL_ERROR[IPOPT]: An unknown internal error occurred. Please contact the IPOPT authors through the mailing list.";
00454                                 osresult->setSolutionStatus(solIdx,  "error", solutionDescription);
00455                         break;
00456                         default:
00457                                 solutionDescription = "OTHER[IPOPT]: other unknown solution status from Ipopt solver";
00458                                 osresult->setSolutionStatus(solIdx,  "other", solutionDescription);
00459                 }
00460                 osresult->setGeneralStatusType("success");
00461                 delete osrlwriter;
00462                 delete[] mdObjValues;
00463                 osrlwriter = NULL;
00464 
00465         }
00466         catch(const ErrorClass& eclass){
00467                 osresult->setGeneralMessage( eclass.errormsg);
00468                 osresult->setGeneralStatusType( "error");
00469                 std::string osrl = osrlwriter->writeOSrL( osresult);
00470                 delete osrlwriter;
00471                 osrlwriter = NULL;
00472                 throw ErrorClass(  osrl) ;
00473                 delete[] mdObjValues;
00474                 mdObjValues = NULL;
00475         }
00477 }
00478 
00479 
00480 void IpoptSolver::buildSolverInstance() throw (ErrorClass) {
00481         try{
00482                 
00483                 if(osil.length() == 0 && osinstance == NULL) throw ErrorClass("there is no instance");
00484                 if(osinstance == NULL){
00485                         m_osilreader = new OSiLReader();
00486                         osinstance = m_osilreader->readOSiL( osil);
00487                 }
00488                 // Create a new instance of your nlp 
00489                 nlp = new IpoptProblem( osinstance, osresult);
00490                 app = new IpoptApplication();
00491                 this->bCallbuildSolverInstance = true;
00492         }
00493         catch(const ErrorClass& eclass){
00494                 std::cout << "THERE IS AN ERROR" << std::endl;
00495                 osresult->setGeneralMessage( eclass.errormsg);
00496                 osresult->setGeneralStatusType( "error");
00497                 osrl = osrlwriter->writeOSrL( osresult);
00498                 throw ErrorClass( osrl) ;
00499         }                               
00500 }//end buildSolverInstance() 
00501 
00502 
00503 //void IpoptSolver::solve() throw (ErrorClass) {
00504 void IpoptSolver::solve() throw (ErrorClass) {
00505         if( this->bCallbuildSolverInstance == false) buildSolverInstance();
00506         try{
00507                 clock_t start, finish;
00508                 double duration;
00509                 start = clock();
00510                 //OSiLWriter osilwriter;
00511                 //cout << osilwriter.writeOSiL( osinstance) << endl;
00512                 if(osinstance->getVariableNumber() <= 0)throw ErrorClass("Ipopt requires decision variables");
00513                 finish = clock();
00514                 duration = (double) (finish - start) / CLOCKS_PER_SEC;
00515                 cout << "Parsing took (seconds): "<< duration << endl;
00516                 //dataEchoCheck();
00517                 /***************now the ipopt invokation*********************/
00518                 app->Options()->SetNumericValue("tol", 1e-9);
00519                 app->Options()->SetIntegerValue("print_level", 0);
00520                 app->Options()->SetIntegerValue("max_iter", 20000);
00521                 app->Options()->SetStringValue("mu_strategy", "adaptive");
00522                 app->Options()->SetStringValue("output_file", "ipopt.out");
00523                 app->Options()->SetStringValue("check_derivatives_for_naninf", "yes");
00524                 // see if we have a linear program
00525                 if(osinstance->getObjectiveNumber() <= 0) throw ErrorClass("Ipopt NEEDS AN OBJECTIVE FUNCTION");
00526                 if( (osinstance->getNumberOfNonlinearExpressions() == 0) && (osinstance->getNumberOfQuadraticTerms() == 0) ) 
00527                         app->Options()->SetStringValue("hessian_approximation", "limited-memory");
00528                 // if it is a max problem call scaling and set to -1
00529                 if( osinstance->instanceData->objectives->obj[ 0]->maxOrMin.compare("min") != 0){
00530                         app->Options()->SetStringValue("nlp_scaling_method", "user-scaling");
00531                 }
00532                 // Intialize the IpoptApplication and process the options
00533                 std::cout << "Call Ipopt Initialize" << std::endl;
00534                 app->Initialize();
00535                 std::cout << "Finished Ipopt Initialize" << std::endl;
00536                 //nlp->osinstance = this->osinstance;
00537                 // Ask Ipopt to solve the problem
00538                 std::cout << "Call Ipopt Optimize" << std::endl;
00539                 ApplicationReturnStatus status = app->OptimizeTNLP( nlp);
00540                 std::cout << "Finish Ipopt Optimize" << std::endl;
00541                 osrl = osrlwriter->writeOSrL( osresult);
00542                 std::cout << "Finish writing the osrl" << std::endl;
00543                 //if (status != Solve_Succeeded) {
00544                 if (status < -2) {
00545                         throw ErrorClass("Ipopt FAILED TO SOLVE THE PROBLEM: " + ipoptErrorMsg);
00546                 }       
00547         }
00548         catch(const ErrorClass& eclass){
00549                 osresult->setGeneralMessage( eclass.errormsg);
00550                 osresult->setGeneralStatusType( "error");
00551                 osrl = osrlwriter->writeOSrL( osresult);
00552                 throw ErrorClass( osrl) ;
00553         }
00554 }//solve
00555 
00556 
00557 void IpoptSolver::dataEchoCheck(){
00558         
00559         int i;
00560         
00561         // print out problem parameters
00562         cout << "This is problem:  " << osinstance->getInstanceName() << endl;
00563         cout << "The problem source is:  " << osinstance->getInstanceSource() << endl;
00564         cout << "The problem description is:  " << osinstance->getInstanceDescription() << endl;
00565         cout << "number of variables = " << osinstance->getVariableNumber() << endl;
00566         cout << "number of Rows = " << osinstance->getConstraintNumber() << endl;
00567 
00568         // print out the variable information
00569         if(osinstance->getVariableNumber() > 0){
00570                 for(i = 0; i < osinstance->getVariableNumber(); i++){
00571                         if(osinstance->getVariableNames() != NULL) cout << "variable Names  " << osinstance->getVariableNames()[ i]  << endl;
00572                         if(osinstance->getVariableTypes() != NULL) cout << "variable Types  " << osinstance->getVariableTypes()[ i]  << endl;
00573                         if(osinstance->getVariableLowerBounds() != NULL) cout << "variable Lower Bounds  " << osinstance->getVariableLowerBounds()[ i]  << endl;
00574                         if(osinstance->getVariableUpperBounds() != NULL) cout << "variable Upper Bounds  " <<  osinstance->getVariableUpperBounds()[i] << endl;
00575                 }
00576         }
00577         
00578         // print out objective function information
00579         if(osinstance->getVariableNumber() > 0 || osinstance->instanceData->objectives->obj != NULL || osinstance->instanceData->objectives->numberOfObjectives > 0){
00580                 if( osinstance->getObjectiveMaxOrMins()[0] == "min")  cout <<  "problem is a minimization" << endl;
00581                 else cout <<  "problem is a maximization" << endl;
00582                 for(i = 0; i < osinstance->getVariableNumber(); i++){
00583                         cout << "OBJ COEFFICIENT =  " <<  osinstance->getDenseObjectiveCoefficients()[0][i] << endl;
00584                 }
00585         }
00586         // print out constraint information
00587         if(osinstance->getConstraintNumber() > 0){
00588                 for(i = 0; i < osinstance->getConstraintNumber(); i++){
00589                         if(osinstance->getConstraintNames() != NULL) cout << "row name = " << osinstance->getConstraintNames()[i] <<  endl;
00590                         if(osinstance->getConstraintLowerBounds() != NULL) cout << "row lower bound = " << osinstance->getConstraintLowerBounds()[i] <<  endl;
00591                         if(osinstance->getConstraintUpperBounds() != NULL) cout << "row upper bound = " << osinstance->getConstraintUpperBounds()[i] <<  endl; 
00592                 }
00593         }
00594         
00595         // print out linear constraint data
00596         cout << endl;
00597         cout << "number of nonzeros =  " << osinstance->getLinearConstraintCoefficientNumber() << endl;
00598         for(i = 0; i <= osinstance->getVariableNumber(); i++){
00599                 cout << "Start Value =  " << osinstance->getLinearConstraintCoefficientsInColumnMajor()->starts[ i] << endl;
00600         }
00601         cout << endl;
00602         for(i = 0; i < osinstance->getLinearConstraintCoefficientNumber(); i++){
00603                 cout << "Index Value =  " << osinstance->getLinearConstraintCoefficientsInColumnMajor()->indexes[i] << endl;
00604                 cout << "Nonzero Value =  " << osinstance->getLinearConstraintCoefficientsInColumnMajor()->values[i] << endl;
00605         }
00606 
00607         // print out quadratic data
00608         cout << "number of qterms =  " <<  osinstance->getNumberOfQuadraticTerms() << endl;
00609         for(int i = 0; i <  osinstance->getNumberOfQuadraticTerms(); i++){
00610                 cout << "Row Index = " <<  osinstance->getQuadraticTerms()->rowIndexes[i] << endl;
00611                 cout << "Var Index 1 = " << osinstance->getQuadraticTerms()->varOneIndexes[ i] << endl;
00612                 cout << "Var Index 2 = " << osinstance->getQuadraticTerms()->varTwoIndexes[ i] << endl;
00613                 cout << "Coefficient = " << osinstance->getQuadraticTerms()->coefficients[ i] << endl;
00614         }
00615 } // end dataEchoCheck
00616 
00617 
00618 IpoptProblem::IpoptProblem(OSInstance *osinstance_,  OSResult *osresult_) {
00619         osinstance = osinstance_;
00620         osresult = osresult_;
00621 }
00622 
00623 IpoptProblem::~IpoptProblem() {
00624 
00625 }
00626 
00627 
00628 

Generated on Sat Mar 29 22:38:02 2008 by  doxygen 1.5.3