20 class SlicedTube :
public SlicedTubeBase
24 explicit SlicedTube(
const std::shared_ptr<TDomain>& tdomain,
26 : SlicedTubeBase(tdomain)
28 for(
auto it = _tdomain->begin(); it != _tdomain->end(); ++it)
31 std::make_shared<Slice<T>>(*this, it, codomain)
35 explicit SlicedTube(
const std::shared_ptr<TDomain>& tdomain,
36 const AnalyticFunction<
typename ExprType<T>::Type>& f)
37 : SlicedTubeBase(tdomain)
39 assert_release(f.args().size() == 1
40 &&
"function's inputs must be limited to one system variable");
42 for(
auto it = _tdomain->begin(); it != _tdomain->end(); ++it)
45 std::make_shared<Slice<T>>(*
this, it, f.eval((Interval)*it))
50 requires std::is_same_v<typename Wrapper<V>::Domain,T>
51 explicit SlicedTube(
const std::shared_ptr<TDomain>& tdomain,
52 const SampledTraj<V>& f)
53 : SlicedTubeBase(tdomain)
55 for(
auto it = _tdomain->begin(); it != _tdomain->end(); ++it)
58 std::make_shared<Slice<T>>(*this, it, f((Interval)*it))
62 SlicedTube(
const SlicedTube<T>& x)
63 : SlicedTubeBase(x.tdomain())
65 for(
auto it = _tdomain->begin(); it != _tdomain->end(); ++it)
68 std::make_shared<Slice<T>>(*x.slice(it), *this)
72 inline SlicedTube& operator=(
const SlicedTube& x)
74 assert_release(_tdomain == x._tdomain);
76 for(
auto it = _tdomain->begin(); it != _tdomain->end(); ++it)
77 slice(it)->set(x.slice(it)->codomain(),
false);
82 inline Index size()
const
84 return first_slice()->size();
87 virtual std::pair<Index,Index> shape()
const
89 if constexpr(std::is_same_v<typename ExprType<T>::Type,ScalarType>)
92 return {first_slice()->codomain().rows(),first_slice()->codomain().cols()};
95 inline double volume()
const
98 for(
const auto& s : *
this)
100 volume += s.volume();
104 inline std::shared_ptr<Slice<T>> first_slice()
106 return std::const_pointer_cast<Slice<T>>(
107 static_cast<const SlicedTube&
>(*this).first_slice());
110 inline std::shared_ptr<const Slice<T>> first_slice()
const
112 return std::static_pointer_cast<const Slice<T>>(
113 this->SlicedTubeBase::first_slice());
116 inline std::shared_ptr<Slice<T>> last_slice()
118 return std::const_pointer_cast<Slice<T>>(
119 static_cast<const SlicedTube&
>(*this).last_slice());
122 inline std::shared_ptr<const Slice<T>> last_slice()
const
124 return std::static_pointer_cast<const Slice<T>>(
125 this->SlicedTubeBase::last_slice());
128 inline std::shared_ptr<Slice<T>> slice(
const std::list<TSlice>::iterator& it)
130 return std::const_pointer_cast<Slice<T>>(
131 static_cast<const SlicedTube&
>(*this).slice(it));
134 inline std::shared_ptr<const Slice<T>> slice(
const std::list<TSlice>::const_iterator& it)
const
136 return std::static_pointer_cast<const Slice<T>>(
137 it->slices().at(
this));
140 inline std::shared_ptr<Slice<T>> slice(
const std::list<TSlice>::reverse_iterator& it)
142 return std::const_pointer_cast<Slice<T>>(
143 static_cast<const SlicedTube&
>(*this).slice(it));
146 inline std::shared_ptr<const Slice<T>> slice(
const std::list<TSlice>::const_reverse_iterator& it)
const
148 return std::static_pointer_cast<const Slice<T>>(
149 it->slices().at(
this));
152 inline std::shared_ptr<Slice<T>> slice(std::shared_ptr<TSlice> ptr)
155 auto it = std::find_if(_tdomain->begin(), _tdomain->end(),
156 [&](TSlice& t){ return &t == ptr.get(); });
160 inline bool is_empty()
const
164 for(
const auto& s : *
this)
165 if(s.is_gate() && s.is_empty())
167 for(
const auto& s : *
this)
168 if(!s.is_gate() && s.is_empty())
173 inline bool is_unbounded()
const
175 for(
const auto& s : *
this)
181 inline T codomain()
const
183 T x = first_slice()->codomain();
184 for(
const auto& s : *
this)
189 template<
typename Func>
190 inline T eval_common(
const Interval& t,
const Func& apply_eval)
const
193 return all_reals_value();
195 auto it = _tdomain->tslice(t.lb());
196 assert(it != _tdomain->end());
197 T codomain = apply_eval(it, t & *it);
199 while(it != std::next(_tdomain->tslice(t.ub())))
201 if(it->lb() == t.ub())
break;
202 codomain |= apply_eval(it, t & *it);
209 T operator()(
const Interval& t)
const
211 return eval_common(t,
212 [
this](
auto it,
const Interval& t_) {
213 return slice(it)->operator()(t_);
217 T operator()(
const Interval& t,
const SlicedTube<T>& v)
const
218 requires (std::is_same_v<T,Interval> || std::is_same_v<T,IntervalVector>)
220 return eval_common(t,
221 [
this,&v](
auto it,
const Interval& t_) {
222 return slice(it)->operator()(t_, *v.slice(it));
226 std::pair<T,T> enclosed_bounds(
const Interval& t)
const
228 auto x = this->empty_value();
229 auto bounds = std::make_pair(x,x);
231 if(t.lb() < _tdomain->t0_tf().lb() || t.ub() > _tdomain->t0_tf().ub())
233 x.init(Interval(-oo,0));
235 x.init(Interval(0,oo));
239 Interval t_inter = t & _tdomain->t0_tf();
240 auto si = slice(_tdomain->tslice(t_inter.lb()));
242 while(si && si->t0_tf().lb() <= t_inter.ub())
244 auto slice_bounds = si->enclosed_bounds(t_inter & si->t0_tf());
245 bounds.first |= slice_bounds.first;
246 bounds.second |= slice_bounds.second;
247 si = si->next_slice();
253 inline void set(
const T& codomain)
255 assert_release(codomain.size() == this->size());
259 s.set(codomain,
false);
263 s.set(codomain,
false);
266 inline void set(
const T& codomain,
double t)
268 assert_release(codomain.size() == this->size());
269 slice(_tdomain->sample(t,
true))->set(codomain);
272 inline void set(
const T& codomain,
const Interval& t)
274 auto it_lb = _tdomain->sample(t.lb(), t.is_degenerated());
275 std::list<TSlice>::iterator it_ub;
277 if(!t.is_degenerated())
279 it_ub = _tdomain->sample(t.ub(),
false);
281 if(it_ub->lb() == t.ub())
284 if(it_lb->ub() == t.lb())
293 slice(it_lb)->set(codomain);
294 }
while(it_lb != it_ub && (++it_lb) != _tdomain->end());
297 inline void set_ith_slice(
const T& codomain, Index i)
300 for(
auto& si : *
this)
309 requires (std::is_same_v<typename Wrapper<V>::Domain,T> || std::is_same_v<V,double>)
310 const SlicedTube<T>& inflate(
const V&
rad)
315 s.set(T(s.codomain()).inflate(
rad),
false);
319 s.set(T(s.codomain()).inflate(
rad),
false);
325 requires (std::is_same_v<typename Wrapper<V>::Domain,T> || std::is_same_v<V,double>)
326 const SlicedTube<T>& inflate(
const SampledTraj<V>&
rad)
329 assert_release(tdomain()->t0_tf().
is_subset(
rad.tdomain()));
333 s.set(T(s.codomain()).inflate(
rad(s.t0_tf()).ub()),
false);
337 s.set(T(s.codomain()).inflate(
rad(s.t0_tf()).ub()),
false);
342 SlicedTube<Interval> operator[](Index i)
const
344 assert_release(i >= 0 && i < size());
345 SlicedTube<Interval> xi(tdomain(), Interval());
346 for(
auto it = tdomain()->begin() ; it != tdomain()->end() ; it++)
347 xi.slice(it)->codomain() = slice(it)->codomain()[i];
351 SlicedTube<IntervalVector> subvector(Index i, Index j)
const
353 assert_release(i >= 0 && i <= j && j < size());
355 for(
auto it = tdomain()->begin() ; it != tdomain()->end() ; it++)
356 xij.slice(it)->codomain() = slice(it)->codomain().subvector(i,j);
360 inline bool operator==(
const SlicedTube& x)
const
362 if(!TDomain::are_same(tdomain(), x.tdomain()))
365 auto it_this = _tdomain->begin();
366 auto it_x = x.tdomain()->cbegin();
368 while(it_this != _tdomain->end())
369 if(*slice(it_this++) != *x.slice(it_x++))
375 inline SlicedTube operator&=(
const SlicedTube& x)
377 assert(TDomain::are_same(tdomain(), x.tdomain()));
378 auto it_this = _tdomain->begin();
379 auto it_x = x.tdomain()->cbegin();
381 while(it_this != _tdomain->end())
383 auto s = slice(it_this);
384 s->set(s->codomain() & x.slice(it_x)->codomain());
388 assert(it_x == x.tdomain()->cend());
392 friend inline std::ostream& operator<<(std::ostream& os,
const SlicedTube<T>& x)
395 <<
"↦" << x.codomain()
396 <<
", " << x.nb_slices()
397 <<
" slice" << (x.nb_slices() > 1 ?
"s" :
"")
402 AnalyticFunction<typename ExprType<T>::Type> as_function()
const
406 AnalyticExprWrapper<typename ExprType<T>::Type>(
407 std::make_shared<AnalyticOperationExpr<
408 TubeOp<SlicedTube<T>>,
typename ExprType<T>::Type,ScalarType>>(*
this,t))
412 template<
typename Func>
413 Interval invert_common(
const T& y,
const Interval& t,
const Func& apply_invert)
const
415 assert_release(y.size() == this->size());
417 if(t.is_empty() || y.is_empty())
420 else if(!t.is_subset(_tdomain->t0_tf()))
424 Interval t_ = t & _tdomain->t0_tf();
426 assert(!t_.is_empty() && !t_.is_unbounded());
427 for(
auto it = _tdomain->tslice(t_.lb()) ;
428 it != _tdomain->end() && it->lb() < t_.ub() ; it++)
433 assert(it != _tdomain->end());
434 invert |= apply_invert(it, t & *it);
440 template<
typename Func>
441 void invert_common_subsets(
const T& y, std::vector<Interval> &v_t,
442 const Interval& t,
const Func& apply_invert)
const
444 assert_release(y.size() == this->size());
447 if(t.is_empty() || y.is_empty())
450 else if(!t.is_subset(_tdomain->t0_tf()))
452 v_t.push_back(Interval());
457 Interval t_ = t & _tdomain->t0_tf();
459 assert(!t_.is_empty() && !t_.is_unbounded());
460 for(
auto it = _tdomain->tslice(t_.lb()) ;
461 it != _tdomain->end() && it->lb() < t_.ub() ; it++)
466 assert(it != _tdomain->end());
467 Interval local_invert = apply_invert(it, t & *it);
469 if(local_invert.is_empty() && !invert.is_empty())
471 v_t.push_back(invert);
476 invert |= local_invert;
479 if(!invert.is_empty())
480 v_t.push_back(invert);
492 Interval invert(
const T& y,
const Interval& t = Interval())
const
494 return invert_common(y, t,
495 [
this,&y](
auto it,
const Interval& t_) {
496 return slice(it)->invert(y, t_);
508 void invert(
const T& y, std::vector<Interval> &v_t,
const Interval& t = Interval())
const
510 return invert_common_subsets(y, v_t, t,
511 [
this,&y](
auto it,
const Interval& t_) {
512 return slice(it)->invert(y, t_);
527 Interval invert(
const T& y,
const SlicedTube<T>& v,
const Interval& t = Interval()) const
528 requires (std::is_same_v<T,Interval> || std::is_same_v<T,
IntervalVector>)
530 return invert_common(y, t,
531 [
this,&v,&y](
auto it,
const Interval& t_) {
532 return slice(it)->invert(y, *v.slice(it), t_);
547 void invert(
const T& y, std::vector<Interval> &v_t,
const SlicedTube<T>& v,
const Interval& t = Interval()) const
548 requires (std::is_same_v<T,Interval> || std::is_same_v<T,
IntervalVector>)
550 return invert_common_subsets(y, v_t, t,
551 [
this,&v,&y](
auto it,
const Interval& t_) {
552 return slice(it)->invert(y, *v.slice(it), t_);
558 T integral(
const Interval& t)
const;
559 T integral(
const Interval& t1,
const Interval& t2)
const;
560 std::pair<T,T> partial_integral(
const Interval& t)
const;
561 std::pair<T,T> partial_integral(
const Interval& t1,
const Interval& t2)
const;
563 inline SlicedTube<T> primitive()
const
565 auto x0 = all_reals_value();
567 return primitive(x0);
570 inline SlicedTube<T> primitive(
const T& x0)
const
572 auto x = all_reals_value();
573 auto p = SlicedTube<T>(this->tdomain(), x);
574 p.set(x0, this->tdomain()->t0_tf().
lb());
580 inline T all_reals_value()
const
582 T x = first_slice()->codomain();
587 inline T empty_value()
const
589 T x = first_slice()->codomain();
597 using base_container = std::list<TSlice>;
599 struct iterator :
public base_container::iterator
603 iterator(SlicedTube& x, base_container::iterator it)
604 : base_container::iterator(it), _x(x) { }
606 std::shared_ptr<Slice<T>> operator->()
608 return _x.slice(*
this);
611 Slice<T>& operator*()
613 return *operator->();
621 iterator begin() {
return { *
this, _tdomain->begin() }; }
622 iterator end() {
return { *
this, _tdomain->end() }; }
624 struct reverse_iterator :
public base_container::reverse_iterator
628 reverse_iterator(SlicedTube& x, base_container::reverse_iterator it)
629 : base_container::reverse_iterator(it), _x(x) { }
631 std::shared_ptr<Slice<T>> operator->()
633 return _x.slice(*
this);
636 Slice<T>& operator*()
638 return *operator->();
646 reverse_iterator rbegin() {
return { *
this, _tdomain->rbegin() }; }
647 reverse_iterator rend() {
return { *
this, _tdomain->rend() }; }
649 struct const_iterator :
public base_container::const_iterator
653 const_iterator(
const SlicedTube& x, base_container::const_iterator it)
654 : base_container::const_iterator(it), _x(x) { }
656 std::shared_ptr<const Slice<T>> operator->()
658 return _x.slice(*
this);
661 const Slice<T>& operator*()
663 return *operator->();
668 const SlicedTube& _x;
671 const_iterator begin()
const {
return { *
this, _tdomain->cbegin() }; }
672 const_iterator end()
const {
return { *
this, _tdomain->cend() }; }
674 struct const_reverse_iterator :
public base_container::const_reverse_iterator
678 const_reverse_iterator(
const SlicedTube& x, base_container::const_reverse_iterator it)
679 : base_container::const_reverse_iterator(it), _x(x) { }
681 std::shared_ptr<const Slice<T>> operator->()
683 return _x.slice(*
this);
686 const Slice<T>& operator*()
688 return *operator->();
693 const SlicedTube& _x;
696 const_reverse_iterator rbegin()
const {
return { *
this, _tdomain->crbegin() }; }
697 const_reverse_iterator rend()
const {
return { *
this, _tdomain->crend() }; }
704 SlicedTube(
const std::shared_ptr<TDomain>& tdomain,
const AnalyticFunction<T>& f) ->
705 SlicedTube<typename Wrapper<T>::Domain>;
708 SlicedTube(
const std::shared_ptr<TDomain>& tdomain,
const SampledTraj<T>& f) ->
709 SlicedTube<typename Wrapper<T>::Domain>;
714 template<
typename... X>
715 inline void CtcBase<X...>::contract_tube(SlicedTube<X>&... x)
const
717 auto tdomain = std::get<0>(std::make_tuple(x...));
718 for(
auto it = tdomain.begin() ; it != tdomain.end() ; it++)
719 contract((x.slice(it)->codomain())...);
static Interval empty()
Provides an empty interval.
Definition codac2_Interval_impl.h:568
auto lb() const
Returns a matrix containing the lower bounds of each interval element.
Definition codac2_MatrixBase_addons_IntervalMatrixBase.h:91
bool is_superset(const Matrix< codac2::Interval, RowsAtCompileTime, ColsAtCompileTime > &x) const
Checks whether this matrix is a superset of another interval matrix.
Definition codac2_MatrixBase_addons_IntervalMatrixBase.h:832
bool is_subset(const Matrix< codac2::Interval, RowsAtCompileTime, ColsAtCompileTime > &x) const
Checks whether this matrix is a subset of another interval matrix.
Definition codac2_MatrixBase_addons_IntervalMatrixBase.h:652
auto rad() const
Returns a matrix containing the radii of each interval element.
Definition codac2_MatrixBase_addons_IntervalMatrixBase.h:197
Definition codac2_OctaSym.h:21
Eigen::Matrix< Interval,-1, 1 > IntervalVector
Alias for a dynamic-size column vector of intervals.
Definition codac2_IntervalVector.h:25