- 3.0.2 optimal control module.
CostEvaluatorFull.h
Go to the documentation of this file.
1 /**********************************************************************************************************************
2 This file is part of the Control Toolbox (https://github.com/ethz-adrl/control-toolbox), copyright by ETH Zurich.
3 Licensed under the BSD-2 license (see LICENSE file in main directory)
4 **********************************************************************************************************************/
5 
6 #pragma once
7 
8 #include <omp.h>
9 #include <math.h>
10 #include <cmath>
11 #include <functional>
12 
14 
18 
19 
20 namespace ct {
21 namespace optcon {
22 
31 template <size_t STATE_DIM, size_t CONTROL_DIM, typename SCALAR = double>
33 {
34 public:
35  EIGEN_MAKE_ALIGNED_OPERATOR_NEW
36 
38 
39  typedef typename DIMENSIONS::state_vector_t state_vector_t;
40  typedef typename DIMENSIONS::control_vector_t control_vector_t;
41  typedef typename DIMENSIONS::state_vector_array_t state_vector_array_t;
42  typedef typename DIMENSIONS::control_vector_array_t control_vector_array_t;
43  typedef typename DIMENSIONS::time_array_t time_array_t;
44 
45  CostEvaluatorFull() = delete;
46 
58  std::shared_ptr<SplinerBase<control_vector_t, SCALAR>> controlSpliner,
59  std::vector<std::shared_ptr<ShotContainer<STATE_DIM, CONTROL_DIM, SCALAR>>> shotInt,
60  DmsSettings settings)
61  : costFct_(costFct), w_(w), controlSpliner_(controlSpliner), shotContainers_(shotInt), settings_(settings)
62  {
63  }
64 
68  ~CostEvaluatorFull() override = default;
69 
70  SCALAR eval() override
71  {
72  SCALAR cost = SCALAR(0.0);
73 
74 #pragma omp parallel for num_threads(settings_.nThreads_)
75  for (auto shotContainer = shotContainers_.begin(); shotContainer < shotContainers_.end(); ++shotContainer)
76  {
77  (*shotContainer)->integrateCost();
78  }
79 
80  for (auto shotContainer : shotContainers_)
81  cost += shotContainer->getCostIntegrated();
82 
83  costFct_->setCurrentStateAndControl(w_->getOptimizedState(settings_.N_), control_vector_t::Zero());
84  cost += costFct_->evaluateTerminal();
85  return cost;
86  }
87 
88  void evalGradient(size_t grad_length, Eigen::Map<Eigen::Matrix<SCALAR, Eigen::Dynamic, 1>>& grad) override
89  {
90  grad.setZero();
91 
92  assert(shotContainers_.size() == settings_.N_);
93 
94 // go through all shots, integrate the state trajectories and evaluate cost accordingly
95 // intermediate costs
96 #pragma omp parallel for num_threads(settings_.nThreads_)
97  for (auto shotContainer = shotContainers_.begin(); shotContainer < shotContainers_.end(); ++shotContainer)
98  {
99  (*shotContainer)->integrateCostSensitivities();
100  }
101 
102  for (size_t shotNr = 0; shotNr < shotContainers_.size(); ++shotNr)
103  {
104  switch (settings_.splineType_)
105  {
107  {
108  grad.segment(w_->getStateIndex(shotNr), STATE_DIM) += shotContainers_[shotNr]->getdLdSiIntegrated();
109  grad.segment(w_->getControlIndex(shotNr), CONTROL_DIM) +=
110  shotContainers_[shotNr]->getdLdQiIntegrated();
111  break;
112  }
114  {
115  grad.segment(w_->getStateIndex(shotNr), STATE_DIM) += shotContainers_[shotNr]->getdLdSiIntegrated();
116  grad.segment(w_->getControlIndex(shotNr), CONTROL_DIM) +=
117  shotContainers_[shotNr]->getdLdQiIntegrated();
118  grad.segment(w_->getControlIndex(shotNr + 1), CONTROL_DIM) +=
119  shotContainers_[shotNr]->getdLdQip1Integrated();
120  break;
121  }
122  default:
123  throw(std::runtime_error(
124  " cost gradient not yet implemented for this type of interpolation. Exiting"));
125  }
126 
127  // H-part.
128  // if(settings_.objectiveType_ == DmsSettings::OPTIMIZE_GRID)
129  // {
130  // costFct_->setCurrentStateAndControl(shotContainers_[shotNr]->getStateIntegrated(),
131  // controlSpliner_->evalSpline(shotContainers_[shotNr]->getIntegrationTimeFinal(), shotNr));
132  // grad(w_->getTimeSegmentIndex(shotNr)) = costFct_->evaluateIntermediate() + shotContainers_[shotNr]->getdLdHiIntegrated();
133  // }
134  }
135 
136  /* gradient of terminal cost */
137  costFct_->setCurrentStateAndControl(w_->getOptimizedState(settings_.N_), control_vector_t::Zero());
138  grad.segment(w_->getStateIndex(settings_.N_), STATE_DIM) +=
139  costFct_->stateDerivativeTerminal(); // * dXdSi.back();
140  }
141 
142 private:
143  std::shared_ptr<ct::optcon::CostFunctionQuadratic<STATE_DIM, CONTROL_DIM, SCALAR>> costFct_;
144  std::shared_ptr<OptVectorDms<STATE_DIM, CONTROL_DIM, SCALAR>> w_;
145  std::shared_ptr<SplinerBase<control_vector_t, SCALAR>> controlSpliner_;
146  std::vector<std::shared_ptr<ShotContainer<STATE_DIM, CONTROL_DIM, SCALAR>>> shotContainers_;
147 
148  const DmsSettings settings_;
149 };
150 
151 } // namespace optcon
152 } // namespace ct
Implements an abstract base class which evaluates the cost function and its gradient in the NLP...
Definition: DiscreteCostEvaluatorBase.h:22
DIMENSIONS::state_vector_array_t state_vector_array_t
Definition: CostEvaluatorFull.h:41
Definition: DmsSettings.h:26
Performs the full cost integration over the shots.
Definition: CostEvaluatorFull.h:32
Abstract base class for the control input splining between the DMS shots.
Definition: SplinerBase.h:20
CostEvaluatorFull(std::shared_ptr< ct::optcon::CostFunctionQuadratic< STATE_DIM, CONTROL_DIM, SCALAR >> costFct, std::shared_ptr< OptVectorDms< STATE_DIM, CONTROL_DIM, SCALAR >> w, std::shared_ptr< SplinerBase< control_vector_t, SCALAR >> controlSpliner, std::vector< std::shared_ptr< ShotContainer< STATE_DIM, CONTROL_DIM, SCALAR >>> shotInt, DmsSettings settings)
Custom constructor.
Definition: CostEvaluatorFull.h:56
This class is a wrapper around the NLP Optvector. It wraps the Vectors from the NLP solvers into stat...
Definition: OptVectorDms.h:37
DIMENSIONS::time_array_t time_array_t
Definition: CostEvaluatorFull.h:43
This class performs the state and the sensitivity integration on a shot.
Definition: ShotContainer.h:32
DIMENSIONS::state_vector_t state_vector_t
Definition: CostEvaluatorFull.h:39
size_t N_
Definition: DmsSettings.h:51
Describes a cost function with a quadratic approximation, i.e. one that can compute first and second ...
Definition: CostFunctionQuadratic.hpp:29
Definition: DmsSettings.h:26
CppAD::AD< CppAD::cg::CG< double > > SCALAR
EIGEN_MAKE_ALIGNED_OPERATOR_NEW typedef DmsDimensions< STATE_DIM, CONTROL_DIM, SCALAR > DIMENSIONS
Definition: CostEvaluatorFull.h:37
SCALAR eval() override
Evaluates the cost function.
Definition: CostEvaluatorFull.h:70
Defines basic types used in the DMS algorithm.
Definition: DmsDimensions.h:18
void evalGradient(size_t grad_length, Eigen::Map< Eigen::Matrix< SCALAR, Eigen::Dynamic, 1 >> &grad) override
Evaluates the cost gradient.
Definition: CostEvaluatorFull.h:88
SplineType_t splineType_
Definition: DmsSettings.h:54
Defines the DMS settings.
Definition: DmsSettings.h:23
DIMENSIONS::control_vector_t control_vector_t
Definition: CostEvaluatorFull.h:40
~CostEvaluatorFull() override=default
The destructor.
DIMENSIONS::control_vector_array_t control_vector_array_t
Definition: CostEvaluatorFull.h:42