codac 2.0.0
Loading...
Searching...
No Matches
codac2::SlicedTube< T > Class Template Reference

Tube represented over a sliced temporal domain. More...

#include <codac2_SlicedTube.h>

Inheritance diagram for codac2::SlicedTube< T >:
Collaboration diagram for codac2::SlicedTube< T >:

Classes

struct  iterator
 Mutable iterator over the slices of this tube. More...
struct  reverse_iterator
 Mutable reverse iterator over the slices of this tube. More...
struct  const_iterator
 Constant iterator over the slices of this tube. More...
struct  const_reverse_iterator
 Constant reverse iterator over the slices of this tube. More...

Public Member Functions

 SlicedTube (const std::shared_ptr< TDomain > &tdomain, const T &codomain)
 Creates a sliced tube with constant codomain over all temporal slices.
 SlicedTube (const std::shared_ptr< TDomain > &tdomain, const AnalyticFunction< typename ExprType< T >::Type > &f)
 Creates a sliced tube by evaluating an analytic function on each temporal slice.
template<typename V>
 SlicedTube (const std::shared_ptr< TDomain > &tdomain, const AnalyticTraj< V > &x)
 Creates a sliced tube by evaluating an analytic trajectory on each temporal slice.
template<typename V>
 SlicedTube (const std::shared_ptr< TDomain > &tdomain, const SampledTraj< V > &x)
 Creates a sliced tube from a sampled trajectory.
 SlicedTube (const SlicedTube< T > &x)
 Copy constructor.
SlicedTubeoperator= (const SlicedTube &x)
 Assignment operator.
Index size () const
 Returns the codomain dimension.
virtual std::pair< Index, Index > shape () const
 Returns the matrix shape of this tube codomain.
double volume () const
 Returns the volume of this tube.
std::shared_ptr< Slice< T > > first_slice ()
 Returns the first slice of this tube.
std::shared_ptr< const Slice< T > > first_slice () const
 Returns the first slice of this tube.
std::shared_ptr< Slice< T > > last_slice ()
 Returns the last slice of this tube.
std::shared_ptr< const Slice< T > > last_slice () const
 Returns the last slice of this tube.
std::shared_ptr< Slice< T > > slice (const std::list< TSlice >::iterator &it)
 Returns the slice attached to a temporal iterator.
std::shared_ptr< const Slice< T > > slice (const std::list< TSlice >::const_iterator &it) const
 Returns the slice attached to a temporal iterator.
std::shared_ptr< Slice< T > > slice (const std::list< TSlice >::reverse_iterator &it)
 Returns the slice attached to a reverse temporal iterator.
std::shared_ptr< const Slice< T > > slice (const std::list< TSlice >::const_reverse_iterator &it) const
 Returns the slice attached to a reverse temporal iterator.
std::shared_ptr< Slice< T > > slice (std::shared_ptr< TSlice > ptr)
 Returns the slice attached to a temporal slice pointer.
bool is_empty () const
 Tests whether this tube is empty.
bool is_unbounded () const
 Tests whether this tube is unbounded.
codomain () const
 Returns the global codomain of this tube.
operator() (const Interval &t) const
 Evaluates this tube over a temporal interval.
operator() (const Interval &t, const SlicedTube< T > &v) const
 Returns the optimal evaluation over a temporal interval using a derivative tube.
std::pair< T, T > enclosed_bounds (const Interval &t) const
 Returns enclosed lower and upper bounds over a temporal interval.
void set (const T &codomain)
 Sets all codomains of this tube to the same value.
void set (const T &codomain, double t)
 Sets the codomain at one temporal instant.
void set (const T &codomain, const Interval &t)
 Sets the codomain over a temporal interval.
void set_ith_slice (const T &codomain, Index i)
 Sets the codomain of the \(i\)-th stored slice.
template<typename V>
const SlicedTube< T > & inflate (const V &rad)
 Inflates this tube by a constant radius.
template<typename V>
const SlicedTube< T > & inflate (const SampledTraj< V > &rad)
 Inflates this tube by a time-varying sampled radius.
SlicedTube< Intervaloperator[] (Index i) const
 Extracts one scalar component of this tube.
SlicedTube< IntervalVectorsubvector (Index i, Index j) const
 Extracts a subvector of this tube.
bool operator== (const SlicedTube &x) const
 Compares two sliced tubes.
SlicedTubeoperator&= (const SlicedTube &x)
 Intersects this tube with another one.
AnalyticFunction< typename ExprType< T >::Type > as_function () const
 Returns an analytic wrapper of this tube, allowing to evaluate this tube as an operator in an analytic expression.
Interval invert (const T &y) const
 Returns the interval inversion \(\left[[x]^{-1}([y])\right]\).
Interval invert (const T &y, const Interval &t) const
 Returns the interval inversion \(\left[[x]^{-1}([y])\right]\).
void invert (const T &y, std::vector< Interval > &v_t) const
 Computes the set of continuous values of the inversion \([x]^{-1}([y])\).
void invert (const T &y, std::vector< Interval > &v_t, const Interval &t) const
 Computes the set of continuous values of the inversion \([x]^{-1}([y])\).
Interval invert (const T &y, const SlicedTube< T > &v) const
 Returns the optimal interval inversion \(\left[[x]^{-1}([y])\right]\).
Interval invert (const T &y, const SlicedTube< T > &v, const Interval &t) const
 Returns the optimal interval inversion \(\left[[x]^{-1}([y])\right]\).
void invert (const T &y, std::vector< Interval > &v_t, const SlicedTube< T > &v) const
 Computes the set of continuous values of the optimal inversion \([x]^{-1}([y])\).
void invert (const T &y, std::vector< Interval > &v_t, const SlicedTube< T > &v, const Interval &t) const
 Computes the set of continuous values of the optimal inversion \([x]^{-1}([y])\).
integral (const Interval &t) const
 Returns an enclosure of the integrals of this tube from \(t_0\) to \([t]\).
integral (const Interval &t1, const Interval &t2) const
 Returns an enclosure of the integrals of this tube between the time intervals \([t_1]\) and \([t_2]\).
std::pair< T, T > partial_integral (const Interval &t) const
 Returns lower and upper enclosures of the integrals of this tube \([x](\cdot)=[x^-(\cdot),x^+(\cdot)]\) from \(t_0\) to \([t]\).
std::pair< T, T > partial_integral (const Interval &t1, const Interval &t2) const
 Returns lower and upper enclosures of the integrals of this tube \([x](\cdot)=[x^-(\cdot),x^+(\cdot)]\) between \([t_1]\) and \([t_2]\).
SlicedTube< T > primitive () const
 Returns a primitive of this tube with zero initial condition.
SlicedTube< T > primitive (const T &x0) const
 Returns a primitive of this tube with prescribed initial condition.
all_reals_value () const
 Returns the unbounded codomain value associated with T.
empty_value () const
 Returns the empty codomain value associated with T.
SampledTraj< typename Scalar< T >::Type > mid () const
 Returns the approximated midpoint trajectory of this tube.
iterator begin ()
 Returns an iterator to the first temporal slice.
iterator end ()
 Returns an iterator past the last temporal slice.
reverse_iterator rbegin ()
 Returns a reverse iterator to the last temporal slice.
reverse_iterator rend ()
 Returns a reverse iterator past the first temporal slice.
const_iterator begin () const
 Returns a constant iterator to the first temporal slice.
const_iterator end () const
 Returns a constant iterator past the last temporal slice.
const_reverse_iterator rbegin () const
 Returns a constant reverse iterator to the last temporal slice.
const_reverse_iterator rend () const
 Returns a constant reverse iterator past the first temporal slice.
Public Member Functions inherited from codac2::SlicedTubeBase
 SlicedTubeBase (const std::shared_ptr< TDomain > &tdomain)
 Creates a sliced tube over a given temporal domain.
 ~SlicedTubeBase ()
 Destroys this sliced tube.
size_t nb_slices () const
 Returns the number of temporal elements of this tube.
std::shared_ptr< const SliceBasefirst_slice () const
 Returns the first slice of this tube.
std::shared_ptr< const SliceBaselast_slice () const
 Returns the last slice of this tube.
Public Member Functions inherited from codac2::TubeBase
 TubeBase (const std::shared_ptr< TDomain > &tdomain)
 Creates a tube over a given temporal domain.
const std::shared_ptr< TDomain > & tdomain () const
 Returns the temporal domain of this tube.
Interval t0_tf () const
 Returns the global temporal interval of this tube.

Additional Inherited Members

Protected Attributes inherited from codac2::TubeBase
const std::shared_ptr< TDomain_tdomain
 Shared temporal domain of this tube.

Detailed Description

template<typename T>
class codac2::SlicedTube< T >

Tube represented over a sliced temporal domain.

A SlicedTube<T> is a tube whose temporal domain is represented by a shared TDomain made of TSlice objects.

For each temporal slice of this partition, the tube stores a Slice<T> describing the codomain of the tube over that temporal support.

The codomain type T is typically Interval or IntervalVector, or any custom domain implemented by the user.

Template Parameters
Tcodomain type of the tube

Constructor & Destructor Documentation

◆ SlicedTube() [1/5]

template<typename T>
codac2::SlicedTube< T >::SlicedTube ( const std::shared_ptr< TDomain > & tdomain,
const T & codomain )
inlineexplicit

Creates a sliced tube with constant codomain over all temporal slices.

Parameters
tdomainshared temporal domain of the tube
codomaincodomain assigned to each slice
50 {
51 for(auto it = _tdomain->begin(); it != _tdomain->end(); ++it)
52 it->_slices.insert({
53 this,
55 });
56 }
SlicedTubeBase(const std::shared_ptr< TDomain > &tdomain)
Creates a sliced tube over a given temporal domain.
Definition codac2_SlicedTubeBase.h:43
Tube represented over a sliced temporal domain.
Definition codac2_SlicedTube.h:38
T codomain() const
Returns the global codomain of this tube.
Definition codac2_SlicedTube.h:347
const std::shared_ptr< TDomain > _tdomain
Shared temporal domain of this tube.
Definition codac2_TubeBase.h:71
const std::shared_ptr< TDomain > & tdomain() const
Returns the temporal domain of this tube.
Definition codac2_TubeBase.h:49

◆ SlicedTube() [2/5]

template<typename T>
codac2::SlicedTube< T >::SlicedTube ( const std::shared_ptr< TDomain > & tdomain,
const AnalyticFunction< typename ExprType< T >::Type > & f )
inlineexplicit

Creates a sliced tube by evaluating an analytic function on each temporal slice.

The function is evaluated on each temporal interval of the associated TDomain.

Parameters
tdomainshared temporal domain of the tube
fanalytic function of one scalar variable
70 {
71 assert_release(f.args().size() == 1
72 && "function's inputs must be limited to one system variable");
73
74 for(auto it = _tdomain->begin(); it != _tdomain->end(); ++it)
75 it->_slices.insert({
76 this,
77 std::make_shared<Slice<T>>(*this, it, f.eval((Interval)*it))
78 });
79 }
Index size() const
Returns the codomain dimension.
Definition codac2_SlicedTube.h:163

◆ SlicedTube() [3/5]

template<typename T>
template<typename V>
codac2::SlicedTube< T >::SlicedTube ( const std::shared_ptr< TDomain > & tdomain,
const AnalyticTraj< V > & x )
inlineexplicit

Creates a sliced tube by evaluating an analytic trajectory on each temporal slice.

The function is evaluated on each temporal interval of the associated TDomain.

Parameters
tdomainshared temporal domain of the tube
xanalytic trajectory
95 {
96 for(auto it = _tdomain->begin(); it != _tdomain->end(); ++it)
97 it->_slices.insert({
98 this,
99 std::make_shared<Slice<T>>(*this, it, x.eval((Interval)*it))
100 });
101 }

◆ SlicedTube() [4/5]

template<typename T>
template<typename V>
codac2::SlicedTube< T >::SlicedTube ( const std::shared_ptr< TDomain > & tdomain,
const SampledTraj< V > & x )
inlineexplicit

Creates a sliced tube from a sampled trajectory.

Template Parameters
Vsampled value type
Parameters
tdomainshared temporal domain of the tube
xsampled trajectory evaluated on each temporal slice
115 {
116 for(auto it = _tdomain->begin(); it != _tdomain->end(); ++it)
117 it->_slices.insert({
118 this,
120 });
121 }

◆ SlicedTube() [5/5]

template<typename T>
codac2::SlicedTube< T >::SlicedTube ( const SlicedTube< T > & x)
inline

Copy constructor.

The copied tube shares the same TDomain as the source tube.

Parameters
xtube to be copied
132 {
133 for(auto it = _tdomain->begin(); it != _tdomain->end(); ++it)
134 it->_slices.insert({
135 this,
137 });
138 }
std::shared_ptr< Slice< T > > slice(const std::list< TSlice >::iterator &it)
Returns the slice attached to a temporal iterator.
Definition codac2_SlicedTube.h:249

Member Function Documentation

◆ operator=()

template<typename T>
SlicedTube & codac2::SlicedTube< T >::operator= ( const SlicedTube< T > & x)
inline

Assignment operator.

Both tubes are expected to share the same temporal domain object.

Parameters
xsource tube
Returns
reference to this tube
149 {
151
152 for(auto it = _tdomain->begin(); it != _tdomain->end(); ++it)
153 slice(it)->set(x.slice(it)->codomain(), false);
154
155 return *this;
156 }

◆ size()

template<typename T>
Index codac2::SlicedTube< T >::size ( ) const
inline

Returns the codomain dimension.

Returns
dimension of the codomain type T
164 {
165 return first_slice()->size();
166 }
std::shared_ptr< Slice< T > > first_slice()
Returns the first slice of this tube.
Definition codac2_SlicedTube.h:204

◆ shape()

template<typename T>
virtual std::pair< Index, Index > codac2::SlicedTube< T >::shape ( ) const
inlinevirtual

Returns the matrix shape of this tube codomain.

For scalar codomains, the returned shape is {1,1}.

Returns
pair (rows, cols)
176 {
178 return {1,1};
179 else
180 return {first_slice()->codomain().rows(),first_slice()->codomain().cols()};
181 }

◆ volume()

template<typename T>
double codac2::SlicedTube< T >::volume ( ) const
inline

Returns the volume of this tube.

Only non-gate slices contribute to the volume.

Returns
sum of the volumes of all non-gate slices
191 {
192 double volume = 0.;
193 for(const auto& s : *this)
194 if(!s.is_gate())
195 volume += s.volume();
196 return volume;
197 }
double volume() const
Returns the volume of this tube.
Definition codac2_SlicedTube.h:190

◆ first_slice() [1/2]

template<typename T>
std::shared_ptr< Slice< T > > codac2::SlicedTube< T >::first_slice ( )
inline

Returns the first slice of this tube.

Returns
shared pointer to the first slice
205 {
207 static_cast<const SlicedTube&>(*this).first_slice());
208 }
SlicedTube(const std::shared_ptr< TDomain > &tdomain, const T &codomain)
Creates a sliced tube with constant codomain over all temporal slices.
Definition codac2_SlicedTube.h:47

◆ first_slice() [2/2]

template<typename T>
std::shared_ptr< const Slice< T > > codac2::SlicedTube< T >::first_slice ( ) const
inline

Returns the first slice of this tube.

Returns
shared pointer to the first slice
216 {
219 }
std::shared_ptr< const SliceBase > first_slice() const
Returns the first slice of this tube.
Definition codac2_SlicedTubeBase.h:81

◆ last_slice() [1/2]

template<typename T>
std::shared_ptr< Slice< T > > codac2::SlicedTube< T >::last_slice ( )
inline

Returns the last slice of this tube.

Returns
shared pointer to the last slice
227 {
229 static_cast<const SlicedTube&>(*this).last_slice());
230 }

◆ last_slice() [2/2]

template<typename T>
std::shared_ptr< const Slice< T > > codac2::SlicedTube< T >::last_slice ( ) const
inline

Returns the last slice of this tube.

Returns
shared pointer to the last slice
238 {
241 }
std::shared_ptr< const SliceBase > last_slice() const
Returns the last slice of this tube.
Definition codac2_SlicedTubeBase.h:97

◆ slice() [1/5]

template<typename T>
std::shared_ptr< Slice< T > > codac2::SlicedTube< T >::slice ( const std::list< TSlice >::iterator & it)
inline

Returns the slice attached to a temporal iterator.

Parameters
ititerator to a TSlice of the shared TDomain
Returns
shared pointer to the corresponding slice
250 {
252 static_cast<const SlicedTube&>(*this).slice(it));
253 }

◆ slice() [2/5]

template<typename T>
std::shared_ptr< const Slice< T > > codac2::SlicedTube< T >::slice ( const std::list< TSlice >::const_iterator & it) const
inline

Returns the slice attached to a temporal iterator.

Parameters
itconstant iterator to a TSlice of the shared TDomain
Returns
shared pointer to the corresponding slice
262 {
264 it->slices().at(this));
265 }

◆ slice() [3/5]

template<typename T>
std::shared_ptr< Slice< T > > codac2::SlicedTube< T >::slice ( const std::list< TSlice >::reverse_iterator & it)
inline

Returns the slice attached to a reverse temporal iterator.

Parameters
itreverse iterator to a TSlice of the shared TDomain
Returns
shared pointer to the corresponding slice
274 {
276 static_cast<const SlicedTube&>(*this).slice(it));
277 }

◆ slice() [4/5]

template<typename T>
std::shared_ptr< const Slice< T > > codac2::SlicedTube< T >::slice ( const std::list< TSlice >::const_reverse_iterator & it) const
inline

Returns the slice attached to a reverse temporal iterator.

Parameters
itconstant reverse iterator to a TSlice of the shared TDomain
Returns
shared pointer to the corresponding slice
286 {
288 it->slices().at(this));
289 }

◆ slice() [5/5]

template<typename T>
std::shared_ptr< Slice< T > > codac2::SlicedTube< T >::slice ( std::shared_ptr< TSlice > ptr)
inline

Returns the slice attached to a temporal slice pointer.

This overload is mainly intended for language bindings.

Parameters
ptrpointer to a temporal slice of the shared TDomain
Returns
shared pointer to the corresponding slice
300 {
301 // Used in Python binding
302 auto it = std::find_if(_tdomain->begin(), _tdomain->end(),
303 [&](TSlice& t){ return &t == ptr.get(); });
304 assert(it != _tdomain->end());
305 return slice(it);
306 }

◆ is_empty()

template<typename T>
bool codac2::SlicedTube< T >::is_empty ( ) const
inline

Tests whether this tube is empty.

A fast evaluation is done by considering gates first, then envelopes, which allows to quickly identify an empty set

Returns
true if at least one slice is empty, false otherwise
317 {
318 for(const auto& s : *this)
319 if(s.is_gate() && s.is_empty())
320 return true;
321 for(const auto& s : *this)
322 if(!s.is_gate() && s.is_empty())
323 return true;
324 return false;
325 }
bool is_empty() const
Tests whether this tube is empty.
Definition codac2_SlicedTube.h:316

◆ is_unbounded()

template<typename T>
bool codac2::SlicedTube< T >::is_unbounded ( ) const
inline

Tests whether this tube is unbounded.

Returns
true if at least one slice is unbounded, false otherwise
333 {
334 for(const auto& s : *this)
335 if(s.is_unbounded())
336 return true;
337 return false;
338 }
bool is_unbounded() const
Tests whether this tube is unbounded.
Definition codac2_SlicedTube.h:332

◆ codomain()

template<typename T>
T codac2::SlicedTube< T >::codomain ( ) const
inline

Returns the global codomain of this tube.

The returned codomain is the union of the codomains of all slices.

Returns
global codomain enclosure
348 {
349 T x = first_slice()->codomain();
350 for(const auto& s : *this)
351 x |= s.codomain();
352 return x;
353 }

◆ operator()() [1/2]

template<typename T>
T codac2::SlicedTube< T >::operator() ( const Interval & t) const
inline

Evaluates this tube over a temporal interval.

Parameters
ttemporal interval
Returns
enclosure of \([x]([t])\)
382 {
383 return eval_common(t,
384 [this](auto it, const Interval& t_) {
385 return slice(it)->operator()(t_);
386 });
387 }

◆ operator()() [2/2]

template<typename T>
T codac2::SlicedTube< T >::operator() ( const Interval & t,
const SlicedTube< T > & v ) const
inline

Returns the optimal evaluation over a temporal interval using a derivative tube.

Parameters
ttemporal interval
vderivative tube such that \(\dot{x}(\cdot)\in[v](\cdot)\)
Returns
enclosure of \([x]([t])\)
398 {
399 return eval_common(t,
400 [this,&v](auto it, const Interval& t_) {
401 return slice(it)->operator()(t_, *v.slice(it));
402 });
403 }

◆ enclosed_bounds()

template<typename T>
std::pair< T, T > codac2::SlicedTube< T >::enclosed_bounds ( const Interval & t) const
inline

Returns enclosed lower and upper bounds over a temporal interval.

Parameters
ttemporal interval
Returns
pair made of enclosed lower and upper bounds
412 {
413 auto x = this->empty_value();
414 auto bounds = std::make_pair(x,x);
415
416 if(t.lb() < _tdomain->t0_tf().lb() || t.ub() > _tdomain->t0_tf().ub())
417 {
418 x.init(Interval(-oo,0));
419 bounds.first |= x;
420 x.init(Interval(0,oo));
421 bounds.second |= x;
422 }
423
424 Interval t_inter = t & _tdomain->t0_tf();
425 auto si = slice(_tdomain->tslice(t_inter.lb()));
426
427 while(si && si->t0_tf().lb() <= t_inter.ub())
428 {
430 bounds.first |= slice_bounds.first;
431 bounds.second |= slice_bounds.second;
432 si = si->next_slice();
433 }
434
435 return bounds;
436 }
double lb() const
Returns the lower bound of this.
Definition codac2_Interval_impl.h:110
std::pair< T, T > enclosed_bounds(const Interval &t) const
Returns enclosed lower and upper bounds over a temporal interval.
Definition codac2_SlicedTube.h:411
T empty_value() const
Returns the empty codomain value associated with T.
Definition codac2_SlicedTube.h:1063
Interval t0_tf() const
Returns the global temporal interval of this tube.
Definition codac2_TubeBase.h:61

◆ set() [1/3]

template<typename T>
void codac2::SlicedTube< T >::set ( const T & codomain)
inline

Sets all codomains of this tube to the same value.

Parameters
codomainnew codomain
444 {
445 assert_release(codomain.size() == this->size());
446
447 for(auto& s : *this)
448 if(!s.is_gate())
449 s.set(codomain, false);
450
451 for(auto& s : *this)
452 if(s.is_gate())
453 s.set(codomain, false);
454 }
void set(const T &codomain)
Sets all codomains of this tube to the same value.
Definition codac2_SlicedTube.h:443

◆ set() [2/3]

template<typename T>
void codac2::SlicedTube< T >::set ( const T & codomain,
double t )
inline

Sets the codomain at one temporal instant.

The temporal domain may be refined so that t becomes an explicit gate.

Parameters
codomainnew codomain
ttemporal instant
465 {
466 assert_release(codomain.size() == this->size());
467 slice(_tdomain->sample(t,true))->set(codomain);
468 }

◆ set() [3/3]

template<typename T>
void codac2::SlicedTube< T >::set ( const T & codomain,
const Interval & t )
inline

Sets the codomain over a temporal interval.

The temporal domain may be refined so that the bounds of t become explicit gates.

Parameters
codomainnew codomain
ttemporal interval
480 {
481 auto it_lb = _tdomain->sample(t.lb(), t.is_degenerated());
483
484 if(!t.is_degenerated())
485 {
486 it_ub = _tdomain->sample(t.ub(), false);
487
488 if(it_ub->lb() == t.ub())
489 it_ub--; // pointing to the tslice [..,t.ub()]
490
491 if(it_lb->ub() == t.lb())
492 it_lb++;
493 }
494
495 else
496 it_ub = it_lb;
497
498 do
499 {
500 slice(it_lb)->set(codomain);
501 } while(it_lb != it_ub && (++it_lb) != _tdomain->end());
502 }

◆ set_ith_slice()

template<typename T>
void codac2::SlicedTube< T >::set_ith_slice ( const T & codomain,
Index i )
inline

Sets the codomain of the \(i\)-th stored slice.

Parameters
codomainnew codomain
islice index
511 {
512 Index j = 0;
513 for(auto& si : *this)
514 if(j++ == i)
515 {
516 si.set(codomain);
517 break;
518 }
519 }

◆ inflate() [1/2]

template<typename T>
template<typename V>
const SlicedTube< T > & codac2::SlicedTube< T >::inflate ( const V & rad)
inline

Inflates this tube by a constant radius.

Template Parameters
Vradius type
Parameters
radinflation radius
Returns
reference to this tube
531 {
532 // todo: faster implementation with iterators
533 for(auto& s : *this)
534 if(!s.is_gate())
535 s.set(T(s.codomain()).inflate(rad), false);
536
537 for(auto& s : *this)
538 if(s.is_gate())
539 s.set(T(s.codomain()).inflate(rad), false);
540
541 return *this;
542 }
const SlicedTube< T > & inflate(const V &rad)
Inflates this tube by a constant radius.
Definition codac2_SlicedTube.h:530

◆ inflate() [2/2]

template<typename T>
template<typename V>
const SlicedTube< T > & codac2::SlicedTube< T >::inflate ( const SampledTraj< V > & rad)
inline

Inflates this tube by a time-varying sampled radius.

Template Parameters
Vsampled radius type
Parameters
radsampled inflation radius
Returns
reference to this tube
554 {
555 // todo: faster implementation with iterators
557
558 for(auto& s : *this)
559 if(!s.is_gate())
560 s.set(T(s.codomain()).inflate(rad(s.t0_tf()).ub()), false);
561
562 for(auto& s : *this)
563 if(s.is_gate())
564 s.set(T(s.codomain()).inflate(rad(s.t0_tf()).ub()), false);
565
566 return *this;
567 }
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:678

◆ operator[]()

template<typename T>
SlicedTube< Interval > codac2::SlicedTube< T >::operator[] ( Index i) const
inline

Extracts one scalar component of this tube.

Parameters
icomponent index
Returns
scalar sliced tube corresponding to the \(i\)-th component
576 {
577 assert_release(i >= 0 && i < size());
579 for(auto it = tdomain()->begin() ; it != tdomain()->end() ; it++)
580 xi.slice(it)->codomain() = slice(it)->codomain()[i];
581 return xi;
582 }
iterator begin()
Returns an iterator to the first temporal slice.
Definition codac2_SlicedTube.h:1120
iterator end()
Returns an iterator past the last temporal slice.
Definition codac2_SlicedTube.h:1125

◆ subvector()

template<typename T>
SlicedTube< IntervalVector > codac2::SlicedTube< T >::subvector ( Index i,
Index j ) const
inline

Extracts a subvector of this tube.

Parameters
ifirst component index
jlast component index
Returns
sliced tube associated with the subvector \([i,j]\)
592 {
593 assert_release(i >= 0 && i <= j && j < size());
595 for(auto it = tdomain()->begin() ; it != tdomain()->end() ; it++)
596 xij.slice(it)->codomain() = slice(it)->codomain().subvector(i,j);
597 return xij;
598 }
Eigen::Matrix< Interval,-1, 1 > IntervalVector
Alias for a dynamic-size column vector of intervals.
Definition codac2_IntervalVector.h:25

◆ operator==()

template<typename T>
bool codac2::SlicedTube< T >::operator== ( const SlicedTube< T > & x) const
inline

Compares two sliced tubes.

Two sliced tubes are equal if they share the same temporal partition and if all corresponding codomains are equal.

Parameters
xtube to compare with
Returns
true if both tubes are equal, false otherwise
610 {
612 return false;
613
614 auto it_this = _tdomain->begin();
615 auto it_x = x.tdomain()->cbegin();
616
617 while(it_this != _tdomain->end())
618 if(*slice(it_this++) != *x.slice(it_x++))
619 return false;
620
621 return true;
622 }

◆ operator&=()

template<typename T>
SlicedTube & codac2::SlicedTube< T >::operator&= ( const SlicedTube< T > & x)
inline

Intersects this tube with another one.

Both tubes are expected to be defined over the same temporal partition.

Parameters
xtube to intersect with
Returns
updated tube
633 {
635 auto it_this = _tdomain->begin();
636 auto it_x = x.tdomain()->cbegin();
637
638 while(it_this != _tdomain->end())
639 {
640 auto s = slice(it_this);
641 s->set(s->codomain() & x.slice(it_x)->codomain());
642 it_this++; it_x++;
643 }
644
645 assert(it_x == x.tdomain()->cend());
646 return *this;
647 }

◆ as_function()

template<typename T>
AnalyticFunction< typename ExprType< T >::Type > codac2::SlicedTube< T >::as_function ( ) const
inline

Returns an analytic wrapper of this tube, allowing to evaluate this tube as an operator in an analytic expression.

Returns
analytic function associated with this tube

◆ invert() [1/4]

template<typename T>
Interval codac2::SlicedTube< T >::invert ( const T & y) const
inline

Returns the interval inversion \(\left[[x]^{-1}([y])\right]\).

If the inversion results in several pre-images, their union is returned.

Parameters
yinterval codomain
Returns
the hull of \([x]^{-1}([y])\)
768 {
769 return invert(y,t0_tf());
770 }
Interval invert(const T &y) const
Returns the interval inversion .
Definition codac2_SlicedTube.h:767

◆ invert() [2/4]

template<typename T>
void codac2::SlicedTube< T >::invert ( const T & y,
std::vector< Interval > & v_t ) const
inline

Computes the set of continuous values of the inversion \([x]^{-1}([y])\).

Parameters
yinterval codomain
v_tvector of the sub-tdomains \([t_k]\) for which \(\forall t\in[t_k] \mid x(t)\in[y], x(\cdot)\in[x](\cdot)\)
798 {
799 invert(y,v_t,t0_tf());
800 }

◆ invert() [3/4]

template<typename T>
Interval codac2::SlicedTube< T >::invert ( const T & y,
const SlicedTube< T > & v ) const
inline

Returns the optimal interval inversion \(\left[[x]^{-1}([y])\right]\).

The knowledge of the derivative tube \([v](\cdot)\) allows a finer inversion. If the inversion results in several pre-images, their union is returned.

Parameters
yinterval codomain
vderivative tube such that \(\dot{x}(\cdot)\in[v](\cdot)\)
Returns
hull of \([x]^{-1}([y])\)
831 {
832 return invert(y,v,t0_tf());
833 }

◆ invert() [4/4]

template<typename T>
void codac2::SlicedTube< T >::invert ( const T & y,
std::vector< Interval > & v_t,
const SlicedTube< T > & v ) const
inline

Computes the set of continuous values of the optimal inversion \([x]^{-1}([y])\).

The knowledge of the derivative tube \([v](\cdot)\) allows finer inversions.

Parameters
yinterval codomain
v_tvector of the sub-tdomains \([t_k]\) for which \(\forall t\in[t_k] \mid x(t)\in[y], x(\cdot)\in[x](\cdot), \dot{x}(\cdot)\in[v](\cdot)\)
vderivative tube such that \(\dot{x}(\cdot)\in[v](\cdot)\)
868 {
869 invert(y,v_t,v,t0_tf());
870 }

◆ all_reals_value()

template<typename T>
T codac2::SlicedTube< T >::all_reals_value ( ) const
inline

Returns the unbounded codomain value associated with T.

Returns
unbounded value of type T
1052 {
1053 T x = first_slice()->codomain();
1054 x.init();
1055 return x;
1056 }

◆ empty_value()

template<typename T>
T codac2::SlicedTube< T >::empty_value ( ) const
inline

Returns the empty codomain value associated with T.

Returns
empty value of type T
1064 {
1065 T x = first_slice()->codomain();
1066 x.set_empty();
1067 return x;
1068 }

◆ mid()

template<typename T>
SampledTraj< typename Scalar< T >::Type > codac2::SlicedTube< T >::mid ( ) const
inline

Returns the approximated midpoint trajectory of this tube.

Template Parameters
Vsampled value type
Returns
midpoint sampled trajectory
1077 {
1079 double t0 = _tdomain->t0_tf().lb();
1080 m.set((*this)(t0).mid(), t0);
1081 for(const auto& s : *this)
1082 if(!s.is_gate())
1083 m.set(s.output_gate().mid(),s.t0_tf().ub());
1084 return m;
1085 }
double ub() const
Returns the upper bound of this.
Definition codac2_Interval_impl.h:115
SampledTraj< typename Scalar< T >::Type > mid() const
Returns the approximated midpoint trajectory of this tube.
Definition codac2_SlicedTube.h:1076

The documentation for this class was generated from the following files: