codac 2.0.0
Loading...
Searching...
No Matches
codac2_SlicedTube.h
Go to the documentation of this file.
1
9
10#pragma once
11
15#include "codac2_CtcDeriv.h"
16
17namespace codac2
18{
19 template<typename T>
20 class SlicedTube : public SlicedTubeBase
21 {
22 public:
23
24 explicit SlicedTube(const std::shared_ptr<TDomain>& tdomain,
25 const T& codomain)
26 : SlicedTubeBase(tdomain)
27 {
28 for(auto it = _tdomain->begin(); it != _tdomain->end(); ++it)
29 it->_slices.insert({
30 this,
31 std::make_shared<Slice<T>>(*this, it, codomain)
32 });
33 }
34
35 explicit SlicedTube(const std::shared_ptr<TDomain>& tdomain,
36 const AnalyticFunction<typename ExprType<T>::Type>& f)
37 : SlicedTubeBase(tdomain)
38 {
39 assert_release(f.args().size() == 1
40 && "function's inputs must be limited to one system variable");
41
42 for(auto it = _tdomain->begin(); it != _tdomain->end(); ++it)
43 it->_slices.insert({
44 this,
45 std::make_shared<Slice<T>>(*this, it, f.eval((Interval)*it))
46 });
47 }
48
49 template<typename V>
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)
54 {
55 for(auto it = _tdomain->begin(); it != _tdomain->end(); ++it)
56 it->_slices.insert({
57 this,
58 std::make_shared<Slice<T>>(*this, it, f((Interval)*it))
59 });
60 }
61
62 SlicedTube(const SlicedTube<T>& x)
63 : SlicedTubeBase(x.tdomain())
64 {
65 for(auto it = _tdomain->begin(); it != _tdomain->end(); ++it)
66 it->_slices.insert({
67 this,
68 std::make_shared<Slice<T>>(*x.slice(it), *this)
69 });
70 }
71
72 inline SlicedTube& operator=(const SlicedTube& x)
73 {
74 assert_release(_tdomain == x._tdomain);
75
76 for(auto it = _tdomain->begin(); it != _tdomain->end(); ++it)
77 slice(it)->set(x.slice(it)->codomain(), false);
78
79 return *this;
80 }
81
82 inline Index size() const
83 {
84 return first_slice()->size();
85 }
86
87 virtual std::pair<Index,Index> shape() const
88 {
89 if constexpr(std::is_same_v<typename ExprType<T>::Type,ScalarType>)
90 return {1,1};
91 else
92 return {first_slice()->codomain().rows(),first_slice()->codomain().cols()};
93 }
94
95 inline double volume() const
96 {
97 double volume = 0.;
98 for(const auto& s : *this)
99 if(!s.is_gate())
100 volume += s.volume();
101 return volume;
102 }
103
104 inline std::shared_ptr<Slice<T>> first_slice()
105 {
106 return std::const_pointer_cast<Slice<T>>(
107 static_cast<const SlicedTube&>(*this).first_slice());
108 }
109
110 inline std::shared_ptr<const Slice<T>> first_slice() const
111 {
112 return std::static_pointer_cast<const Slice<T>>(
113 this->SlicedTubeBase::first_slice());
114 }
115
116 inline std::shared_ptr<Slice<T>> last_slice()
117 {
118 return std::const_pointer_cast<Slice<T>>(
119 static_cast<const SlicedTube&>(*this).last_slice());
120 }
121
122 inline std::shared_ptr<const Slice<T>> last_slice() const
123 {
124 return std::static_pointer_cast<const Slice<T>>(
125 this->SlicedTubeBase::last_slice());
126 }
127
128 inline std::shared_ptr<Slice<T>> slice(const std::list<TSlice>::iterator& it)
129 {
130 return std::const_pointer_cast<Slice<T>>(
131 static_cast<const SlicedTube&>(*this).slice(it));
132 }
133
134 inline std::shared_ptr<const Slice<T>> slice(const std::list<TSlice>::const_iterator& it) const
135 {
136 return std::static_pointer_cast<const Slice<T>>(
137 it->slices().at(this));
138 }
139
140 inline std::shared_ptr<Slice<T>> slice(const std::list<TSlice>::reverse_iterator& it)
141 {
142 return std::const_pointer_cast<Slice<T>>(
143 static_cast<const SlicedTube&>(*this).slice(it));
144 }
145
146 inline std::shared_ptr<const Slice<T>> slice(const std::list<TSlice>::const_reverse_iterator& it) const
147 {
148 return std::static_pointer_cast<const Slice<T>>(
149 it->slices().at(this));
150 }
151
152 inline std::shared_ptr<Slice<T>> slice(std::shared_ptr<TSlice> ptr)
153 {
154 // Used in Python binding
155 auto it = std::find_if(_tdomain->begin(), _tdomain->end(),
156 [&](TSlice& t){ return &t == ptr.get(); });
157 return slice(it);
158 }
159
160 inline bool is_empty() const
161 {
162 // Fast evaluation by considering gates first, then envelopes,
163 // which allows to quickly identify an empty set
164 for(const auto& s : *this)
165 if(s.is_gate() && s.is_empty())
166 return true;
167 for(const auto& s : *this)
168 if(!s.is_gate() && s.is_empty())
169 return true;
170 return false;
171 }
172
173 inline bool is_unbounded() const
174 {
175 for(const auto& s : *this)
176 if(s.is_unbounded())
177 return true;
178 return false;
179 }
180
181 inline T codomain() const
182 {
183 T x = first_slice()->codomain();
184 for(const auto& s : *this)
185 x |= s.codomain();
186 return x;
187 }
188
189 template<typename Func>
190 inline T eval_common(const Interval& t, const Func& apply_eval) const
191 {
192 if(!tdomain()->t0_tf().is_superset(t))
193 return all_reals_value();
194
195 auto it = _tdomain->tslice(t.lb());
196 assert(it != _tdomain->end());
197 T codomain = apply_eval(it, t & *it);
198
199 while(it != std::next(_tdomain->tslice(t.ub())))
200 {
201 if(it->lb() == t.ub()) break;
202 codomain |= apply_eval(it, t & *it);
203 it++;
204 }
205
206 return codomain;
207 }
208
209 T operator()(const Interval& t) const
210 {
211 return eval_common(t,
212 [this](auto it, const Interval& t_) {
213 return slice(it)->operator()(t_);
214 });
215 }
216
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>)
219 {
220 return eval_common(t,
221 [this,&v](auto it, const Interval& t_) {
222 return slice(it)->operator()(t_, *v.slice(it));
223 });
224 }
225
226 std::pair<T,T> enclosed_bounds(const Interval& t) const
227 {
228 auto x = this->empty_value();
229 auto bounds = std::make_pair(x,x);
230
231 if(t.lb() < _tdomain->t0_tf().lb() || t.ub() > _tdomain->t0_tf().ub())
232 {
233 x.init(Interval(-oo,0));
234 bounds.first |= x;
235 x.init(Interval(0,oo));
236 bounds.second |= x;
237 }
238
239 Interval t_inter = t & _tdomain->t0_tf();
240 auto si = slice(_tdomain->tslice(t_inter.lb()));
241
242 while(si && si->t0_tf().lb() <= t_inter.ub())
243 {
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();
248 }
249
250 return bounds;
251 }
252
253 inline void set(const T& codomain)
254 {
255 assert_release(codomain.size() == this->size());
256
257 for(auto& s : *this)
258 if(!s.is_gate())
259 s.set(codomain, false);
260
261 for(auto& s : *this)
262 if(s.is_gate())
263 s.set(codomain, false);
264 }
265
266 inline void set(const T& codomain, double t)
267 {
268 assert_release(codomain.size() == this->size());
269 slice(_tdomain->sample(t,true))->set(codomain);
270 }
271
272 inline void set(const T& codomain, const Interval& t)
273 {
274 auto it_lb = _tdomain->sample(t.lb(), t.is_degenerated());
275 std::list<TSlice>::iterator it_ub;
276
277 if(!t.is_degenerated())
278 {
279 it_ub = _tdomain->sample(t.ub(), false);
280
281 if(it_ub->lb() == t.ub())
282 it_ub--; // pointing to the tslice [..,t.ub()]
283
284 if(it_lb->ub() == t.lb())
285 it_lb++;
286 }
287
288 else
289 it_ub = it_lb;
290
291 do
292 {
293 slice(it_lb)->set(codomain);
294 } while(it_lb != it_ub && (++it_lb) != _tdomain->end());
295 }
296
297 inline void set_ith_slice(const T& codomain, Index i)
298 {
299 Index j = 0;
300 for(auto& si : *this)
301 if(j++ == i)
302 {
303 si.set(codomain);
304 break;
305 }
306 }
307
308 template<typename V>
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)
311 {
312 // todo: faster implementation with iterators
313 for(auto& s : *this)
314 if(!s.is_gate())
315 s.set(T(s.codomain()).inflate(rad), false);
316
317 for(auto& s : *this)
318 if(s.is_gate())
319 s.set(T(s.codomain()).inflate(rad), false);
320
321 return *this;
322 }
323
324 template<typename V>
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)
327 {
328 // todo: faster implementation with iterators
329 assert_release(tdomain()->t0_tf().is_subset(rad.tdomain()));
330
331 for(auto& s : *this)
332 if(!s.is_gate())
333 s.set(T(s.codomain()).inflate(rad(s.t0_tf()).ub()), false);
334
335 for(auto& s : *this)
336 if(s.is_gate())
337 s.set(T(s.codomain()).inflate(rad(s.t0_tf()).ub()), false);
338
339 return *this;
340 }
341
342 SlicedTube<Interval> operator[](Index i) const
343 {
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];
348 return xi;
349 }
350
351 SlicedTube<IntervalVector> subvector(Index i, Index j) const
352 {
353 assert_release(i >= 0 && i <= j && j < size());
354 SlicedTube<IntervalVector> xij(tdomain(), IntervalVector(j-i+1));
355 for(auto it = tdomain()->begin() ; it != tdomain()->end() ; it++)
356 xij.slice(it)->codomain() = slice(it)->codomain().subvector(i,j);
357 return xij;
358 }
359
360 inline bool operator==(const SlicedTube& x) const
361 {
362 if(!TDomain::are_same(tdomain(), x.tdomain()))
363 return false;
364
365 auto it_this = _tdomain->begin();
366 auto it_x = x.tdomain()->cbegin();
367
368 while(it_this != _tdomain->end())
369 if(*slice(it_this++) != *x.slice(it_x++))
370 return false;
371
372 return true;
373 }
374
375 inline SlicedTube operator&=(const SlicedTube& x)
376 {
377 assert(TDomain::are_same(tdomain(), x.tdomain()));
378 auto it_this = _tdomain->begin();
379 auto it_x = x.tdomain()->cbegin();
380
381 while(it_this != _tdomain->end())
382 {
383 auto s = slice(it_this);
384 s->set(s->codomain() & x.slice(it_x)->codomain());
385 it_this++; it_x++;
386 }
387
388 assert(it_x == x.tdomain()->cend());
389 return *this;
390 }
391
392 friend inline std::ostream& operator<<(std::ostream& os, const SlicedTube<T>& x)
393 {
394 os << x.t0_tf()
395 << "↦" << x.codomain()
396 << ", " << x.nb_slices()
397 << " slice" << (x.nb_slices() > 1 ? "s" : "")
398 << std::flush;
399 return os;
400 }
401
402 AnalyticFunction<typename ExprType<T>::Type> as_function() const
403 {
404 ScalarVar t;
405 return {{t},
406 AnalyticExprWrapper<typename ExprType<T>::Type>(
407 std::make_shared<AnalyticOperationExpr<
408 TubeOp<SlicedTube<T>>,typename ExprType<T>::Type,ScalarType>>(*this,t))
409 };
410 }
411
412 template<typename Func>
413 Interval invert_common(const T& y, const Interval& t, const Func& apply_invert) const
414 {
415 assert_release(y.size() == this->size());
416
417 if(t.is_empty() || y.is_empty())
418 return Interval::empty();
419
420 else if(!t.is_subset(_tdomain->t0_tf()))
421 return Interval();
422
423 Interval invert = Interval::empty();
424 Interval t_ = t & _tdomain->t0_tf();
425
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++)
429 {
430 if(it->is_gate())
431 continue;
432
433 assert(it != _tdomain->end());
434 invert |= apply_invert(it, t & *it);
435 }
436
437 return invert;
438 }
439
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
443 {
444 assert_release(y.size() == this->size());
445 v_t.clear();
446
447 if(t.is_empty() || y.is_empty())
448 return;
449
450 else if(!t.is_subset(_tdomain->t0_tf()))
451 {
452 v_t.push_back(Interval());
453 return;
454 }
455
456 Interval invert = Interval::empty();
457 Interval t_ = t & _tdomain->t0_tf();
458
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++)
462 {
463 if(it->is_gate())
464 continue;
465
466 assert(it != _tdomain->end());
467 Interval local_invert = apply_invert(it, t & *it);
468
469 if(local_invert.is_empty() && !invert.is_empty())
470 {
471 v_t.push_back(invert);
472 invert.set_empty();
473 }
474
475 else
476 invert |= local_invert;
477 }
478
479 if(!invert.is_empty())
480 v_t.push_back(invert);
481 }
482
492 Interval invert(const T& y, const Interval& t = Interval()) const
493 {
494 return invert_common(y, t,
495 [this,&y](auto it, const Interval& t_) {
496 return slice(it)->invert(y, t_);
497 });
498 }
499
508 void invert(const T& y, std::vector<Interval> &v_t, const Interval& t = Interval()) const
509 {
510 return invert_common_subsets(y, v_t, t,
511 [this,&y](auto it, const Interval& t_) {
512 return slice(it)->invert(y, t_);
513 });
514 }
515
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>)
529 {
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_);
533 });
534 }
535
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>)
549 {
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_);
553 });
554 }
555
556 // Integral related methods
557
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;
562
563 inline SlicedTube<T> primitive() const
564 {
565 auto x0 = all_reals_value();
566 x0.init(0.);
567 return primitive(x0);
568 }
569
570 inline SlicedTube<T> primitive(const T& x0) const
571 {
572 auto x = all_reals_value();
573 auto p = SlicedTube<T>(this->tdomain(), x);
574 p.set(x0, this->tdomain()->t0_tf().lb()); // may create an unwanted gate
575 CtcDeriv c;
576 c.contract(p,*this);
577 return p;
578 }
579
580 inline T all_reals_value() const
581 {
582 T x = first_slice()->codomain();
583 x.init();
584 return x;
585 }
586
587 inline T empty_value() const
588 {
589 T x = first_slice()->codomain();
590 x.set_empty();
591 return x;
592 }
593
594
595 public:
596
597 using base_container = std::list<TSlice>;
598
599 struct iterator : public base_container::iterator
600 {
601 public:
602
603 iterator(SlicedTube& x, base_container::iterator it)
604 : base_container::iterator(it), _x(x) { }
605
606 std::shared_ptr<Slice<T>> operator->()
607 {
608 return _x.slice(*this);
609 }
610
611 Slice<T>& operator*()
612 {
613 return *operator->();
614 }
615
616 protected:
617
618 SlicedTube& _x;
619 };
620
621 iterator begin() { return { *this, _tdomain->begin() }; }
622 iterator end() { return { *this, _tdomain->end() }; }
623
624 struct reverse_iterator : public base_container::reverse_iterator
625 {
626 public:
627
628 reverse_iterator(SlicedTube& x, base_container::reverse_iterator it)
629 : base_container::reverse_iterator(it), _x(x) { }
630
631 std::shared_ptr<Slice<T>> operator->()
632 {
633 return _x.slice(*this);
634 }
635
636 Slice<T>& operator*()
637 {
638 return *operator->();
639 }
640
641 protected:
642
643 SlicedTube& _x;
644 };
645
646 reverse_iterator rbegin() { return { *this, _tdomain->rbegin() }; }
647 reverse_iterator rend() { return { *this, _tdomain->rend() }; }
648
649 struct const_iterator : public base_container::const_iterator
650 {
651 public:
652
653 const_iterator(const SlicedTube& x, base_container::const_iterator it)
654 : base_container::const_iterator(it), _x(x) { }
655
656 std::shared_ptr<const Slice<T>> operator->()
657 {
658 return _x.slice(*this);
659 }
660
661 const Slice<T>& operator*()
662 {
663 return *operator->();
664 }
665
666 protected:
667
668 const SlicedTube& _x;
669 };
670
671 const_iterator begin() const { return { *this, _tdomain->cbegin() }; }
672 const_iterator end() const { return { *this, _tdomain->cend() }; }
673
674 struct const_reverse_iterator : public base_container::const_reverse_iterator
675 {
676 public:
677
678 const_reverse_iterator(const SlicedTube& x, base_container::const_reverse_iterator it)
679 : base_container::const_reverse_iterator(it), _x(x) { }
680
681 std::shared_ptr<const Slice<T>> operator->()
682 {
683 return _x.slice(*this);
684 }
685
686 const Slice<T>& operator*()
687 {
688 return *operator->();
689 }
690
691 protected:
692
693 const SlicedTube& _x;
694 };
695
696 const_reverse_iterator rbegin() const { return { *this, _tdomain->crbegin() }; }
697 const_reverse_iterator rend() const { return { *this, _tdomain->crend() }; }
698 };
699
700
701 // Template deduction guide:
702
703 template<typename T>
704 SlicedTube(const std::shared_ptr<TDomain>& tdomain, const AnalyticFunction<T>& f) ->
705 SlicedTube<typename Wrapper<T>::Domain>;
706
707 template<typename T>
708 SlicedTube(const std::shared_ptr<TDomain>& tdomain, const SampledTraj<T>& f) ->
709 SlicedTube<typename Wrapper<T>::Domain>;
710
711
712 // Ctc
713
714 template<typename... X>
715 inline void CtcBase<X...>::contract_tube(SlicedTube<X>&... x) const
716 {
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())...);
720 }
721}
722
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