DirectTrajectoryOptimization  v0.3
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Pages
base_solver_interface.hpp
Go to the documentation of this file.
1 /***********************************************************************************
2 Copyright (c) 2017, Diego Pardo. All rights reserved.
3 
4 Redistribution and use in source and binary forms, with or without modification,
5 are permitted provided that the following conditions are met:
6  * Redistributions of source code must retain the above copyright notice,
7  this list of conditions and the following disclaimer.
8  * Redistributions in binary form must reproduce the above copyright notice,
9  this list of conditions and the following disclaimer in the documentation
10  and/or other materials provided with the distribution.
11  * Neither the name of ETH ZURICH nor the names of its contributors may be used
12  to endorse or promote products derived from this software without specific
13  prior written permission.
14 
15 THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" AND ANY
16 EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED WARRANTIES
17 OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT
18 SHALL ETH ZURICH BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY,
19 OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE
20 GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
21 CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
22 (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE,
23 EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
24 ***************************************************************************************/
25 
31 #include <Eigen/Dense>
32 #include <unsupported/Eigen/NumericalDiff>
34 
35 #ifndef INCLUDE_SOLVER_INTERFACES_BASE_SOLVER_INTERFACE_HPP_
36 #define INCLUDE_SOLVER_INTERFACES_BASE_SOLVER_INTERFACE_HPP_
37 
38 namespace DirectTrajectoryOptimization {
39 
46 template <typename DIMENSIONS>
48 
49 public:
50 
51  EIGEN_MAKE_ALIGNED_OPERATOR_NEW
52 
53  typedef typename DIMENSIONS::state_vector_t state_vector_t;
54  typedef typename DIMENSIONS::control_vector_t control_vector_t;
55  typedef typename DIMENSIONS::state_vector_array_t state_vector_array_t;
56  typedef typename DIMENSIONS::control_vector_array_t control_vector_array_t;
57 
58  virtual ~DTOSolverInterface() {};
59 
67  virtual void InitializeSolver(
68  const std::string & problemName,
69  const std::string & outputFile,
70  const std::string & paramFile /*= "" */,
71  const bool &computeJacobian) = 0;
72 
77  virtual void Solve(bool warminit = false) = 0;
78 
86  virtual void GetDTOSolution(
87  state_vector_array_t &yTraj,
88  control_vector_array_t &uTraj,
89  Eigen::VectorXd &hTraj) = 0;
90 
99  virtual void GetDTOSolution(
100  state_vector_array_t & yTraj,
101  control_vector_array_t & uTraj,
102  Eigen::VectorXd & hTraj,
103  std::vector<int> & phaseId) = 0;
104 
114  virtual void GetDTOSolution(
115  state_vector_array_t &yTraj,
116  state_vector_array_t &ydotTraj,
117  control_vector_array_t &uTraj,
118  Eigen::VectorXd &hTraj,
119  std::vector<int> &phaseId) { std::cout << "not here" << std::endl; getchar();};
120 
121 
131  Eigen::VectorXi & _x_state,
132  Eigen::VectorXd & _x_mul,
133  Eigen::VectorXi & _f_state,
134  Eigen::VectorXd & _f_mul) {}
135 
136 
144  const state_vector_array_t & y_sol,
145  const control_vector_array_t & u_sol,
146  const Eigen::VectorXd & h_sol) = 0;
147 
148 
157  const Eigen::VectorXi & x_state,
158  const Eigen::VectorXd & x_mul,
159  const Eigen::VectorXi & f_state,
160  const Eigen::VectorXd & f_mul) {};
161 
162  virtual void SetVerifyLevel(const int & level) { }
163  virtual void UsePrintFile(const bool & flag) { }
164 
165 
169  virtual void SetupSolverParameters() = 0 ;
170 
175  virtual void GetFVector(Eigen::VectorXd &fVec) = 0;
176 
181  void SetVerbosity(bool verb) { _verbose = verb; }
182 
183  void TotallyQuiet(bool quiet) { _totally_quiet_ = quiet; }
184 
185  void SaveBasisFile(const bool & flag) {
186  _save_this_solution_ = true;
187  new_basis_file_param = 30;
188  }
189 
190  void LoadBasisFile(const bool & flag , int name_number ) {
191  _load_from_solution_ = true;
192  old_basis_file_param = name_number;
193  }
194 
195 
196  //member variables
197  bool _totally_quiet_ = false;
198  bool _verbose = true;
199  bool _save_this_solution_ = false;
200  bool _load_from_solution_ = false;
201  int new_basis_file_param = 0;
202  int old_basis_file_param = 0;
203 
204 protected:
205 
212 
213 };
214 } // namespace DirectTrajectoryOptimization
215 #endif /* INCLUDE_SOLVER_INTERFACES_BASE_SOLVER_INTERFACE_HPP_ */
virtual void GetDTOSolution(state_vector_array_t &yTraj, control_vector_array_t &uTraj, Eigen::VectorXd &hTraj)=0
Method to get the current values of the solution.
virtual void SetupSolverParameters()=0
Setting up solver parameters.
void SetVerbosity(bool verb)
Set verbosity of the solver.
Definition: base_solver_interface.hpp:181
virtual void GetSolverSolutionVariables(Eigen::VectorXi &_x_state, Eigen::VectorXd &_x_mul, Eigen::VectorXi &_f_state, Eigen::VectorXd &_f_mul)
Method to get final values of internal solver variables.
Definition: base_solver_interface.hpp:130
virtual void InitializeSolver(const std::string &problemName, const std::string &outputFile, const std::string &paramFile, const bool &computeJacobian)=0
DTO initializes the solver using this method.
This class serves as a base interface for inheriting classes.
Definition: base_solver_interface.hpp:47
virtual void Solve(bool warminit=false)=0
DTO request the solver to proceed with the NLP.
virtual void WarmSolverInitialization(const Eigen::VectorXi &x_state, const Eigen::VectorXd &x_mul, const Eigen::VectorXi &f_state, const Eigen::VectorXd &f_mul)
Warm Initialization of NLP Solver.
Definition: base_solver_interface.hpp:156
DTOSolverInterface()
Construct an instance of the DTOSolverInterface class.
Definition: base_solver_interface.hpp:211
virtual void InitializeDecisionVariablesFromTrajectory(const state_vector_array_t &y_sol, const control_vector_array_t &u_sol, const Eigen::VectorXd &h_sol)=0
Initial solution point for the solver.
virtual void GetFVector(Eigen::VectorXd &fVec)=0
Returns the vector of NLP-constraints.
virtual void GetDTOSolution(state_vector_array_t &yTraj, state_vector_array_t &ydotTraj, control_vector_array_t &uTraj, Eigen::VectorXd &hTraj, std::vector< int > &phaseId)
Method to get the current values of the solution.
Definition: base_solver_interface.hpp:114
Base class for direct methods.