TVM  0.9.4
BasicLinearFunction.h
Go to the documentation of this file.
1 
3 #pragma once
4 
8 #include <tvm/utils/AffineExpr.h>
9 
10 #include <utility>
11 
12 namespace tvm
13 {
14 
15 namespace function
16 {
17 
22 {
23 public:
27  BasicLinearFunction(const std::vector<MatrixConstRef> & A, const std::vector<VariablePtr> & x);
28 
32  BasicLinearFunction(const std::vector<MatrixConstRef> & A,
33  const std::vector<VariablePtr> & x,
34  const VectorConstRef & b);
35 
41 
46  BasicLinearFunction(int m, const std::vector<VariablePtr> & x);
47 
51  template<typename Derived>
53 
57  template<typename CstDerived, typename... Derived>
59 
62  virtual void A(const MatrixConstRef & A, const Variable & x, const tvm::internal::MatrixProperties & p = {});
64  virtual void A(const MatrixConstRef & A, const tvm::internal::MatrixProperties & p = {});
66  virtual void b(const VectorConstRef & b, const tvm::internal::MatrixProperties & p = {});
67 
68  using LinearFunction::b;
69 
70 private:
71  template<typename Derived>
72  void add(const Eigen::MatrixBase<Derived> & A, VariablePtr x);
73 
74  template<typename Tuple, size_t... Indices>
75  void add(Tuple && tuple, std::index_sequence<Indices...>);
76 
77  template<typename Derived>
78  void add(const utils::LinearExpr<Derived> & lin);
79 };
80 
81 namespace internal
82 {
83 template<typename Tuple, size_t... Indices>
84 void addVar(tvm::internal::VariableCountingVector & v, Tuple && tuple, std::index_sequence<Indices...>)
85 { (v.add(std::get<Indices>(std::forward<Tuple>(tuple)).variable()), ...); }
86 } // namespace internal
87 
88 template<typename Derived>
90 : LinearFunction(static_cast<int>(lin.matrix().rows()))
91 {
92  add(lin);
93  this->b(Eigen::VectorXd::Zero(size()),
96 }
97 
98 template<typename CstDerived, typename... Derived>
100 : LinearFunction(static_cast<int>(std::get<0>(aff.linear()).matrix().rows()))
101 {
102  using Indices = std::make_index_sequence<sizeof...(Derived)>;
104  internal::addVar(v, aff.linear(), Indices{});
105  const auto & vars = v.variables();
106  const auto & simple = v.simple(); // simple[i]==true means vars[i] that appears in v based on a single first add to v
107  for(int i = 0; i < vars.numberOfVariables(); ++i)
108  {
109  if(!simple[i]) // if vars[i] is the results of adding several potentially overlapping variables to v...
110  {
111  addVariable(vars[i], true); // ... we add the variable now and initialize its jacobian to 0 below.
112  // This will force using += to fill the part of the jacobian corresponding to
113  // a subvariable in method add. Simple variables will be added directly in add.
114  }
115  }
116  for(auto & J : jacobian_)
117  {
118  J.second.setZero();
119  J.second.properties({tvm::internal::MatrixProperties::Constness(true)});
120  }
121  add(aff.linear(), Indices{});
122  if constexpr(std::is_same_v<CstDerived, utils::internal::NoConstant>)
123  {
124  this->b(Eigen::VectorXd::Zero(size()),
126  }
127  else
128  {
129  this->b(aff.constant(), {tvm::internal::MatrixProperties::Constness(true)});
130  }
132 }
133 
134 template<typename Derived>
135 inline void BasicLinearFunction::add(const Eigen::MatrixBase<Derived> & A, VariablePtr x)
136 {
137  if(!x->space().isEuclidean() && x->isBasePrimitive())
138  throw std::runtime_error("We allow linear function only on Euclidean variables.");
139  if(A.rows() != size())
140  throw std::runtime_error("Matrix A doesn't have coherent row size.");
141  if(A.cols() != x->size())
142  throw std::runtime_error("Matrix A doesn't have its column size coherent with its corresponding variable.");
143 
144  if(variables().contains(*x))
145  {
146  auto J = jacobian_.at(x.get(), tvm::utils::internal::with_sub{});
147  J += A; // this erase the properties of J;
148  J.properties({tvm::internal::MatrixProperties::Constness(true)}); // We restore only this one
149  }
150  else
151  {
153  addVariable(x, true);
154  jacobian_.at(x.get()) = A;
155  Shape shape = Shape::GENERAL;
156  // TODO: generalize test to matrices with fixed size, to triangular matrices
157  // and to richer diagonal case
158  if constexpr(std::is_same_v<Derived, std::remove_const_t<tvm::utils::internal::IdentityType>>)
159  {
160  shape = Shape::IDENTITY;
161  }
162  if constexpr(std::is_same_v<Derived, std::remove_const_t<tvm::utils::internal::MinusIdentityType>>)
163  {
164  shape = Shape::MINUS_IDENTITY;
165  }
166  if constexpr(std::is_same_v<Derived, std::remove_const_t<tvm::utils::internal::MultIdentityType>>)
167  {
168  shape = Shape::MULTIPLE_OF_IDENTITY;
169  }
170  if constexpr(std::is_same_v<Derived, std::remove_const_t<decltype(Eigen::VectorXd(1).asDiagonal())>>)
171  {
172  shape = Shape::DIAGONAL;
173  }
174  jacobian_.at(x.get()).properties({tvm::internal::MatrixProperties::Constness(true), shape});
175  }
176 }
177 
178 template<typename Tuple, size_t... Indices>
179 inline void BasicLinearFunction::add(Tuple && tuple, std::index_sequence<Indices...>)
180 { (add(std::get<Indices>(std::forward<Tuple>(tuple))), ...); }
181 
182 template<typename Derived>
183 inline void BasicLinearFunction::add(const utils::LinearExpr<Derived> & lin)
184 { add(lin.matrix(), lin.variable()); }
185 
186 } // namespace function
187 
188 } // namespace tvm
#define TVM_DLLAPI
Definition: api.h:35
Definition: Variable.h:49
Definition: BasicLinearFunction.h:22
virtual void b(const VectorConstRef &b, const tvm::internal::MatrixProperties &p={})
BasicLinearFunction(int m, VariablePtr x)
BasicLinearFunction(const std::vector< MatrixConstRef > &A, const std::vector< VariablePtr > &x)
BasicLinearFunction(const MatrixConstRef &A, VariablePtr x)
BasicLinearFunction(int m, const std::vector< VariablePtr > &x)
virtual void A(const MatrixConstRef &A, const Variable &x, const tvm::internal::MatrixProperties &p={})
BasicLinearFunction(const std::vector< MatrixConstRef > &A, const std::vector< VariablePtr > &x, const VectorConstRef &b)
BasicLinearFunction(const MatrixConstRef &A, VariablePtr x, const VectorConstRef &b)
virtual void A(const MatrixConstRef &A, const tvm::internal::MatrixProperties &p={})
Definition: LinearFunction.h:62
const tvm::internal::VectorWithProperties & b() const
Definition: LinearFunction.h:70
int size() const
Definition: FirstOrderProvider.h:200
utils::internal::MapWithVariableAsKey< MatrixWithProperties, slice_matrix, true > jacobian_
Definition: FirstOrderProvider.h:179
const VariableVector & variables() const
Definition: FirstOrderProvider.h:205
void addVariable(VariablePtr v, bool linear)
Definition: MatrixProperties.h:17
Shape
Definition: MatrixProperties.h:21
@ ZERO
Definition: MatrixProperties.h:39
Definition: VariableCountingVector.h:30
const VariableVector & variables() const
const std::vector< uint8_t > simple() const
Definition: AffineExpr.h:71
const CstDerived & constant() const
Definition: AffineExpr.h:98
const std::tuple< LinearExpr< Derived >... > & linear() const
Definition: AffineExpr.h:97
Definition: AffineExpr.h:26
Definition: probe.h:44
void addVar(tvm::internal::VariableCountingVector &v, Tuple &&tuple, std::index_sequence< Indices... >)
Definition: BasicLinearFunction.h:84
@ DIAGONAL
Definition: CompiledAssignment.h:38
@ IDENTITY
Definition: CompiledAssignment.h:47
@ GENERAL
Definition: CompiledAssignment.h:49
Definition: Clock.h:12
std::shared_ptr< Variable > VariablePtr
Definition: defs.h:65
Eigen::Ref< const Eigen::VectorXd > VectorConstRef
Definition: defs.h:50
Eigen::Ref< const Eigen::MatrixXd > MatrixConstRef
Definition: defs.h:48
Definition: MatrixProperties.h:65