33 #ifndef _InterPhaseBase_HPP_
34 #define _InterPhaseBase_HPP_
39 namespace DirectTrajectoryOptimization{
42 template <
class DIMENSIONS>
43 class InterPhaseBase :
public NumDiffDerivativesBase,
44 public std::enable_shared_from_this<InterPhaseBase<DIMENSIONS> > {
48 EIGEN_MAKE_ALIGNED_OPERATOR_NEW
50 typedef typename DIMENSIONS::state_vector_t state_vector_t;
51 typedef typename DIMENSIONS::control_vector_t control_vector_t;
53 static const int s_size =
static_cast<int>(DIMENSIONS::DimensionsSize::STATE_SIZE);
54 static const int c_size =
static_cast<int>(DIMENSIONS::DimensionsSize::CONTROL_SIZE);
56 int complete_vector_size = 0;
57 int size_of_jacobian = 0;
58 int total_inputs_jacobian=0;
60 std::shared_ptr<InterPhaseBase<DIMENSIONS> > my_pointer;
63 InterPhaseBase(
const bool & no_ipc) {
66 complete_vector_size = s_size;
67 ipc_lb = Eigen::VectorXd::Zero(complete_vector_size);
68 ipc_ub = Eigen::VectorXd::Zero(complete_vector_size);
73 virtual ~InterPhaseBase(){ };
74 int GetSize() {
return complete_vector_size;}
75 int GetNnzJacobian() {
return size_of_jacobian;}
76 Eigen::VectorXd getLowerBound(){
return ipc_lb;}
77 Eigen::VectorXd getUpperBound(){
return ipc_ub;}
79 void setSize(
const int & soc);
80 void setBounds(
const Eigen::VectorXd & lb ,
const Eigen::VectorXd & ub);
83 void initialize_num_diff();
84 void fx(
const Eigen::VectorXd & in, Eigen::VectorXd & out);
87 virtual Eigen::VectorXd evalConstraintsFct(
const state_vector_t& x1,
const state_vector_t & x2,
88 const control_vector_t & u1,
const control_vector_t & u2);
90 Eigen::VectorXd evalConstraintsFctDerivatives(
const state_vector_t & x1,
const state_vector_t & x2,
91 const control_vector_t & u1,
const control_vector_t & u2);
94 Eigen::VectorXd ipc_lb;
95 Eigen::VectorXd ipc_ub;
97 bool soc_ready_ =
false;
98 bool bounds_ready_ =
false;
102 template <
class DIMENSIONS>
103 void InterPhaseBase<DIMENSIONS>::setSize(
const int & soc) {
105 std::cout <<
"Interphase size must be greater than 0" << std::endl;
106 std::exit(EXIT_FAILURE);
108 complete_vector_size = soc;
112 template <
class DIMENSIONS>
113 void InterPhaseBase<DIMENSIONS>::setBounds(
const Eigen::VectorXd & lb ,
const Eigen::VectorXd & ub) {
114 if(!(lb.size() == complete_vector_size) || !(ub.size() == complete_vector_size)) {
115 std::cout <<
"Interphase bounds size must be equal to SoG" << std::endl;
116 std::exit(EXIT_FAILURE);
120 bounds_ready_ =
true;
123 template <
class DIMENSIONS>
124 void InterPhaseBase<DIMENSIONS>::initialize() {
126 if(!soc_ready_ || !bounds_ready_) {
127 std::cout <<
"Interphase is not ready!" << std::endl;
128 std::exit(EXIT_FAILURE);
132 total_inputs_jacobian = 2*s_size + 2*c_size;
133 size_of_jacobian = complete_vector_size*total_inputs_jacobian;
134 initialize_num_diff();
137 template <
class DIMENSIONS>
138 Eigen::VectorXd InterPhaseBase<DIMENSIONS>::evalConstraintsFct(
const state_vector_t& x1,
const state_vector_t & x2,
139 const control_vector_t & u1,
const control_vector_t & u2) {
144 template <
class DIMENSIONS>
145 Eigen::VectorXd InterPhaseBase<DIMENSIONS>::evalConstraintsFctDerivatives(
146 const state_vector_t & x1,
const state_vector_t & x2 ,
147 const control_vector_t &u1,
const control_vector_t &u2) {
149 Eigen::VectorXd in_vector(total_inputs_jacobian);
151 in_vector.segment<s_size>(0) = x1;
152 in_vector.segment<c_size>(s_size) = u1;
153 in_vector.segment<s_size>(s_size+c_size)= x2;
154 in_vector.segment<c_size>(s_size+s_size+c_size) = u2;
156 this->numDiff->df(in_vector,this->mJacobian);
159 Eigen::Matrix<double,Eigen::Dynamic,Eigen::Dynamic,Eigen::RowMajor> temp(this->mJacobian);
160 return (Eigen::Map<Eigen::VectorXd>(temp.data(),size_of_jacobian));
163 template <
class DIMENSIONS>
164 void InterPhaseBase<DIMENSIONS>::fx(
const Eigen::VectorXd & inputs, Eigen::VectorXd & outputs) {
166 state_vector_t x1 = inputs.segment<s_size>(0);
167 control_vector_t u1 = inputs.segment<c_size>(s_size);
168 state_vector_t x2 = inputs.segment<s_size>(s_size+c_size);
169 control_vector_t u2 = inputs.segment<c_size>(s_size+s_size+c_size);
171 outputs = evalConstraintsFct(x1,x2,u1,u2);
174 template <
class DIMENSIONS>
175 void InterPhaseBase<DIMENSIONS>::initialize_num_diff() {
177 this->my_pointer = std::enable_shared_from_this<InterPhaseBase<DIMENSIONS>>::shared_from_this();
178 this->mJacobian.resize(complete_vector_size,total_inputs_jacobian);
179 this->numdifoperator = std::shared_ptr<FunctionOperator>(
new FunctionOperator(total_inputs_jacobian, complete_vector_size , this->my_pointer));
180 this->numDiff = std::shared_ptr<Eigen::NumericalDiff<FunctionOperator> >(
new Eigen::NumericalDiff<FunctionOperator>(*this->numdifoperator,std::numeric_limits<double>::epsilon()));
Base class for vector constraints and numdiff derivatives.