31 template <
size_t STATE_DIM,
size_t CONTROL_DIM,
typename SCALAR =
double>
35 EIGEN_MAKE_ALIGNED_OPERATOR_NEW
61 : costFct_(costFct), w_(w), controlSpliner_(controlSpliner), shotContainers_(shotInt), settings_(settings)
74 #pragma omp parallel for num_threads(settings_.nThreads_) 75 for (
auto shotContainer = shotContainers_.begin(); shotContainer < shotContainers_.end(); ++shotContainer)
77 (*shotContainer)->integrateCost();
80 for (
auto shotContainer : shotContainers_)
81 cost += shotContainer->getCostIntegrated();
83 costFct_->setCurrentStateAndControl(w_->getOptimizedState(settings_.
N_), control_vector_t::Zero());
84 cost += costFct_->evaluateTerminal();
88 void evalGradient(
size_t grad_length, Eigen::Map<Eigen::Matrix<SCALAR, Eigen::Dynamic, 1>>& grad)
override 92 assert(shotContainers_.size() == settings_.
N_);
96 #pragma omp parallel for num_threads(settings_.nThreads_) 97 for (
auto shotContainer = shotContainers_.begin(); shotContainer < shotContainers_.end(); ++shotContainer)
99 (*shotContainer)->integrateCostSensitivities();
102 for (
size_t shotNr = 0; shotNr < shotContainers_.size(); ++shotNr)
108 grad.segment(w_->getStateIndex(shotNr), STATE_DIM) += shotContainers_[shotNr]->getdLdSiIntegrated();
109 grad.segment(w_->getControlIndex(shotNr), CONTROL_DIM) +=
110 shotContainers_[shotNr]->getdLdQiIntegrated();
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();
123 throw(std::runtime_error(
124 " cost gradient not yet implemented for this type of interpolation. Exiting"));
137 costFct_->setCurrentStateAndControl(w_->getOptimizedState(settings_.
N_), control_vector_t::Zero());
138 grad.segment(w_->getStateIndex(settings_.
N_), STATE_DIM) +=
139 costFct_->stateDerivativeTerminal();
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_;
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
CostEvaluatorFull()=delete
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