optimization_problem_interface.h
Go to the documentation of this file.
1 /*********************************************************************
2  *
3  * Software License Agreement
4  *
5  * Copyright (c) 2020,
6  * TU Dortmund - Institute of Control Theory and Systems Engineering.
7  * All rights reserved.
8  *
9  * This program is free software: you can redistribute it and/or modify
10  * it under the terms of the GNU General Public License as published by
11  * the Free Software Foundation, either version 3 of the License, or
12  * (at your option) any later version.
13  *
14  * This program is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17  * GNU General Public License for more details.
18  *
19  * You should have received a copy of the GNU General Public License
20  * along with this program. If not, see <https://www.gnu.org/licenses/>.
21  *
22  * Authors: Christoph Rösmann
23  *********************************************************************/
24 
25 #ifndef SRC_OPTIMIZATION_INCLUDE_CORBO_OPTIMIZATION_OPTIMIZATION_PROBLEM_INTERFACE_H_
26 #define SRC_OPTIMIZATION_INCLUDE_CORBO_OPTIMIZATION_OPTIMIZATION_PROBLEM_INTERFACE_H_
27 
28 #include <corbo-core/macros.h>
29 #include <corbo-core/types.h>
30 
31 #include <Eigen/Sparse>
32 
33 #include <functional>
34 #include <memory>
35 
36 namespace corbo {
37 
71 {
72  public:
73  using Ptr = std::shared_ptr<OptimizationProblemInterface>;
74  using UPtr = std::unique_ptr<OptimizationProblemInterface>;
75 
77 
78  virtual void clear() {}
79 
82 
84  virtual int getNonLsqObjectiveDimension() = 0;
86  virtual int getLsqObjectiveDimension() = 0;
88  virtual int getObjectiveDimension() = 0;
89 
91  virtual int getEqualityDimension() = 0;
93  virtual int getInequalityDimension() = 0;
95  virtual int getParameterDimension() = 0;
96 
98 
101 
111  virtual void applyIncrement(const Eigen::Ref<const Eigen::VectorXd>& increment);
112 
119  virtual void applyIncrement(int idx, double increment);
120 
126  virtual double getParameterValue(int idx) = 0;
127 
133  virtual void setParameterValue(int idx, double x) = 0;
134 
140 
147 
149  virtual void backupParameters() = 0;
151  virtual void restoreBackupParameters(bool keep_backup) = 0;
153  virtual void discardBackupParameters(bool all = false) = 0;
154 
156 
159 
173  virtual bool isLeastSquaresProblem() const = 0;
174 
176 
179 
185  virtual double getLowerBound(int idx) = 0;
186 
192  virtual double getUpperBound(int idx) = 0;
193 
199  virtual void setLowerBound(int idx, double lb) = 0;
200 
206  virtual void setUpperBound(int idx, double ub) = 0;
207 
214 
221 
223 
226 
233 
234  virtual double computeValueNonLsqObjective() = 0;
235 
236  virtual double computeValueObjective();
237 
243  virtual void computeValuesEquality(Eigen::Ref<Eigen::VectorXd> values) = 0;
244 
250  virtual void computeValuesInequality(Eigen::Ref<Eigen::VectorXd> values) = 0;
251 
252  virtual void computeValues(double& non_lsq_obj_value, Eigen::Ref<Eigen::VectorXd> lsq_obj_values, Eigen::Ref<Eigen::VectorXd> eq_values,
253  Eigen::Ref<Eigen::VectorXd> ineq_values)
254  {
255  non_lsq_obj_value = computeValueNonLsqObjective();
256  computeValuesLsqObjective(lsq_obj_values);
257  computeValuesEquality(eq_values);
258  computeValuesInequality(ineq_values);
259  }
260 
262 
265 
274  virtual int finiteCombinedBoundsDimension();
275 
285  virtual int finiteBoundsDimension();
286 
288 
300  virtual void computeValuesActiveInequality(Eigen::Ref<Eigen::VectorXd> values, double weight = 1.0);
301 
314 
329 
342  virtual void getParametersAndBoundsFinite(Eigen::Ref<Eigen::VectorXd> lb_finite_bounds, Eigen::Ref<Eigen::VectorXd> ub_finite_bounds,
343  Eigen::Ref<Eigen::VectorXd> x_finite_bounds);
344 
346 
348 
355  virtual void computeDenseJacobianLsqObjective(Eigen::Ref<Eigen::MatrixXd> jacobian, const double* multipliers = nullptr);
358  virtual void computeSparseJacobianLsqObjectiveValues(Eigen::Ref<Eigen::VectorXd> values, const double* multipliers = nullptr);
359  virtual void computeSparseJacobianLsqObjective(Eigen::SparseMatrix<double>& jacobian, const double* multipliers = nullptr);
360 
367  virtual void computeDenseJacobianEqualities(Eigen::Ref<Eigen::MatrixXd> jacobian, const double* multipliers = nullptr);
370  virtual void computeSparseJacobianEqualitiesValues(Eigen::Ref<Eigen::VectorXd> values, const double* multipliers = nullptr);
371  virtual void computeSparseJacobianEqualities(Eigen::SparseMatrix<double>& jacobian, const double* multipliers = nullptr);
372 
379  virtual void computeDenseJacobianInequalities(Eigen::Ref<Eigen::MatrixXd> jacobian, const double* multipliers = nullptr);
382  virtual void computeSparseJacobianInequalitiesValues(Eigen::Ref<Eigen::VectorXd> values, const double* multipliers = nullptr);
383  virtual void computeSparseJacobianInequalities(Eigen::SparseMatrix<double>& jacobian, const double* multipliers = nullptr);
384 
391  virtual void computeDenseJacobianActiveInequalities(Eigen::Ref<Eigen::MatrixXd> jacobian, double weight = 1.0);
392  virtual void computeSparseJacobianActiveInequalitiesValues(Eigen::Ref<Eigen::VectorXd> values, double weight = 1.0);
393  virtual void computeSparseJacobianActiveInequalities(Eigen::SparseMatrix<double>& jacobian, double weight = 1.0);
394 
402  virtual void computeDenseJacobianFiniteCombinedBounds(Eigen::Ref<Eigen::MatrixXd> jacobian, double weight = 1.0);
405  virtual void computeSparseJacobianFiniteCombinedBoundsValues(Eigen::Ref<Eigen::VectorXd> values, double weight = 1.0);
406  virtual void computeSparseJacobianFiniteCombinedBounds(Eigen::SparseMatrix<double>& jacobian, double weight = 1.0);
407 
415 
433  virtual void computeDenseJacobians(Eigen::Ref<Eigen::VectorXd> gradient_non_lsq_obj, Eigen::Ref<Eigen::MatrixXd> jacobian_lsq_obj,
434  Eigen::Ref<Eigen::MatrixXd> jacobian_eq, Eigen::Ref<Eigen::MatrixXd> jacobian_ineq,
435  const double* multipliers_lsq_obj = nullptr, const double* multipliers_eq = nullptr,
436  const double* multipliers_ineq = nullptr, bool active_ineq = false, double active_ineq_weight = 1.0);
437  virtual void computeSparseJacobiansNNZ(int& nnz_lsq_obj, int& nnz_eq, int& nnz_ineq);
442  Eigen::Ref<Eigen::VectorXd> values_ineq, const double* multipliers_obj = nullptr,
443  const double* multipliers_eq = nullptr, const double* multipliers_ineq = nullptr,
444  bool active_ineq = false, double active_ineq_weight = 1.0);
445  virtual void computeSparseJacobians(Eigen::SparseMatrix<double>& jacobian_lsq_obj, Eigen::SparseMatrix<double>& jacobian_eq,
446  Eigen::SparseMatrix<double>& jacobian_ineq, const double* multipliers_lsq_obj = nullptr,
447  const double* multipliers_eq = nullptr, const double* multipliers_ineq = nullptr, bool active_ineq = false,
448  double active_ineq_weight = 1.0);
449 
450  virtual void computeCombinedSparseJacobian(Eigen::SparseMatrix<double>& jacobian, bool objective_lsq, bool equality, bool inequality,
451  bool finite_combined_bounds, bool active_ineq = false, double weight_eq = 1.0,
452  double weight_ineq = 1.0, double weight_bounds = 1.0, const Eigen::VectorXd* values = nullptr,
453  const Eigen::VectorXi* col_nnz = nullptr);
454 
455  virtual int computeCombinedSparseJacobiansNNZ(bool objective_lsq = true, bool equality = true, bool inequality = true);
457  bool objective_lsq = true, bool equality = true, bool inequality = true);
458  virtual void computeCombinedSparseJacobiansValues(Eigen::Ref<Eigen::VectorXd> values, bool objective_lsq = true, bool equality = true,
459  bool inequality = true, const double* multipliers_obj = nullptr,
460  const double* multipliers_eq = nullptr, const double* multipliers_ineq = nullptr);
461 
462  // useful for IPOPT (w/ hessian-approx)
464  Eigen::Ref<Eigen::VectorXd> jac_values, bool equality = true,
465  bool inequality = true, const double* multipliers_eq = nullptr,
466  const double* multipliers_ineq = nullptr);
467 
479  const double* multipliers = nullptr, bool jacob_scaled = true);
480 
481  virtual void computeDenseHessianObjective(Eigen::Ref<Eigen::MatrixXd> hessian, double multiplier = 1.0);
482  virtual int computeSparseHessianObjectiveNNZ(bool lower_part_only = false);
484  bool lower_part_only = false);
485  virtual void computeSparseHessianObjectiveValues(Eigen::Ref<Eigen::VectorXd> values, double multiplier = 1.0, bool lower_part_only = false);
486  virtual void computeSparseHessianObjective(Eigen::SparseMatrix<double>& hessian, double multiplier = 1.0);
487 
488  virtual void computeSparseHessianObjectiveLL(Eigen::SparseMatrix<double, Eigen::ColMajor, long long>& H, const Eigen::VectorXi* col_nnz = nullptr,
489  bool upper_part_only = false);
490  virtual void computeSparseHessianObjectiveNNZperCol(Eigen::Ref<Eigen::VectorXi> col_nnz, bool upper_part_only = false);
491 
492  virtual void computeDenseHessianEqualities(Eigen::Ref<Eigen::MatrixXd> hessian, const double* multipliers = nullptr);
493  virtual int computeSparseHessianEqualitiesNNZ(bool lower_part_only = false);
495  bool lower_part_only = false);
496  virtual void computeSparseHessianEqualitiesValues(Eigen::Ref<Eigen::VectorXd> values, const double* multipliers = nullptr,
497  bool lower_part_only = false);
498  virtual void computeSparseHessianEqualities(Eigen::SparseMatrix<double>& hessian, const double* multipliers = nullptr);
499 
500  virtual void computeDenseHessianInequalities(Eigen::Ref<Eigen::MatrixXd> hessian, const double* multipliers = nullptr);
501  virtual int computeSparseHessianInequalitiesNNZ(bool lower_part_only = false);
503  bool lower_part_only = false);
504  virtual void computeSparseHessianInequalitiesValues(Eigen::Ref<Eigen::VectorXd> values, const double* multipliers = nullptr,
505  bool lower_part_only = false);
506  virtual void computeSparseHessianInequalities(Eigen::SparseMatrix<double>& hessian, const double* multipliers = nullptr);
507 
508  // virtual void computeDenseJacobianHessianObjective(Eigen::Ref<Eigen::MatrixXd> jacobian, Eigen::Ref<Eigen::MatrixXd> hessian,
509  // const double* multipliers = nullptr);
510  // virtual void computeDenseJacobianHessianEqualities(Eigen::Ref<Eigen::MatrixXd> jacobian, Eigen::Ref<Eigen::MatrixXd> hessian,
511  // const double* multipliers = nullptr);
512  // virtual void computeDenseJacobianHessianInequalities(Eigen::Ref<Eigen::MatrixXd> jacobian, Eigen::Ref<Eigen::MatrixXd> hessian,
513  // const double* multipliers = nullptr);
514 
516  Eigen::Ref<Eigen::MatrixXd> hessian_ineq, double multiplier_obj = 1.0, const double* multipliers_eq = nullptr,
517  const double* multipliers_ineq = nullptr);
519  Eigen::SparseMatrix<double>& hessian_ineq, double multiplier_obj = 1.0, const double* multipliers_eq = nullptr,
520  const double* multipliers_ineq = nullptr);
521  virtual void computeSparseHessiansNNZ(int& nnz_obj, int& nnz_eq, int& nnz_ineq, bool lower_part_only = false);
525  bool lower_part_only = false);
527  Eigen::Ref<Eigen::VectorXd> values_ineq, double multiplier_obj = 1.0,
528  const double* multipliers_eq = nullptr, const double* multipliers_ineq = nullptr,
529  bool lower_part_only = false);
530 
535  const double* multipliers_ineq, const Eigen::VectorXi* col_nnz = nullptr,
536  bool upper_part_only = false);
538  bool upper_part_only); // TODO(roesmann): lower_part_only? does it make sense for NNZ?
539 
542 
570  const Eigen::VectorXi* col_nnz = nullptr);
571  virtual void computeSparseJacobianTwoSideBoundedLinearFormNNZPerColumn(Eigen::Ref<Eigen::VectorXi> col_nnz, bool include_finite_bounds);
572 
573  virtual int computeSparseJacobianTwoSideBoundedLinearFormNNZ(bool include_finite_bounds);
575  bool include_finite_bounds);
576  virtual void computeSparseJacobianTwoSideBoundedLinearFormValues(Eigen::Ref<Eigen::VectorXd> values, bool include_finite_bounds);
577 
586  Eigen::SparseMatrix<double, Eigen::ColMajor, long long>& H, const double* multipliers_eq, const double* multipliers_ineq,
587  Eigen::SparseMatrix<double, Eigen::ColMajor, long long>& A, bool include_finite_bounds, const Eigen::VectorXi* col_nnz_H = nullptr,
588  const Eigen::VectorXi* col_nnz_A = nullptr, bool upper_part_only_H = false);
589 
615  bool include_finite_bounds);
616 
618 
621  bool include_finite_bounds, const Eigen::VectorXi* col_nnz_H,
622  const Eigen::VectorXi* col_nnz_A, bool upper_part_only_H);
623 
625  virtual bool checkIfAllUnfixedParam(std::function<bool(double, int)> fun);
626 
627  protected:
629 };
630 
631 } // namespace corbo
632 
633 #endif // SRC_OPTIMIZATION_INCLUDE_CORBO_OPTIMIZATION_OPTIMIZATION_PROBLEM_INTERFACE_H_
virtual void computeSparseJacobianTwoSideBoundedLinearFormValues(Eigen::Ref< Eigen::VectorXd > values, bool include_finite_bounds)
virtual double getUpperBound(int idx)=0
Return specific upper bound of a parameter.
virtual void computeDenseHessians(Eigen::Ref< Eigen::MatrixXd > hessian_obj, Eigen::Ref< Eigen::MatrixXd > hessian_eq, Eigen::Ref< Eigen::MatrixXd > hessian_ineq, double multiplier_obj=1.0, const double *multipliers_eq=nullptr, const double *multipliers_ineq=nullptr)
virtual int finiteCombinedBoundsDimension()
Dimension of the set of finite bounds (combined such that each ub and lb component define a single di...
virtual void computeSparseHessiansStructure(Eigen::Ref< Eigen::VectorXi > i_row_obj, Eigen::Ref< Eigen::VectorXi > j_col_obj, Eigen::Ref< Eigen::VectorXi > i_row_eq, Eigen::Ref< Eigen::VectorXi > j_col_eq, Eigen::Ref< Eigen::VectorXi > i_row_ineq, Eigen::Ref< Eigen::VectorXi > j_col_ineq, bool lower_part_only=false)
Scalar * x
virtual void computeSparseJacobianInequalities(Eigen::SparseMatrix< double > &jacobian, const double *multipliers=nullptr)
virtual void computeSparseJacobianLsqObjectiveStructure(Eigen::Ref< Eigen::VectorXi > i_row, Eigen::Ref< Eigen::VectorXi > j_col)
virtual void computeBoundsForTwoSideBoundedLinearForm(Eigen::Ref< Eigen::VectorXd > lbA, Eigen::Ref< Eigen::VectorXd > ubA, bool include_finite_bounds)
Compute lower and upper bounds lbA and ubA for the linear form lbA <= A x <= ubA. ...
virtual void setBounds(const Eigen::Ref< const Eigen::VectorXd > &lb, const Eigen::Ref< const Eigen::VectorXd > &ub)
Set lower and upper bound vector.
virtual void backupParameters()=0
Restore parameter set from the last backup and keep backup if desired.
virtual int getInequalityDimension()=0
Total dimension of general inequality constraints.
virtual void computeDenseJacobianEqualities(Eigen::Ref< Eigen::MatrixXd > jacobian, const double *multipliers=nullptr)
Compute the equality constraint Jacobian Jceq(x) for the current parameter set.
virtual void setParameterValue(int idx, double x)=0
Set specific value of the parameter vector.
virtual int finiteBoundsDimension()
Dimension of the set of finite bounds (individual bounds ub and lb)
virtual void computeSparseHessiansNNZ(int &nnz_obj, int &nnz_eq, int &nnz_ineq, bool lower_part_only=false)
virtual void computeSparseHessians(Eigen::SparseMatrix< double > &hessian_obj, Eigen::SparseMatrix< double > &hessian_eq, Eigen::SparseMatrix< double > &hessian_ineq, double multiplier_obj=1.0, const double *multipliers_eq=nullptr, const double *multipliers_ineq=nullptr)
virtual void computeSparseJacobianTwoSideBoundedLinearFormStructure(Eigen::Ref< Eigen::VectorXi > i_row, Eigen::Ref< Eigen::VectorXi > j_col, bool include_finite_bounds)
virtual int computeSparseHessianObjectiveNNZ(bool lower_part_only=false)
virtual void computeSparseJacobianFiniteCombinedBounds(Eigen::SparseMatrix< double > &jacobian, double weight=1.0)
virtual void computeDenseJacobianFiniteCombinedBounds(Eigen::Ref< Eigen::MatrixXd > jacobian, double weight=1.0)
Compute the Jacobian for finite combined bounds.
virtual void computeDenseJacobianActiveInequalities(Eigen::Ref< Eigen::MatrixXd > jacobian, double weight=1.0)
Compute the Jacobian Jc(x) with non-zeros for active constraints c(x)>= 0 and zeros for inactive ones...
virtual void computeDistanceFiniteCombinedBounds(Eigen::Ref< Eigen::VectorXd > values)
Compute the distance to finite bound values (combined lower and upper)
virtual void computeSparseHessianEqualitiesValues(Eigen::Ref< Eigen::VectorXd > values, const double *multipliers=nullptr, bool lower_part_only=false)
virtual void restoreBackupParameters(bool keep_backup)=0
Discard last backup (or all)
virtual void computeCombinedSparseJacobiansValues(Eigen::Ref< Eigen::VectorXd > values, bool objective_lsq=true, bool equality=true, bool inequality=true, const double *multipliers_obj=nullptr, const double *multipliers_eq=nullptr, const double *multipliers_ineq=nullptr)
virtual void computeSparseJacobianEqualities(Eigen::SparseMatrix< double > &jacobian, const double *multipliers=nullptr)
virtual void computeSparseJacobianFiniteCombinedBoundsStructure(Eigen::Ref< Eigen::VectorXi > i_row, Eigen::Ref< Eigen::VectorXi > j_col)
virtual int getEqualityDimension()=0
Total dimension of equality constraints.
virtual void computeSparseJacobianTwoSideBoundedLinearFormAndHessianObjective(Eigen::SparseMatrix< double, Eigen::ColMajor, long long > &H, Eigen::SparseMatrix< double, Eigen::ColMajor, long long > &A, bool include_finite_bounds, const Eigen::VectorXi *col_nnz_H, const Eigen::VectorXi *col_nnz_A, bool upper_part_only_H)
virtual void computeSparseJacobianLsqObjectiveValues(Eigen::Ref< Eigen::VectorXd > values, const double *multipliers=nullptr)
std::shared_ptr< OptimizationProblemInterface > Ptr
virtual void computeSparseJacobianEqualitiesValues(Eigen::Ref< Eigen::VectorXd > values, const double *multipliers=nullptr)
virtual void computeSparseJacobianLsqObjective(Eigen::SparseMatrix< double > &jacobian, const double *multipliers=nullptr)
virtual void computeDenseJacobianFiniteCombinedBoundsIdentity(Eigen::Ref< Eigen::MatrixXd > jacobian)
Compute the Jacobian for finite combined bounds.
virtual void applyIncrement(const Eigen::Ref< const Eigen::VectorXd > &increment)
Apply increment to the current parameter set.
virtual void computeSparseHessianObjectiveStructure(Eigen::Ref< Eigen::VectorXi > i_row, Eigen::Ref< Eigen::VectorXi > j_col, bool lower_part_only=false)
MatrixType A(a, *n, *n, *lda)
virtual void computeSparseHessianInequalitiesValues(Eigen::Ref< Eigen::VectorXd > values, const double *multipliers=nullptr, bool lower_part_only=false)
virtual void setLowerBound(int idx, double lb)=0
Set specific lower bound of a parameter.
virtual int computeCombinedSparseJacobiansNNZ(bool objective_lsq=true, bool equality=true, bool inequality=true)
virtual int getLsqObjectiveDimension()=0
Total dimension of least-squares objective function terms.
virtual double getParameterValue(int idx)=0
Return specific value of the parameter vector.
virtual void computeSparseHessianObjective(Eigen::SparseMatrix< double > &hessian, double multiplier=1.0)
virtual int computeSparseHessianInequalitiesNNZ(bool lower_part_only=false)
virtual void setUpperBound(int idx, double ub)=0
Set specific upper bound of a parameter.
virtual void computeSparseHessianLagrangianNNZperCol(Eigen::Ref< Eigen::VectorXi > col_nnz, bool upper_part_only)
virtual void computeLowerAndUpperBoundDiff(Eigen::Ref< Eigen::VectorXd > lb_minus_x, Eigen::Ref< Eigen::VectorXd > ub_minus_x)
Compute the distance between parameters and bounds.
virtual void computeValues(double &non_lsq_obj_value, Eigen::Ref< Eigen::VectorXd > lsq_obj_values, Eigen::Ref< Eigen::VectorXd > eq_values, Eigen::Ref< Eigen::VectorXd > ineq_values)
virtual int computeSparseHessianEqualitiesNNZ(bool lower_part_only=false)
virtual void computeValuesActiveInequality(Eigen::Ref< Eigen::VectorXd > values, double weight=1.0)
Compute the values of the active inequality constraints (elementwise max(0, c(x))) ...
virtual void computeSparseJacobianTwoSideBoundedLinearFormAndHessianLagrangian(Eigen::SparseMatrix< double, Eigen::ColMajor, long long > &H, const double *multipliers_eq, const double *multipliers_ineq, Eigen::SparseMatrix< double, Eigen::ColMajor, long long > &A, bool include_finite_bounds, const Eigen::VectorXi *col_nnz_H=nullptr, const Eigen::VectorXi *col_nnz_A=nullptr, bool upper_part_only_H=false)
Compute the Jacobian and Hessian of the lagrangian.
virtual void computeSparseJacobianInequalitiesStructure(Eigen::Ref< Eigen::VectorXi > i_row, Eigen::Ref< Eigen::VectorXi > j_col)
virtual void computeSparseJacobiansNNZ(int &nnz_lsq_obj, int &nnz_eq, int &nnz_ineq)
std::unique_ptr< OptimizationProblemInterface > UPtr
virtual void computeSparseJacobianEqualitiesStructure(Eigen::Ref< Eigen::VectorXi > i_row, Eigen::Ref< Eigen::VectorXi > j_col)
virtual void computeSparseHessianObjectiveValues(Eigen::Ref< Eigen::VectorXd > values, double multiplier=1.0, bool lower_part_only=false)
virtual void computeSparseJacobians(Eigen::SparseMatrix< double > &jacobian_lsq_obj, Eigen::SparseMatrix< double > &jacobian_eq, Eigen::SparseMatrix< double > &jacobian_ineq, const double *multipliers_lsq_obj=nullptr, const double *multipliers_eq=nullptr, const double *multipliers_ineq=nullptr, bool active_ineq=false, double active_ineq_weight=1.0)
Generic interface for optimization problem definitions.
virtual void getParameterVector(Eigen::Ref< Eigen::VectorXd > x)
Return deep copy of the complete parameter vector.
virtual double computeValueNonLsqObjective()=0
virtual void computeSparseHessianObjectiveLL(Eigen::SparseMatrix< double, Eigen::ColMajor, long long > &H, const Eigen::VectorXi *col_nnz=nullptr, bool upper_part_only=false)
virtual void computeValuesEquality(Eigen::Ref< Eigen::VectorXd > values)=0
Compute the equality constraint values ceq(x) for the current parameter set.
virtual void computeSparseJacobianInequalitiesValues(Eigen::Ref< Eigen::VectorXd > values, const double *multipliers=nullptr)
virtual void getBounds(Eigen::Ref< Eigen::VectorXd > lb, Eigen::Ref< Eigen::VectorXd > ub)
Get lower and upper bound vector.
virtual double getLowerBound(int idx)=0
Return specific lower bound value of a parameter.
virtual void computeValuesInequality(Eigen::Ref< Eigen::VectorXd > values)=0
Compute the inequality constraint values c(x) for the current parameter set.
virtual int computeSparseJacobianTwoSideBoundedLinearFormNNZ(bool include_finite_bounds)
A matrix or vector expression mapping an existing expression.
Definition: Ref.h:192
virtual void computeGradientNonLsqObjective(Eigen::Ref< Eigen::VectorXd > gradient)
virtual void computeSparseJacobianFiniteCombinedBoundsValues(Eigen::Ref< Eigen::VectorXd > values, double weight=1.0)
virtual void computeGradientObjectiveAndCombinedSparseJacobiansValues(Eigen::Ref< Eigen::VectorXd > gradient, Eigen::Ref< Eigen::VectorXd > jac_values, bool equality=true, bool inequality=true, const double *multipliers_eq=nullptr, const double *multipliers_ineq=nullptr)
virtual void computeSparseJacobianTwoSideBoundedLinearFormNNZPerColumn(Eigen::Ref< Eigen::VectorXi > col_nnz, bool include_finite_bounds)
virtual void computeSparseJacobianActiveInequalities(Eigen::SparseMatrix< double > &jacobian, double weight=1.0)
virtual int getNonLsqObjectiveDimension()=0
Get dimension of the objective (should be zero or one)
virtual void computeDenseHessianInequalities(Eigen::Ref< Eigen::MatrixXd > hessian, const double *multipliers=nullptr)
virtual void computeSparseJacobiansStructure(Eigen::Ref< Eigen::VectorXi > i_row_obj, Eigen::Ref< Eigen::VectorXi > j_col_obj, Eigen::Ref< Eigen::VectorXi > i_row_eq, Eigen::Ref< Eigen::VectorXi > j_col_eq, Eigen::Ref< Eigen::VectorXi > i_row_ineq, Eigen::Ref< Eigen::VectorXi > j_col_ineq)
virtual void computeSparseJacobiansValues(Eigen::Ref< Eigen::VectorXd > values_obj, Eigen::Ref< Eigen::VectorXd > values_eq, Eigen::Ref< Eigen::VectorXd > values_ineq, const double *multipliers_obj=nullptr, const double *multipliers_eq=nullptr, const double *multipliers_ineq=nullptr, bool active_ineq=false, double active_ineq_weight=1.0)
virtual void computeSparseHessiansValues(Eigen::Ref< Eigen::VectorXd > values_obj, Eigen::Ref< Eigen::VectorXd > values_eq, Eigen::Ref< Eigen::VectorXd > values_ineq, double multiplier_obj=1.0, const double *multipliers_eq=nullptr, const double *multipliers_ineq=nullptr, bool lower_part_only=false)
virtual void computeDenseHessianEqualities(Eigen::Ref< Eigen::MatrixXd > hessian, const double *multipliers=nullptr)
virtual bool isLeastSquaresProblem() const =0
Check if the underlying problem is defined in the least squares form.
virtual void computeSparseHessianLagrangian(Eigen::SparseMatrix< double, Eigen::ColMajor, long long > &H, const double *multipliers_eq, const double *multipliers_ineq, const Eigen::VectorXi *col_nnz=nullptr, bool upper_part_only=false)
Compute the hessian of the lagrangian L(x) = f(x) + lambda1 * c(x) + lambda2 * ceq(x) ...
virtual void getParametersAndBoundsFinite(Eigen::Ref< Eigen::VectorXd > lb_finite_bounds, Eigen::Ref< Eigen::VectorXd > ub_finite_bounds, Eigen::Ref< Eigen::VectorXd > x_finite_bounds)
Return bound and parameter vectors only for finite boudns.
virtual void computeValuesLsqObjective(Eigen::Ref< Eigen::VectorXd > values)=0
Compute the objective function values f(x) for the current parameter set.
virtual int getObjectiveDimension()=0
Get dimension of the objective (should be zero or one, includes Lsq objectives if present) ...
virtual void discardBackupParameters(bool all=false)=0
virtual void computeSparseHessianEqualitiesStructure(Eigen::Ref< Eigen::VectorXi > i_row, Eigen::Ref< Eigen::VectorXi > j_col, bool lower_part_only=false)
virtual int getParameterDimension()=0
Effictive dimension of the optimization parameter set (changeable, non-fixed part) ...
virtual void computeCombinedSparseJacobiansStructure(Eigen::Ref< Eigen::VectorXi > i_row, Eigen::Ref< Eigen::VectorXi > j_col, bool objective_lsq=true, bool equality=true, bool inequality=true)
virtual void computeDenseJacobianInequalities(Eigen::Ref< Eigen::MatrixXd > jacobian, const double *multipliers=nullptr)
Compute the inequality constraint Jacobian Jc(x) for the current parameter set.
virtual void computeSparseHessianInequalitiesStructure(Eigen::Ref< Eigen::VectorXi > i_row, Eigen::Ref< Eigen::VectorXi > j_col, bool lower_part_only=false)
virtual void computeSparseHessianObjectiveNNZperCol(Eigen::Ref< Eigen::VectorXi > col_nnz, bool upper_part_only=false)
virtual void setParameterVector(const Eigen::Ref< const Eigen::VectorXd > &x)
Set complete parameter vector.
virtual bool checkIfAllUnfixedParam(std::function< bool(double, int)> fun)
Check if a function taking the parameter value and unfixed-idx is true for all unfixed parameter valu...
virtual void computeDenseHessianObjective(const Eigen::Ref< const Eigen::MatrixXd > &jacobian, Eigen::Ref< Eigen::MatrixXd > hessian, const double *multipliers=nullptr, bool jacob_scaled=true)
Compute the objective Hessian Hf(x) for the current parameter set.
virtual void computeSparseJacobianTwoSideBoundedLinearForm(Eigen::SparseMatrix< double, Eigen::ColMajor, long long > &A, bool include_finite_bounds, const Eigen::VectorXi *col_nnz=nullptr)
Compute the jacobian A for the linear form lbA <= A x <= lbB.
virtual void computeDenseJacobianLsqObjective(Eigen::Ref< Eigen::MatrixXd > jacobian, const double *multipliers=nullptr)
Compute the objective Jacobian Jf(x) for the current parameter set.
virtual void computeSparseHessianEqualities(Eigen::SparseMatrix< double > &hessian, const double *multipliers=nullptr)
virtual void computeGradientObjective(Eigen::Ref< Eigen::VectorXd > gradient)
virtual void computeSparseHessianInequalities(Eigen::SparseMatrix< double > &hessian, const double *multipliers=nullptr)
virtual void computeDenseJacobians(Eigen::Ref< Eigen::VectorXd > gradient_non_lsq_obj, Eigen::Ref< Eigen::MatrixXd > jacobian_lsq_obj, Eigen::Ref< Eigen::MatrixXd > jacobian_eq, Eigen::Ref< Eigen::MatrixXd > jacobian_ineq, const double *multipliers_lsq_obj=nullptr, const double *multipliers_eq=nullptr, const double *multipliers_ineq=nullptr, bool active_ineq=false, double active_ineq_weight=1.0)
Compute the objective and constraint Jacobians at once.
virtual void computeSparseJacobianActiveInequalitiesValues(Eigen::Ref< Eigen::VectorXd > values, double weight=1.0)
virtual void computeCombinedSparseJacobian(Eigen::SparseMatrix< double > &jacobian, bool objective_lsq, bool equality, bool inequality, bool finite_combined_bounds, bool active_ineq=false, double weight_eq=1.0, double weight_ineq=1.0, double weight_bounds=1.0, const Eigen::VectorXd *values=nullptr, const Eigen::VectorXi *col_nnz=nullptr)


control_box_rst
Author(s): Christoph Rösmann
autogenerated on Mon Feb 28 2022 22:07:06