My Project
OSIpoptSolver.h
Go to the documentation of this file.
1/* $Id$ */
15#ifndef IPOPTSOLVER_H
16#define IPOPTSOLVER_H
17
18#include "OSConfig.h"
19#include "OSDefaultSolver.h"
20#include "OSrLWriter.h"
21#include "OSInstance.h"
22#include "OSParameters.h"
23#include "OSnLNode.h"
24#include "OSiLReader.h"
25#include "OSoLReader.h"
26#include "OSInstance.h"
27#include "OSExpressionTree.h"
28#include "OSnLNode.h"
29#include "OSGeneral.h"
30#include "OSFileUtil.h"
31#include "OSErrorClass.h"
32
33#include "OSResult.h"
34#include "OSInstance.h"
35#include "OSOption.h"
36
37#include <cstddef>
38#include <cstdlib>
39#include <cctype>
40#include <cassert>
41#include <stack>
42#include <string>
43#include <iostream>
44#include <vector>
45#include <map>
46
47
48#include "IpTNLP.hpp"
49#include "IpIpoptApplication.hpp"
50
51// for Stefan
52class IpoptProblem : public Ipopt::TNLP
53{
54public:
55
57 IpoptProblem(OSInstance *osinstance_ , OSOption *osoption_, OSResult *osresult_, std::string *ipoptErrorMsg_);
58
60 virtual ~IpoptProblem();
61
62
64
66
68
69 std::string *ipoptErrorMsg;
70
72 virtual bool get_nlp_info(Ipopt::Index& n, Ipopt::Index& m, Ipopt::Index& nnz_jac_g,
73 Ipopt::Index& nnz_h_lag, IndexStyleEnum& index_style);
74
76 virtual bool get_bounds_info(Ipopt::Index n, Ipopt::Number* x_l, Ipopt::Number* x_u,
77 Ipopt::Index m, Ipopt::Number* g_l, Ipopt::Number* g_u);
78
80 virtual bool get_starting_point(Ipopt::Index n, bool init_x, Ipopt::Number* x,
81 bool init_z, Ipopt::Number* z_L, Ipopt::Number* z_U,
82 Ipopt::Index m, bool init_lambda,
83 Ipopt::Number* lambda);
84
86 virtual bool eval_f(Ipopt::Index n, const Ipopt::Number* x, bool new_x, Ipopt::Number& obj_value);
87
89 virtual bool eval_grad_f(Ipopt::Index n, const Ipopt::Number* x, bool new_x, Ipopt::Number* grad_f);
90
92 virtual bool eval_g(Ipopt::Index n, const Ipopt::Number* x, bool new_x, Ipopt::Index m, Ipopt::Number* g);
93
98 virtual bool eval_jac_g(Ipopt::Index n, const Ipopt::Number* x, bool new_x,
99 Ipopt::Index m, Ipopt::Index nele_jac, Ipopt::Index* iRow, Ipopt::Index *jCol,
100 Ipopt::Number* values);
101
106 virtual bool eval_h(Ipopt::Index n, const Ipopt::Number* x, bool new_x,
107 Ipopt::Number obj_factor, Ipopt::Index m, const Ipopt::Number* lambda,
108 bool new_lambda, Ipopt::Index nele_hess, Ipopt::Index* iRow,
109 Ipopt::Index* jCol, Ipopt::Number* values);
110
112
113
114 virtual bool get_scaling_parameters(Ipopt::Number& obj_scaling,
115 bool& use_x_scaling, Ipopt::Index n,
116 Ipopt::Number* x_scaling,
117 bool& use_g_scaling, Ipopt::Index m,
118 Ipopt::Number* g_scaling);
119
123 virtual void finalize_solution(Ipopt::SolverReturn status,
124 Ipopt::Index n, const Ipopt::Number* x, const Ipopt::Number* z_L, const Ipopt::Number* z_U,
125 Ipopt::Index m, const Ipopt::Number* g, const Ipopt::Number* lambda,
126 Ipopt::Number obj_value,
127 const Ipopt::IpoptData* ip_data,
128 Ipopt::IpoptCalculatedQuantities* ip_cq);
130
131
132
133private:
145 // HS071_NLP();
149
150};
151
152
165//class IpoptSolver : public DefaultSolver, public TNLP{
166
168{
169public:
170
172 IpoptSolver();
173
175 ~IpoptSolver();
176
177 Ipopt::SmartPtr<Ipopt::TNLP> nlp;
178
179 Ipopt::SmartPtr<Ipopt::IpoptApplication> app;
180
181
184 virtual void solve();
185
186
191 virtual void buildSolverInstance();
192
197 virtual void setSolverOptions();
198
204 void dataEchoCheck();
205
211
217
218
219private:
221
233 // IpoptSolver();
234 //IpoptSolver(const IpoptSolver&);
235 //IpoptSolver& operator=(const IpoptSolver&);
237 //std::string ipoptErrorMsg;
238 std::string *ipoptErrorMsg;
239};
240
241
242#endif /*IPOPTSOLVER_H*/
This file defines the OSnLNode class along with its derived classes.
The Default Solver Class.
virtual bool get_starting_point(Ipopt::Index n, bool init_x, Ipopt::Number *x, bool init_z, Ipopt::Number *z_L, Ipopt::Number *z_U, Ipopt::Index m, bool init_lambda, Ipopt::Number *lambda)
Method to return the starting point for the algorithm.
OSOption * osoption
OSInstance * osinstance
virtual bool get_scaling_parameters(Ipopt::Number &obj_scaling, bool &use_x_scaling, Ipopt::Index n, Ipopt::Number *x_scaling, bool &use_g_scaling, Ipopt::Index m, Ipopt::Number *g_scaling)
IpoptProblem(const IpoptProblem &)
std::string * ipoptErrorMsg
virtual bool eval_grad_f(Ipopt::Index n, const Ipopt::Number *x, bool new_x, Ipopt::Number *grad_f)
Method to return the gradient of the objective.
IpoptProblem(OSInstance *osinstance_, OSOption *osoption_, OSResult *osresult_, std::string *ipoptErrorMsg_)
the IpoptProblemclass constructor
virtual bool eval_g(Ipopt::Index n, const Ipopt::Number *x, bool new_x, Ipopt::Index m, Ipopt::Number *g)
Method to return the constraint residuals.
virtual ~IpoptProblem()
the IpoptProblem class destructor
IpoptProblem & operator=(const IpoptProblem &)
virtual bool get_nlp_info(Ipopt::Index &n, Ipopt::Index &m, Ipopt::Index &nnz_jac_g, Ipopt::Index &nnz_h_lag, IndexStyleEnum &index_style)
IPOpt specific methods for defining the nlp problem.
virtual bool eval_jac_g(Ipopt::Index n, const Ipopt::Number *x, bool new_x, Ipopt::Index m, Ipopt::Index nele_jac, Ipopt::Index *iRow, Ipopt::Index *jCol, Ipopt::Number *values)
Method to return: 1) The structure of the jacobian (if "values" is NULL) 2) The values of the jacobia...
OSResult * osresult
virtual bool eval_f(Ipopt::Index n, const Ipopt::Number *x, bool new_x, Ipopt::Number &obj_value)
Method to return the objective value.
virtual bool get_bounds_info(Ipopt::Index n, Ipopt::Number *x_l, Ipopt::Number *x_u, Ipopt::Index m, Ipopt::Number *g_l, Ipopt::Number *g_u)
Method to return the bounds for my problem.
virtual bool eval_h(Ipopt::Index n, const Ipopt::Number *x, bool new_x, Ipopt::Number obj_factor, Ipopt::Index m, const Ipopt::Number *lambda, bool new_lambda, Ipopt::Index nele_hess, Ipopt::Index *iRow, Ipopt::Index *jCol, Ipopt::Number *values)
Method to return: 1) The structure of the hessian of the lagrangian (if "values" is NULL) 2) The valu...
virtual void finalize_solution(Ipopt::SolverReturn status, Ipopt::Index n, const Ipopt::Number *x, const Ipopt::Number *z_L, const Ipopt::Number *z_U, Ipopt::Index m, const Ipopt::Number *g, const Ipopt::Number *lambda, Ipopt::Number obj_value, const Ipopt::IpoptData *ip_data, Ipopt::IpoptCalculatedQuantities *ip_cq)
This method is called when the algorithm is complete so the TNLP can store/write the solution.
The IpoptSolver class solves problems using Ipopt.
std::string * ipoptErrorMsg
OSoLReader * m_osolreader
m_osolreader is an OSoLReader object used to create an osoption from an osol string if needed
OSrLWriter * osrlwriter
OSiLReader * m_osilreader
m_osilreader is an OSiLReader object used to create an osinstance from an osil string if needed
virtual void setSolverOptions()
The implementation of the virtual functions.
Ipopt::SmartPtr< Ipopt::TNLP > nlp
Ipopt::SmartPtr< Ipopt::IpoptApplication > app
virtual void solve()
solve results in an instance being read into the Ipopt data structures and optimize
IpoptSolver()
the IpoptSolver class constructor
~IpoptSolver()
the IpoptSolver class destructor
void dataEchoCheck()
use this for debugging, print out the instance that the solver thinks it has and compare this with th...
virtual void buildSolverInstance()
The implementation of the virtual functions.
The in-memory representation of an OSiL instance..
The Option Class.
Definition OSOption.h:3565
The Result Class.
Definition OSResult.h:2549
Used to read an OSiL string.
Definition OSiLReader.h:38
Used to read an OSoL string.
Definition OSoLReader.h:38
Take an OSResult object and write a string that validates against OSrL.
Definition OSrLWriter.h:31