codac 1.5.6
Loading...
Searching...
No Matches
codac::Tube Class Reference

One dimensional tube \([x](\cdot)\), defined as an interval of scalar trajectories. More...

#include <codac_Tube.h>

Inheritance diagram for codac::Tube:
Collaboration diagram for codac::Tube:

Public Member Functions

Definition
 Tube (const Interval &tdomain, const Interval &codomain=Interval::ALL_REALS)
 Creates a scalar tube \([x](\cdot)\) made of one slice.
 
 Tube (const Interval &tdomain, double timestep, const Interval &codomain=Interval::ALL_REALS)
 Creates a scalar tube \([x](\cdot)\) with some temporal discretization.
 
 Tube (const Interval &tdomain, double timestep, const TFnc &f, int f_image_id=0)
 Creates a scalar tube \([x](\cdot)\) from a TFnc object and with some temporal discretization.
 
 Tube (const std::vector< Interval > &v_tdomains, const std::vector< Interval > &v_codomains)
 Creates a tube \([x](\cdot)\) from a list of \(k\) boxes \(\big([t_1]\times[x_1],\dots,[t_k]\times[x_k]\big)\).
 
 Tube (const Tube &x)
 Creates a copy of a scalar tube \([x](\cdot)\), with the same time discretization.
 
 Tube (const Tube &x, const TFnc &f, int f_image_id=0)
 Creates a copy of a scalar tube \([x](\cdot)\), with the same time discretization but a specific codomain defined by a TFnc object.
 
 Tube (const Trajectory &traj, double timestep)
 Creates a scalar tube \([x](\cdot)\) enclosing a trajectory \(x(\cdot)\), possibly with some temporal discretization.
 
 Tube (const Trajectory &lb, const Trajectory &ub, double timestep)
 Creates a scalar tube \([x](\cdot)\) defined as an interval of two trajectories \([lb(\cdot),ub(\cdot)]\).
 
 Tube (const std::string &binary_file_name)
 Restore a scalar tube from serialization.
 
 Tube (const std::string &binary_file_name, Trajectory *&traj)
 Restore a scalar tube from serialization, together with a Trajectory object.
 
 ~Tube ()
 Tube destructor.
 
int size () const
 Returns the dimension of the scalar tube (always 1)
 
const Tube primitive (const Interval &c=Interval(0.)) const
 Returns the primitive Tube of this tube.
 
const Tubeoperator= (const Tube &x)
 Returns a copy of a Tube.
 
const Interval tdomain () const
 Returns the temporal definition domain of this tube.
 
const Polygon polygon_envelope () const
 Returns the polygon envelope of this tube.
 
const Trajectory lb () const
 Returns a possible lower bound \(x^{-}(\cdot)\) of the tube.
 
const Trajectory ub () const
 Returns a possible upper bound \(x^{+}(\cdot)\) of the tube.
 
Slices structure
int nb_slices () const
 Returns the number of slices of this tube.
 
Sliceslice (int slice_id)
 Returns a pointer to the ith Slice object of this tube.
 
const Sliceslice (int slice_id) const
 Returns a constant pointer to the ith Slice object of this tube.
 
Sliceslice (double t)
 Returns a pointer to the Slice object of this tube that is defined at \(t\).
 
const Sliceslice (double t) const
 Returns a constant pointer to the Slice object of this tube that is defined at \(t\).
 
Slicefirst_slice ()
 Returns a pointer to the first Slice object of this tube.
 
const Slicefirst_slice () const
 Returns a constant pointer to the first Slice object of this tube.
 
Slicelast_slice ()
 Returns a pointer to the last Slice object of this tube.
 
const Slicelast_slice () const
 Returns a constant pointer to the last Slice object of this tube.
 
Slicewider_slice ()
 Returns a pointer to the Slice defined over the wider temporal domain.
 
const Slicewider_slice () const
 Returns a constant pointer to the Slice defined over the wider temporal domain.
 
Slicelargest_slice ()
 Returns a pointer to the Slice object of this tube for which the interval value is the most uncertain.
 
const Slicelargest_slice () const
 Returns a constant pointer to the Slice object of this tube for which the interval value is the most uncertain.
 
Slicesteepest_slice ()
 Returns a pointer to the Slice object of this tube for which the difference between the mid values of the input and output gates is the largest.
 
const Slicesteepest_slice () const
 Returns a const pointer to the Slice object of this tube for which the difference between the mid values of the input and output gates is the largest.
 
const Interval slice_tdomain (int slice_id) const
 Returns the temporal definition domain of the ith Slice of this tube.
 
int time_to_index (double t) const
 Returns the Slice index related to the temporal key \(t\).
 
int index (const Slice *slice) const
 Returns the Slice index related to the Slice pointer.
 
void sample (double t)
 Samples this tube at \(t\).
 
void sample (double t, Slice *slice_to_be_sampled)
 Samples this tube at \(t\) from a pointer to the corresponding slice.
 
void sample (double t, const Interval &gate)
 Samples this tube at \(t\) with a specific gate value.
 
void sample (const Tube &x)
 Samples this tube so that it will share the same sampling of \([x](\cdot)\).
 
bool gate_exists (double t) const
 Tests if a gate exists at time \(t\).
 
void remove_gate (double t)
 Removes the gate at \(t\) and merges the two related slices.
 
void merge_similar_slices (double distance_threshold)
 Merges all adjacent slices whose Hausdorff distance is less than the given threshold.
 
Accessing values
const Interval codomain () const
 Returns the interval of feasible values.
 
double volume () const
 Returns the volume of this tube.
 
const Interval operator() (int slice_id) const
 Returns the value of the ith slice.
 
const Interval operator() (double t) const
 Returns the evaluation of this tube at \(t\).
 
const Interval operator() (const Interval &t) const
 Returns the interval evaluation of this tube over \([t]\).
 
const std::pair< Interval, Interval > eval (const Interval &t=Interval::ALL_REALS) const
 Returns the interval evaluations of the bounds of the tube \(\underline{x^-}(\cdot)\) and \(\overline{x^+}(\cdot)\) over \([t]\).
 
const Interval interpol (double t, const Tube &v) const
 Returns the optimal evaluation of this tube at \(t\), based on the derivative information \(\dot{x}(\cdot)\).
 
const Interval interpol (const Interval &t, const Tube &v) const
 Returns the optimal interval evaluation of this tube over \([t]\), based on the derivative information \(\dot{x}(\cdot)\).
 
const Interval invert (const Interval &y, const Interval &search_tdomain=Interval::ALL_REALS) const
 Returns the interval inversion \([x]^{-1}([y])\).
 
void invert (const Interval &y, std::vector< Interval > &v_t, const Interval &search_tdomain=Interval::ALL_REALS) const
 Computes the set of continuous values of the inversion \([x]^{-1}([y])\).
 
const Interval invert (const Interval &y, const Tube &v, const Interval &search_tdomain=Interval::ALL_REALS) const
 Returns the optimal interval inversion \([x]^{-1}([y])\).
 
void invert (const Interval &y, std::vector< Interval > &v_t, const Tube &v, const Interval &search_tdomain=Interval::ALL_REALS) const
 Computes the set of continuous values of the optimal inversion \([x]^{-1}([y])\).
 
double max_diam () const
 Returns the diameter of the interval value \([x](t)\) that is the more uncertain.
 
double max_gate_diam (double &t) const
 Returns the diameter of the gate of this tube that is the more uncertain.
 
const Trajectory diam (bool gates_thicknesses=false) const
 Returns the diameters of the tube as a trajectory.
 
const Trajectory diam (const Tube &v) const
 Returns the diameters of the tube as a trajectory.
 
Tests
bool operator== (const Tube &x) const
 Returns true if this tube is equal to \([x](\cdot)\).
 
bool operator!= (const Tube &x) const
 Returns true if this tube is different from \([x](\cdot)\).
 
bool is_subset (const Tube &x) const
 Returns true if this tube is a subset of \([x](\cdot)\).
 
bool is_strict_subset (const Tube &x) const
 Returns true if this tube is a subset of \([x](\cdot)\), and not \([x](\cdot)\) itself.
 
bool is_interior_subset (const Tube &x) const
 Returns true if this tube is a subset of the interior of \([x](\cdot)\).
 
bool is_strict_interior_subset (const Tube &x) const
 Returns true if this tube is a subset of the interior of \([x](\cdot)\), and not \([x](\cdot)\) itself.
 
bool is_superset (const Tube &x) const
 Returns true if this tube is a superset of \([x](\cdot)\).
 
bool is_strict_superset (const Tube &x) const
 Returns true if this tube is a superset of \([x](\cdot)\), and not \([x](\cdot)\) itself.
 
bool is_empty () const
 Returns true if this tube is empty.
 
const BoolInterval contains (const Trajectory &x) const
 Returns true if this tube contains the trajectory \(x(\cdot)\).
 
bool overlaps (const Tube &x, float ratio=1.) const
 Returns true if this tube overlaps the tube \([x](\cdot)\).
 
Setting values
const Tubeset (const Interval &y)
 Sets a constant interval value for this tube: \(\forall t, [x](t)=[y]\).
 
const Tubeset (const Interval &y, int slice_id)
 Sets the interval value of the ith slice of this tube.
 
const Tubeset (const Interval &y, double t)
 Sets the interval value of this tube at \(t\): \([x](t)=[y]\).
 
const Tubeset (const Interval &y, const Interval &t)
 Sets the interval value of this tube over \([t]\): \(\forall t\in[t], [x](t)=[y]\).
 
const Tubeset_empty ()
 Sets this tube to the empty set.
 
const Tubeinflate (double rad)
 Inflates this tube by adding \([-rad,+rad]\) to all its codomain components.
 
const Tubeinflate (const Trajectory &rad)
 Inflates this tube by adding non-constant uncertainties defined in a trajectory.
 
Tubetruncate_tdomain (const Interval &tdomain)
 Truncates the tdomain of \([x](\cdot)\).
 
void shift_tdomain (double a)
 Shifts the tdomain \([t_0,t_f]\) of \([x](\cdot)\).
 
Bisection
const std::pair< Tube, Tubebisect (double t, float ratio=0.49) const
 Bisects this tube.
 
Assignments operators
const Tubeoperator+= (const Interval &x)
 Operates +=.
 
const Tubeoperator+= (const Trajectory &x)
 Operates +=.
 
const Tubeoperator+= (const Tube &x)
 Operates +=.
 
const Tubeoperator-= (const Interval &x)
 Operates -=.
 
const Tubeoperator-= (const Trajectory &x)
 Operates -=.
 
const Tubeoperator-= (const Tube &x)
 Operates -=.
 
const Tubeoperator*= (const Interval &x)
 Operates *=.
 
const Tubeoperator*= (const Trajectory &x)
 Operates *=.
 
const Tubeoperator*= (const Tube &x)
 Operates *=.
 
const Tubeoperator/= (const Interval &x)
 Operates /=.
 
const Tubeoperator/= (const Trajectory &x)
 Operates /=.
 
const Tubeoperator/= (const Tube &x)
 Operates /=.
 
const Tubeoperator|= (const Interval &x)
 Operates |=.
 
const Tubeoperator|= (const Trajectory &x)
 Operates |=.
 
const Tubeoperator|= (const Tube &x)
 Operates |=.
 
const Tubeoperator&= (const Interval &x)
 Operates &=.
 
const Tubeoperator&= (const Trajectory &x)
 Operates &=.
 
const Tubeoperator&= (const Tube &x)
 Operates &=.
 
Integration
const Interval integral (double t) const
 Computes the interval integral \(\int_0^t[x](\tau)d\tau\).
 
const Interval integral (const Interval &t) const
 Computes the interval integral \(\int_0^{[t]}[x](\tau)d\tau\).
 
const Interval integral (const Interval &t1, const Interval &t2) const
 Computes the interval integral \(\int_{[t_1]}^{[t_2]}[x](\tau)d\tau\).
 
const std::pair< Interval, Interval > partial_integral (const Interval &t) const
 Computes the partial interval integral \(\int_{0}^{[t]}[x](\tau)d\tau\).
 
const std::pair< Interval, Interval > partial_integral (const Interval &t1, const Interval &t2) const
 Computes the partial interval integral \(\int_{[t_1]}^{[t_2]}[x](\tau)d\tau\).
 
Serialization
void serialize (const std::string &binary_file_name="x.tube", int version_number=SERIALIZATION_VERSION) const
 Serializes this tube.
 
void serialize (const std::string &binary_file_name, const Trajectory &traj, int version_number=SERIALIZATION_VERSION) const
 Serializes this tube together with a Trajectory object.
 
- Public Member Functions inherited from codac::DynamicalItem
virtual ~DynamicalItem ()
 DynamicalItem destructor.
 

Static Public Member Functions

static bool same_slicing (const Tube &x1, const Tube &x2)
 Tests whether the two Tube objects are sharing the same slicing.
 
static void enable_syntheses (bool enable=true)
 Enables the computation of a synthesis tree for any Tube object.
 
static const Tube hull (const std::list< Tube > &l_tubes)
 Computes the hull of several tubes.
 
- Static Public Member Functions inherited from codac::DynamicalItem
static bool valid_tdomain (const Interval &tdomain)
 Verifies that this interval is a feasible tdomain.
 

Protected Member Functions

 Tube ()
 Creates a not-defined scalar tube.
 
const IntervalVector codomain_box () const
 Returns the box \([x]([t_0,t_f])\).
 
void deserialize (const std::string &binary_file_name, Trajectory *&traj)
 Restores a scalar tube from serialization, together with a Trajectory object.
 
void create_synthesis_tree () const
 Creates the synthesis tree associated to the values of this tube.
 
void delete_synthesis_tree () const
 Deletes the synthesis tree of this tube.
 
void create_polynomial_synthesis (double eps) const
 Creates the synthesis tree associated to the values of this tube.
 
void delete_polynomial_synthesis () const
 Deletes the synthesis tree of this tube.
 

Protected Attributes

Slicem_first_slice = nullptr
 pointer to the first Slice object of this tube
 
TubeTreeSynthesis * m_synthesis_tree = nullptr
 pointer to the optional synthesis tree
 
TubePolynomialSynthesis * m_polynomial_synthesis = nullptr
 pointer to the optional synthesis tree
 
SynthesisMode m_synthesis_mode = SynthesisMode::NONE
 enables of the use of a synthesis tree
 
Interval m_tdomain
 redundant information for fast evaluations
 

String

const std::string class_name () const
 Returns the name of this class.
 
void enable_synthesis (SynthesisMode mode=SynthesisMode::BINARY_TREE, double eps=1.e-3) const
 Enables the computation of a synthesis tree.
 

Detailed Description

One dimensional tube \([x](\cdot)\), defined as an interval of scalar trajectories.

Note
Use TubeVector for the multi-dimensional case

Constructor & Destructor Documentation

◆ Tube() [1/11]

codac::Tube::Tube ( const Interval & tdomain,
const Interval & codomain = Interval::ALL_REALS )
explicit

Creates a scalar tube \([x](\cdot)\) made of one slice.

Parameters
tdomaintemporal domain \([t_0,t_f]\)
codomainInterval value of the slice (all reals \([-\infty,\infty]\) by default)

◆ Tube() [2/11]

codac::Tube::Tube ( const Interval & tdomain,
double timestep,
const Interval & codomain = Interval::ALL_REALS )
explicit

Creates a scalar tube \([x](\cdot)\) with some temporal discretization.

Parameters
tdomaintemporal domain \([t_0,t_f]\)
timestepsampling value \(\delta\) for the temporal discretization (double)
codomainInterval value of the slices (all reals \([-\infty,\infty]\) by default)

◆ Tube() [3/11]

codac::Tube::Tube ( const Interval & tdomain,
double timestep,
const TFnc & f,
int f_image_id = 0 )
explicit

Creates a scalar tube \([x](\cdot)\) from a TFnc object and with some temporal discretization.

Note
Due to the slicing implementation of the tube, a wrapping effect will occur to reliably enclose the TFnc object
Parameters
tdomaintemporal domain \([t_0,t_f]\)
timestepsampling value \(\delta\) for the temporal discretization (double)
fTFnc object that will be enclosed by the tube: \(\forall t\in[t_0,t_f], [f](t)\subseteq[x](t)\)
f_image_idcomponent index of the interval function \([f]\) (that is possibly multidimensional, first component by default)

◆ Tube() [4/11]

codac::Tube::Tube ( const std::vector< Interval > & v_tdomains,
const std::vector< Interval > & v_codomains )
explicit

Creates a tube \([x](\cdot)\) from a list of \(k\) boxes \(\big([t_1]\times[x_1],\dots,[t_k]\times[x_k]\big)\).

Note
The slicing will be based on the vector of temporal domains.
The \([t_i]\)'s must cover continuously the tdomain of \([x](\cdot)\).
Parameters
v_tdomainsvector of temporal domains \([t_i]\)
v_codomainsvector of codomains \([x_i]\) related to the \([t_i]\)'s

◆ Tube() [5/11]

codac::Tube::Tube ( const Tube & x)

Creates a copy of a scalar tube \([x](\cdot)\), with the same time discretization.

Parameters
xTube to be duplicated

◆ Tube() [6/11]

codac::Tube::Tube ( const Tube & x,
const TFnc & f,
int f_image_id = 0 )
explicit

Creates a copy of a scalar tube \([x](\cdot)\), with the same time discretization but a specific codomain defined by a TFnc object.

Note
Due to the slicing implementation of the tube, a wrapping effect will occur to reliably enclose the TFnc object
Parameters
xTube from which the sampling will be duplicated
fTFnc object that will be enclosed by the tube: \(\forall t\in[t_0,t_f], [f](t)\subseteq[x](t)\)
f_image_idcomponent index of the interval function \([f]\) (that is possibly multidimensional, first component by default)

◆ Tube() [7/11]

codac::Tube::Tube ( const Trajectory & traj,
double timestep )
explicit

Creates a scalar tube \([x](\cdot)\) enclosing a trajectory \(x(\cdot)\), possibly with some temporal discretization.

Note
Due to the slicing implementation of the tube, a wrapping effect will occur to reliably enclose the Trajectory object
Parameters
trajTrajectory \(x(\cdot)\) to enclose
timestepsampling value \(\delta\) for the temporal discretization (double, no discretization by default: one slice only)

◆ Tube() [8/11]

codac::Tube::Tube ( const Trajectory & lb,
const Trajectory & ub,
double timestep )
explicit

Creates a scalar tube \([x](\cdot)\) defined as an interval of two trajectories \([lb(\cdot),ub(\cdot)]\).

Note
Due to the slicing implementation of the tube, a wrapping effect will occur to reliably enclose the Trajectory object
Parameters
lbTrajectory defining the lower bound \(x^{-}(\cdot)\) of the tube
ubTrajectory defining the upper bound \(x^{+}(\cdot)\) of the tube
timestepsampling value \(\delta\) for the temporal discretization (double, no discretization by default: one slice only)

◆ Tube() [9/11]

codac::Tube::Tube ( const std::string & binary_file_name)
explicit

Restore a scalar tube from serialization.

Note
The Tube must have been serialized beforehand by the appropriate method serialize()
Parameters
binary_file_namepath to the binary file

◆ Tube() [10/11]

codac::Tube::Tube ( const std::string & binary_file_name,
Trajectory *& traj )
explicit

Restore a scalar tube from serialization, together with a Trajectory object.

Note
The Tube and the Trajectory must have been serialized beforehand by the appropriate method serialize()
Parameters
binary_file_namepath to the binary file
traja pointer to the Trajectory object to be instantiated

◆ Tube() [11/11]

codac::Tube::Tube ( )
protected

Creates a not-defined scalar tube.

Note
Constructor necessary for the TubeVector class

Member Function Documentation

◆ size()

int codac::Tube::size ( ) const
virtual

Returns the dimension of the scalar tube (always 1)

Returns
1

Implements codac::DynamicalItem.

◆ primitive()

const Tube codac::Tube::primitive ( const Interval & c = Interval(0.)) const

Returns the primitive Tube of this tube.

Parameters
cthe constant of integration (0. by default)
Returns
a new Tube object with same slicing, enclosing the feasible primitives of this tube

◆ operator=()

const Tube & codac::Tube::operator= ( const Tube & x)

Returns a copy of a Tube.

Parameters
xthe Tube object to be copied
Returns
a new Tube object with same slicing and values

◆ tdomain()

const Interval codac::Tube::tdomain ( ) const
virtual

Returns the temporal definition domain of this tube.

Returns
an Interval object \([t_0,t_f]\)

Implements codac::DynamicalItem.

◆ polygon_envelope()

const Polygon codac::Tube::polygon_envelope ( ) const

Returns the polygon envelope of this tube.

Returns
a Polygon object enclosing the slices

◆ lb()

const Trajectory codac::Tube::lb ( ) const

Returns a possible lower bound \(x^{-}(\cdot)\) of the tube.

Note
The exact lower bound cannot be known. However, the returned trajectory \(x^{-}(\cdot)\) is guaranteed to be enclosed in the tube \([x](\cdot)\).
Returns
a candidate for \(x^{-}(\cdot)\)

◆ ub()

const Trajectory codac::Tube::ub ( ) const

Returns a possible upper bound \(x^{+}(\cdot)\) of the tube.

Note
The exact upper bound cannot be known. However, the returned trajectory \(x^{+}(\cdot)\) is guaranteed to be enclosed in the tube \([x](\cdot)\).
Returns
a candidate for \(x^{+}(\cdot)\)

◆ nb_slices()

int codac::Tube::nb_slices ( ) const

Returns the number of slices of this tube.

Returns
an integer

◆ slice() [1/4]

Slice * codac::Tube::slice ( int slice_id)

Returns a pointer to the ith Slice object of this tube.

Parameters
slice_idthe index of the ith Slice
Returns
a pointer to the corresponding Slice

◆ slice() [2/4]

const Slice * codac::Tube::slice ( int slice_id) const

Returns a constant pointer to the ith Slice object of this tube.

Parameters
slice_idthe index of the ith Slice
Returns
a const pointer to the corresponding Slice

◆ slice() [3/4]

Slice * codac::Tube::slice ( double t)

Returns a pointer to the Slice object of this tube that is defined at \(t\).

Note
If two Slices are defined at \(t\) (common tdomain), then the first Slice is considered
Parameters
tthe temporal key (double, must belong to the Tube's tdomain)
Returns
a pointer to the corresponding Slice

◆ slice() [4/4]

const Slice * codac::Tube::slice ( double t) const

Returns a constant pointer to the Slice object of this tube that is defined at \(t\).

Note
If two Slices are defined at \(t\) (common tdomain), then the first Slice is considered
Parameters
tthe temporal key (double, must belong to the Tube's tdomain)
Returns
a const pointer to the corresponding Slice

◆ first_slice() [1/2]

Slice * codac::Tube::first_slice ( )

Returns a pointer to the first Slice object of this tube.

Returns
a pointer to the corresponding Slice

◆ first_slice() [2/2]

const Slice * codac::Tube::first_slice ( ) const

Returns a constant pointer to the first Slice object of this tube.

Returns
a const pointer to the corresponding Slice

◆ last_slice() [1/2]

Slice * codac::Tube::last_slice ( )

Returns a pointer to the last Slice object of this tube.

Returns
a pointer to the corresponding Slice

◆ last_slice() [2/2]

const Slice * codac::Tube::last_slice ( ) const

Returns a constant pointer to the last Slice object of this tube.

Returns
a const pointer to the corresponding Slice

◆ wider_slice() [1/2]

Slice * codac::Tube::wider_slice ( )

Returns a pointer to the Slice defined over the wider temporal domain.

Note
If two Slice objects have the same tdomain width, then the first one is considered
Returns
a pointer to the corresponding Slice

◆ wider_slice() [2/2]

const Slice * codac::Tube::wider_slice ( ) const

Returns a constant pointer to the Slice defined over the wider temporal domain.

Note
If two Slice objects have the same tdomain width, then the first one is considered
Returns
a const pointer to the corresponding Slice

◆ largest_slice() [1/2]

Slice * codac::Tube::largest_slice ( )

Returns a pointer to the Slice object of this tube for which the interval value is the most uncertain.

Returns
a pointer to the corresponding Slice

◆ largest_slice() [2/2]

const Slice * codac::Tube::largest_slice ( ) const

Returns a constant pointer to the Slice object of this tube for which the interval value is the most uncertain.

Returns
a const pointer to the corresponding Slice

◆ steepest_slice() [1/2]

Slice * codac::Tube::steepest_slice ( )

Returns a pointer to the Slice object of this tube for which the difference between the mid values of the input and output gates is the largest.

Returns
a pointer to the corresponding Slice

◆ steepest_slice() [2/2]

const Slice * codac::Tube::steepest_slice ( ) const

Returns a const pointer to the Slice object of this tube for which the difference between the mid values of the input and output gates is the largest.

Returns
a const pointer to the corresponding Slice

◆ slice_tdomain()

const Interval codac::Tube::slice_tdomain ( int slice_id) const

Returns the temporal definition domain of the ith Slice of this tube.

Parameters
slice_idthe index of the ith Slice
Returns
an Interval object \([t_0^i,t_f^i]\)

◆ time_to_index()

int codac::Tube::time_to_index ( double t) const

Returns the Slice index related to the temporal key \(t\).

Parameters
tthe temporal key (double, must belong to the Tube's tdomain)
Returns
an integer

◆ index()

int codac::Tube::index ( const Slice * slice) const

Returns the Slice index related to the Slice pointer.

Parameters
slicea const pointer to a Slice object of this tube
Returns
an integer

◆ sample() [1/4]

void codac::Tube::sample ( double t)

Samples this tube at \(t\).

Note
Without any effect if two Slice objects are already defined at \(t\) (if the gate \([x](t)\) already exists)
Parameters
tthe temporal key (double, must belong to the Tube's tdomain)

◆ sample() [2/4]

void codac::Tube::sample ( double t,
Slice * slice_to_be_sampled )

Samples this tube at \(t\) from a pointer to the corresponding slice.

Reduces the complexity of related methods by providing a direct access to the Slice object covering the input \(t\).

Note
Without any effect if two Slice objects are already defined at \(t\) (if the gate \([x](t)\) already exists)
Parameters
tthe temporal key (double, must belong to the Tube's tdomain)
slice_to_be_sampleda pointer to the Slice whose tdomain contains \(t\)

◆ sample() [3/4]

void codac::Tube::sample ( double t,
const Interval & gate )

Samples this tube at \(t\) with a specific gate value.

Note
Without any sampling effect if two Slice objects are already defined at \(t\) (if the gate \([x](t)\) already exists)
Parameters
tthe temporal key (double, must belong to the Tube's tdomain)
gatethe Interval value of this tube at \(t\)

◆ sample() [4/4]

void codac::Tube::sample ( const Tube & x)

Samples this tube so that it will share the same sampling of \([x](\cdot)\).

Note
The previous sampling of this tube is preserved
Parameters
xthe Tube from which the new sampling will come from

◆ gate_exists()

bool codac::Tube::gate_exists ( double t) const

Tests if a gate exists at time \(t\).

Parameters
ttime input to test
Returns
true if a gate exists

◆ remove_gate()

void codac::Tube::remove_gate ( double t)

Removes the gate at \(t\) and merges the two related slices.

Parameters
ttime input where the gate to remove is

◆ merge_similar_slices()

void codac::Tube::merge_similar_slices ( double distance_threshold)

Merges all adjacent slices whose Hausdorff distance is less than the given threshold.

Note
Adjacent slices whose Hausdorff distance is larger (or equal) than the given threshold are not merged
Parameters
distance_thresholdthe threshold for the maximum Haussdorf distance between adjacent slices

◆ codomain()

const Interval codac::Tube::codomain ( ) const

Returns the interval of feasible values.

Returns
an Interval object \([x]([t_0,t_f])\)

◆ volume()

double codac::Tube::volume ( ) const

Returns the volume of this tube.

Note
returns POS_INFINITY if the codomain is unbounded
returns 0 if the tube is flat (and so without wrapping effect)
Returns
volume defined as \(w([t_0,t_f])\times w([x]([t_0,t_f]))\)

◆ operator()() [1/3]

const Interval codac::Tube::operator() ( int slice_id) const

Returns the value of the ith slice.

Parameters
slice_idthe index of the ith slice
Returns
Interval value of \([x](i)\)

◆ operator()() [2/3]

const Interval codac::Tube::operator() ( double t) const

Returns the evaluation of this tube at \(t\).

Parameters
tthe temporal key (double, must belong to the Tube's tdomain)
Returns
Interval value of \([x](t)\)

◆ operator()() [3/3]

const Interval codac::Tube::operator() ( const Interval & t) const

Returns the interval evaluation of this tube over \([t]\).

Parameters
tthe subtdomain (Interval, must be a subset of the Tube's tdomain)
Returns
Interval envelope \([x]([t])\)

◆ eval()

const std::pair< Interval, Interval > codac::Tube::eval ( const Interval & t = Interval::ALL_REALS) const

Returns the interval evaluations of the bounds of the tube \(\underline{x^-}(\cdot)\) and \(\overline{x^+}(\cdot)\) over \([t]\).

Parameters
tthe subtdomain (Interval, must be a subset of the Tube's tdomain)
Returns
the pair \(\big([\underline{x^-}]([t]),[\overline{x^+}]([t])\big)\)

◆ interpol() [1/2]

const Interval codac::Tube::interpol ( double t,
const Tube & v ) const

Returns the optimal evaluation of this tube at \(t\), based on the derivative information \(\dot{x}(\cdot)\).

Todo
Change the name of this method?
Parameters
tthe temporal key (double, must belong to the Tube's tdomain)
vthe derivative tube such that \(\dot{x}(\cdot)\in[v](\cdot)\)
Returns
Interval value of \([x](t)\)

◆ interpol() [2/2]

const Interval codac::Tube::interpol ( const Interval & t,
const Tube & v ) const

Returns the optimal interval evaluation of this tube over \([t]\), based on the derivative information \(\dot{x}(\cdot)\).

Todo
Change the name of this method?
Parameters
tthe subtdomain (Interval, must be a subset of the Tube's tdomain)
vthe derivative tube such that \(\dot{x}(\cdot)\in[v](\cdot)\)
Returns
Interval value of \([x]([t])\)

◆ invert() [1/4]

const Interval codac::Tube::invert ( const Interval & y,
const Interval & search_tdomain = Interval::ALL_REALS ) const

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

Note
If the inversion results in several pre-images, their union is returned
Parameters
ythe interval codomain
search_tdomainthe optional temporal domain on which the inversion will be performed
Returns
the hull of \([x]^{-1}([y])\)

◆ invert() [2/4]

void codac::Tube::invert ( const Interval & y,
std::vector< Interval > & v_t,
const Interval & search_tdomain = Interval::ALL_REALS ) const

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

Parameters
ythe interval codomain
v_tthe vector of the sub-tdomains \([t_k]\) for which \(\forall t\in[t_k] \mid x(t)\in[y], x(\cdot)\in[x](\cdot)\)
search_tdomainthe optional temporal domain on which the inversion will be performed

◆ invert() [3/4]

const Interval codac::Tube::invert ( const Interval & y,
const Tube & v,
const Interval & search_tdomain = Interval::ALL_REALS ) const

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

Note
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
ythe interval codomain
vthe derivative tube such that \(\dot{x}(\cdot)\in[v](\cdot)\)
search_tdomainthe optional temporal domain on which the inversion will be performed
Returns
the hull of \([x]^{-1}([y])\)

◆ invert() [4/4]

void codac::Tube::invert ( const Interval & y,
std::vector< Interval > & v_t,
const Tube & v,
const Interval & search_tdomain = Interval::ALL_REALS ) const

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

Note
The knowledge of the derivative tube \([v](\cdot)\) allows finer inversions
Parameters
ythe interval codomain
v_tthe vector of the sub-tdomains \([t_k]\) for which \(\exists t\in[t_k] \mid x(t)\in[y], x(\cdot)\in[x](\cdot), \dot{x}(\cdot)\in[v](\cdot)\)
vthe derivative tube such that \(\dot{x}(\cdot)\in[v](\cdot)\)
search_tdomainthe optional temporal domain on which the inversion will be performed

◆ max_diam()

double codac::Tube::max_diam ( ) const

Returns the diameter of the interval value \([x](t)\) that is the more uncertain.

Returns
the maximal thickness of this tube

◆ max_gate_diam()

double codac::Tube::max_gate_diam ( double & t) const

Returns the diameter of the gate of this tube that is the more uncertain.

Parameters
tthe temporal key of the corresponding uncertain gate
Returns
the maximal thickness of the gate

◆ diam() [1/2]

const Trajectory codac::Tube::diam ( bool gates_thicknesses = false) const

Returns the diameters of the tube as a trajectory.

Note
Without derivative knowledge, and because the tube is made of boxed slices, the trajectory will be discontinuous and so the returned object will not reliably represent the diameters. It can be mainly used for display purposes.
Parameters
gates_thicknessesif true, the diameters of the gates will be evaluated too
Returns
the set of diameters associated to temporal inputs

◆ diam() [2/2]

const Trajectory codac::Tube::diam ( const Tube & v) const

Returns the diameters of the tube as a trajectory.

Note
Because the tube is made of boxed slices, the trajectory may be discontinuous and so the returned object will not reliably represent the diameters. It can be mainly used for display purposes.
Parameters
vthe derivative tube such that \(\dot{x}(\cdot)\in[v](\cdot)\)
Returns
the set of diameters associated to temporal inputs

◆ operator==()

bool codac::Tube::operator== ( const Tube & x) const

Returns true if this tube is equal to \([x](\cdot)\).

Note
Equality is obtained if the tubes share the same bounds, tdomain and sampling
Parameters
xthe Tube object
Returns
true in case of equality

◆ operator!=()

bool codac::Tube::operator!= ( const Tube & x) const

Returns true if this tube is different from \([x](\cdot)\).

Note
The two tubes are different if they do not share the same bounds, tdomain or sampling
Parameters
xthe Tube object
Returns
true in case of difference

◆ is_subset()

bool codac::Tube::is_subset ( const Tube & x) const

Returns true if this tube is a subset of \([x](\cdot)\).

Note
The two tubes may not share the same slicing, but must have the same tdomain
Parameters
xthe Tube object
Returns
true in case of subset

◆ is_strict_subset()

bool codac::Tube::is_strict_subset ( const Tube & x) const

Returns true if this tube is a subset of \([x](\cdot)\), and not \([x](\cdot)\) itself.

Note
The two tubes may not share the same slicing, but must have the same tdomain
Parameters
xthe Tube object
Returns
true in case of strict subset

◆ is_interior_subset()

bool codac::Tube::is_interior_subset ( const Tube & x) const

Returns true if this tube is a subset of the interior of \([x](\cdot)\).

Note
The two tubes may not share the same slicing, but must have the same tdomain
Parameters
xthe Tube object
Returns
true in case of interior subset

◆ is_strict_interior_subset()

bool codac::Tube::is_strict_interior_subset ( const Tube & x) const

Returns true if this tube is a subset of the interior of \([x](\cdot)\), and not \([x](\cdot)\) itself.

Note
The two tubes may not share the same slicing, but must have the same tdomain
Parameters
xthe Tube object
Returns
true in case of strict interior subset

◆ is_superset()

bool codac::Tube::is_superset ( const Tube & x) const

Returns true if this tube is a superset of \([x](\cdot)\).

Note
The two tubes may not share the same slicing, but must have the same tdomain
Parameters
xthe Tube object
Returns
true in case of superset

◆ is_strict_superset()

bool codac::Tube::is_strict_superset ( const Tube & x) const

Returns true if this tube is a superset of \([x](\cdot)\), and not \([x](\cdot)\) itself.

Note
The two tubes may not share the same slicing, but must have the same tdomain
Parameters
xthe Tube object
Returns
true in case of strict superset

◆ is_empty()

bool codac::Tube::is_empty ( ) const

Returns true if this tube is empty.

Note
If \([x](t)=\varnothing\) for some \(t\), then the tube \([x](\cdot)\) is considered empty
Returns
true in case of emptiness

◆ contains()

const BoolInterval codac::Tube::contains ( const Trajectory & x) const

Returns true if this tube contains the trajectory \(x(\cdot)\).

Note
Due to the reliable numerical representation of a trajectory, some wrapping effect may appear for its evaluations (either if it is defined by a map of values or an analytic function). Hence, this "contains" test may not be able to conclude, if the thin envelope of \(x(\cdot)\) overlaps a boundary of the tube.
Parameters
xthe trajectory that might be contained by this tube
Returns
BoolInterval::YES (or BoolInterval::NO) if this tube contains \(x(\cdot)\) (or does not contain) and BoolInterval::MAYBE in case of ambiguity

◆ overlaps()

bool codac::Tube::overlaps ( const Tube & x,
float ratio = 1. ) const

Returns true if this tube overlaps the tube \([x](\cdot)\).

Parameters
xthe other tube
ratioan optional overlapping ratio between 0 and 1 (1 by default). For instance, if ratio=0.3 and there is an overlapping of at least 30%, then the function returns true
Returns
true in case of overlapping with respect to the defined ratio

◆ set() [1/4]

const Tube & codac::Tube::set ( const Interval & y)

Sets a constant interval value for this tube: \(\forall t, [x](t)=[y]\).

Note
The sampling of this tube is preserved
Parameters
yInterval value of the slices
Returns
*this

◆ set() [2/4]

const Tube & codac::Tube::set ( const Interval & y,
int slice_id )

Sets the interval value of the ith slice of this tube.

Parameters
yInterval value of the ith slice
slice_idindex of the ith Slice
Returns
*this

◆ set() [3/4]

const Tube & codac::Tube::set ( const Interval & y,
double t )

Sets the interval value of this tube at \(t\): \([x](t)=[y]\).

Note
It may create a gate (and so one more slice) if the tube is not already sampled at \(t\). Otherwise, it will update the value of the already existing gate.
Parameters
yInterval value of the gate
tthe temporal key (double, must belong to the Tube's tdomain)
Returns
*this

◆ set() [4/4]

const Tube & codac::Tube::set ( const Interval & y,
const Interval & t )

Sets the interval value of this tube over \([t]\): \(\forall t\in[t], [x](t)=[y]\).

Note
It may create two gates (and so further slices) if the tube is not already sampled at \(t^-\) and \(t^+\). This is done to ensure that \(\forall t\in[t], [x](t)=[y]\).
Parameters
yInterval value to be set
tthe subtdomain (Interval, must be a subset of the Tube's tdomain)
Returns
*this

◆ set_empty()

const Tube & codac::Tube::set_empty ( )

Sets this tube to the empty set.

Note
By convention, all slices will be set to the empty set
Returns
*this

◆ inflate() [1/2]

const Tube & codac::Tube::inflate ( double rad)

Inflates this tube by adding \([-rad,+rad]\) to all its codomain components.

Note
All slices and gates will be inflated
Parameters
radhalf of the inflation
Returns
*this

◆ inflate() [2/2]

const Tube & codac::Tube::inflate ( const Trajectory & rad)

Inflates this tube by adding non-constant uncertainties defined in a trajectory.

Note
From the trajectory \(a(\cdot)\), the function will inflate this tube such that \(\forall t, [x](t):=[x](t)+[-a(t),+a(t)]\)
All slices and gates will be inflated
Parameters
radthe Trajectory object defining the non-constant inflation
Returns
*this

◆ truncate_tdomain()

Tube & codac::Tube::truncate_tdomain ( const Interval & tdomain)

Truncates the tdomain of \([x](\cdot)\).

Note
The new tdomain must be a subset of the old one
Parameters
tdomainnew temporal domain \([t_0,t_f]\)
Returns
a reference to this tube

◆ shift_tdomain()

void codac::Tube::shift_tdomain ( double a)

Shifts the tdomain \([t_0,t_f]\) of \([x](\cdot)\).

Parameters
athe offset value so that \([t_0,t_f]:=[t_0+a,t_f+a]\)

◆ bisect()

const std::pair< Tube, Tube > codac::Tube::bisect ( double t,
float ratio = 0.49 ) const

Bisects this tube.

Note
The tube is bisected along the codomain and according to a defined ratio
The bisection is performed on the interval value of the gate \([x](t)\)
If the tube is not already sampled at \(t\), then a sampling is performed
Parameters
tthe temporal key (double, must belong to the Tube's tdomain)
ratiothe bisection ratio (default value: 0.55)
Returns
a pair of two Tube objects resulting from the bisection

◆ operator+=() [1/3]

const Tube & codac::Tube::operator+= ( const Interval & x)

Operates +=.

Parameters
xInterval
Returns
(*this)+=x

◆ operator+=() [2/3]

const Tube & codac::Tube::operator+= ( const Trajectory & x)

Operates +=.

Parameters
xTrajectory
Returns
(*this)+=x

◆ operator+=() [3/3]

const Tube & codac::Tube::operator+= ( const Tube & x)

Operates +=.

Parameters
xTube
Returns
(*this)+=x

◆ operator-=() [1/3]

const Tube & codac::Tube::operator-= ( const Interval & x)

Operates -=.

Parameters
xInterval
Returns
(*this)-=x

◆ operator-=() [2/3]

const Tube & codac::Tube::operator-= ( const Trajectory & x)

Operates -=.

Parameters
xTrajectory
Returns
(*this)-=x

◆ operator-=() [3/3]

const Tube & codac::Tube::operator-= ( const Tube & x)

Operates -=.

Parameters
xTube
Returns
(*this)-=x

◆ operator*=() [1/3]

const Tube & codac::Tube::operator*= ( const Interval & x)

Operates *=.

Parameters
xInterval
Returns
(this)=x

◆ operator*=() [2/3]

const Tube & codac::Tube::operator*= ( const Trajectory & x)

Operates *=.

Parameters
xTrajectory
Returns
(this)=x

◆ operator*=() [3/3]

const Tube & codac::Tube::operator*= ( const Tube & x)

Operates *=.

Parameters
xTube
Returns
(this)=x

◆ operator/=() [1/3]

const Tube & codac::Tube::operator/= ( const Interval & x)

Operates /=.

Parameters
xInterval
Returns
(*this)/=x

◆ operator/=() [2/3]

const Tube & codac::Tube::operator/= ( const Trajectory & x)

Operates /=.

Parameters
xTrajectory
Returns
(*this)/=x

◆ operator/=() [3/3]

const Tube & codac::Tube::operator/= ( const Tube & x)

Operates /=.

Parameters
xTube
Returns
(*this)/=x

◆ operator|=() [1/3]

const Tube & codac::Tube::operator|= ( const Interval & x)

Operates |=.

Parameters
xInterval
Returns
(*this)|=x

◆ operator|=() [2/3]

const Tube & codac::Tube::operator|= ( const Trajectory & x)

Operates |=.

Parameters
xTrajectory
Returns
(*this)|=x

◆ operator|=() [3/3]

const Tube & codac::Tube::operator|= ( const Tube & x)

Operates |=.

Parameters
xTube
Returns
(*this)|=x

◆ operator&=() [1/3]

const Tube & codac::Tube::operator&= ( const Interval & x)

Operates &=.

Parameters
xInterval
Returns
(*this)&=x

◆ operator&=() [2/3]

const Tube & codac::Tube::operator&= ( const Trajectory & x)

Operates &=.

Parameters
xTrajectory
Returns
(*this)&=x

◆ operator&=() [3/3]

const Tube & codac::Tube::operator&= ( const Tube & x)

Operates &=.

Parameters
xTube
Returns
(*this)&=x

◆ class_name()

const std::string codac::Tube::class_name ( ) const
inlinevirtual

Returns the name of this class.

Note
Only used for some generic display method
Returns
the predefined name

Implements codac::DynamicalItem.

1016{ return "Tube"; };

◆ enable_synthesis()

void codac::Tube::enable_synthesis ( SynthesisMode mode = SynthesisMode::BINARY_TREE,
double eps = 1.e-3 ) const

Enables the computation of a synthesis tree.

Note
The synthesis tree speeds up computations such as integrals or evaluations
Parameters
modemode of synthesis
epsprecision of the polynomial approximation, if selected

◆ integral() [1/3]

const Interval codac::Tube::integral ( double t) const

Computes the interval integral \(\int_0^t[x](\tau)d\tau\).

Note
From the monotonicity of the integral operator, \(\int_0^t[x](\tau)d\tau=[\int_0^t x^-(\tau)d\tau,\int_0^t x^+(\tau)d\tau]\)
Parameters
tthe temporal key (double, must belong to the Tube's tdomain)
Returns
the set of feasible integral values

◆ integral() [2/3]

const Interval codac::Tube::integral ( const Interval & t) const

Computes the interval integral \(\int_0^{[t]}[x](\tau)d\tau\).

Note
From the monotonicity of the integral operator, \(\int_0^{[t]}[x](\tau)d\tau=[\int_0^{[t]}x^-(\tau)d\tau,\int_0^{[t]}x^+(\tau)d\tau]\)
Parameters
tthe subtdomain (Interval, must be a subset of the Tube's tdomain)
Returns
the set of feasible integral values

◆ integral() [3/3]

const Interval codac::Tube::integral ( const Interval & t1,
const Interval & t2 ) const

Computes the interval integral \(\int_{[t_1]}^{[t_2]}[x](\tau)d\tau\).

Note
From the monotonicity of the integral operator, \(\int_{[t_1]}^{[t_2]}[x](\tau)d\tau=[\int_{[t_1]}^{[t_2]}x^-(\tau)d\tau,\int_{[t_1]}^{[t_2]}x^+(\tau)d\tau]\)
Parameters
t1lower bound, subset of the Tube's tdomain
t2upper bound, subset of the Tube's tdomain
Returns
the set of feasible integral values

◆ partial_integral() [1/2]

const std::pair< Interval, Interval > codac::Tube::partial_integral ( const Interval & t) const

Computes the partial interval integral \(\int_{0}^{[t]}[x](\tau)d\tau\).

Note
From the monotonicity of the integral operator, \(\int_{0}^{[t]}[x](\tau)d\tau=[\int_{0}^{[t]}x^-(\tau)d\tau,\int_{0}^{[t]}x^+(\tau)d\tau]\)
Parameters
tinterval upper bound, subset of the Tube's tdomain
Returns
the pair \(\big([i^-],[i^+]\big)\), where \([i^-]=\int_{0}^{[t]}x^-(\tau)d\tau\) and \([i^+]=\int_{0}^{[t]}x^+(\tau)d\tau\)

◆ partial_integral() [2/2]

const std::pair< Interval, Interval > codac::Tube::partial_integral ( const Interval & t1,
const Interval & t2 ) const

Computes the partial interval integral \(\int_{[t_1]}^{[t_2]}[x](\tau)d\tau\).

Note
From the monotonicity of the integral operator, \(\int_{[t_1]}^{[t_2]}[x](\tau)d\tau=[\int_{[t_1]}^{[t_2]}x^-(\tau)d\tau,\int_{[t_1]}^{[t_2]}x^+(\tau)d\tau]\)
Parameters
t1interval lower bound, subset of the Tube's tdomain
t2interval upper bound, subset of the Tube's tdomain
Returns
the pair \(\big([i^-],[i^+]\big)\), where \([i^-]=\int_{[t_1]}^{[t_2]}x^-(\tau)d\tau\) and \([i^+]=\int_{[t_1]}^{[t_2]}x^+(\tau)d\tau\)

◆ serialize() [1/2]

void codac::Tube::serialize ( const std::string & binary_file_name = "x.tube",
int version_number = SERIALIZATION_VERSION ) const

Serializes this tube.

Note
The values and sampling (slices and gates) are serialized
Parameters
binary_file_namename of the output file (default value: "x.tube")
version_numberserialization version (used for tests purposes, default value: last version)

◆ serialize() [2/2]

void codac::Tube::serialize ( const std::string & binary_file_name,
const Trajectory & traj,
int version_number = SERIALIZATION_VERSION ) const

Serializes this tube together with a Trajectory object.

Note
The values and sampling (slices and gates) are serialized
The serialization of a Trajectory defined from a TFunction object is not supported
The output file will appear in the executable current directory
Parameters
binary_file_namename of the output file (default value: "x.tube")
trajthe Trajectory object to serialize (for instance, actual but unknown values)
version_numberserialization version (used for tests purposes, default value: last version)

◆ same_slicing()

static bool codac::Tube::same_slicing ( const Tube & x1,
const Tube & x2 )
static

Tests whether the two Tube objects are sharing the same slicing.

Note
If true, it means the two tubes are defined with the same amount of slices and identical sampling
Parameters
x1the first Tube
x2the second Tube
Returns
true in case of same slicing

◆ enable_syntheses()

static void codac::Tube::enable_syntheses ( bool enable = true)
static

Enables the computation of a synthesis tree for any Tube object.

Note
A synthesis tree speeds up computations such as integrals or evaluations
Parameters
enableboolean
Returns
void

◆ hull()

static const Tube codac::Tube::hull ( const std::list< Tube > & l_tubes)
static

Computes the hull of several tubes.

Parameters
l_tubeslist of tubes
Returns
the tube enveloping the other ones

◆ codomain_box()

const IntervalVector codac::Tube::codomain_box ( ) const
protectedvirtual

Returns the box \([x]([t_0,t_f])\).

Note
Used for genericity purposes
Returns
the envelope of codomain values

Implements codac::DynamicalItem.

◆ deserialize()

void codac::Tube::deserialize ( const std::string & binary_file_name,
Trajectory *& traj )
protected

Restores a scalar tube from serialization, together with a Trajectory object.

Note
The Tube and the Trajectory must have been serialized beforehand by the appropriate method serialize()
Parameters
binary_file_namepath to the binary file
traja pointer to the Trajectory object to be instantiated

◆ create_synthesis_tree()

void codac::Tube::create_synthesis_tree ( ) const
protected

Creates the synthesis tree associated to the values of this tube.

Note
The synthesis tree speeds up computations such as integrals or evaluations

◆ create_polynomial_synthesis()

void codac::Tube::create_polynomial_synthesis ( double eps) const
protected

Creates the synthesis tree associated to the values of this tube.

Note
The synthesis tree speeds up computations such as integrals or evaluations
Parameters
epsprecision of the polynomial approximation

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