pub struct ExplicitRungeKutta<E, F, T: Real, Y: State<T>, D: CallBackData, const O: usize, const S: usize, const I: usize> {
pub h0: T,
pub rtol: T,
pub atol: T,
pub h_max: T,
pub h_min: T,
pub max_steps: usize,
pub max_rejects: usize,
pub safety_factor: T,
pub min_scale: T,
pub max_scale: T,
/* private fields */
}Expand description
Runge-Kutta solver that can handle:
- Fixed-step methods with cubic Hermite interpolation
- Adaptive step methods with embedded error estimation and cubic Hermite interpolation
- Advanced methods with dense output interpolation using Butcher tableau coefficients
§Type Parameters
E: Equation type (e.g., Ordinary, Delay, Stochastic)F: Family type (e.g., Adaptive, Fixed, DormandPrince, etc.)T: Real number type (f32, f64)Y: State vector typeD: Callback data typeconst O: Order of the methodconst S: Number of stages in the methodconst I: Total number of stages including interpolation (equal to S for methods without dense output)
Fields§
§h0: T§rtol: T§atol: T§h_max: T§h_min: T§max_steps: usize§max_rejects: usize§safety_factor: T§min_scale: T§max_scale: TImplementations§
Source§impl<E, T: Real, Y: State<T>, D: CallBackData> ExplicitRungeKutta<E, Adaptive, T, Y, D, 5, 6, 6>
impl<E, T: Real, Y: State<T>, D: CallBackData> ExplicitRungeKutta<E, Adaptive, T, Y, D, 5, 6, 6>
Source§impl<E, T: Real, Y: State<T>, D: CallBackData> ExplicitRungeKutta<E, Adaptive, T, Y, D, 5, 6, 6>
impl<E, T: Real, Y: State<T>, D: CallBackData> ExplicitRungeKutta<E, Adaptive, T, Y, D, 5, 6, 6>
Source§impl<E, T: Real, Y: State<T>, D: CallBackData> ExplicitRungeKutta<E, Adaptive, T, Y, D, 6, 9, 10>
impl<E, T: Real, Y: State<T>, D: CallBackData> ExplicitRungeKutta<E, Adaptive, T, Y, D, 6, 9, 10>
Source§impl<E, T: Real, Y: State<T>, D: CallBackData> ExplicitRungeKutta<E, Adaptive, T, Y, D, 6, 9, 12>
impl<E, T: Real, Y: State<T>, D: CallBackData> ExplicitRungeKutta<E, Adaptive, T, Y, D, 6, 9, 12>
Source§impl<E, T: Real, Y: State<T>, D: CallBackData> ExplicitRungeKutta<E, Adaptive, T, Y, D, 7, 10, 13>
impl<E, T: Real, Y: State<T>, D: CallBackData> ExplicitRungeKutta<E, Adaptive, T, Y, D, 7, 10, 13>
Source§impl<E, T: Real, Y: State<T>, D: CallBackData> ExplicitRungeKutta<E, Adaptive, T, Y, D, 7, 10, 16>
impl<E, T: Real, Y: State<T>, D: CallBackData> ExplicitRungeKutta<E, Adaptive, T, Y, D, 7, 10, 16>
Source§impl<E, T: Real, Y: State<T>, D: CallBackData> ExplicitRungeKutta<E, Adaptive, T, Y, D, 8, 13, 17>
impl<E, T: Real, Y: State<T>, D: CallBackData> ExplicitRungeKutta<E, Adaptive, T, Y, D, 8, 13, 17>
Source§impl<E, T: Real, Y: State<T>, D: CallBackData> ExplicitRungeKutta<E, Adaptive, T, Y, D, 8, 13, 21>
impl<E, T: Real, Y: State<T>, D: CallBackData> ExplicitRungeKutta<E, Adaptive, T, Y, D, 8, 13, 21>
Source§impl<E, T: Real, Y: State<T>, D: CallBackData> ExplicitRungeKutta<E, Adaptive, T, Y, D, 9, 16, 21>
impl<E, T: Real, Y: State<T>, D: CallBackData> ExplicitRungeKutta<E, Adaptive, T, Y, D, 9, 16, 21>
Source§impl<E, T: Real, Y: State<T>, D: CallBackData> ExplicitRungeKutta<E, Adaptive, T, Y, D, 9, 16, 26>
impl<E, T: Real, Y: State<T>, D: CallBackData> ExplicitRungeKutta<E, Adaptive, T, Y, D, 9, 16, 26>
Source§impl<E, T: Real, Y: State<T>, D: CallBackData> ExplicitRungeKutta<E, DormandPrince, T, Y, D, 8, 12, 16>
impl<E, T: Real, Y: State<T>, D: CallBackData> ExplicitRungeKutta<E, DormandPrince, T, Y, D, 8, 12, 16>
Source§impl<E, T: Real, Y: State<T>, D: CallBackData> ExplicitRungeKutta<E, DormandPrince, T, Y, D, 5, 7, 7>
impl<E, T: Real, Y: State<T>, D: CallBackData> ExplicitRungeKutta<E, DormandPrince, T, Y, D, 5, 7, 7>
Source§impl<T: Real, Y: State<T>, D: CallBackData, const O: usize, const S: usize, const I: usize> ExplicitRungeKutta<Delay, Fixed, T, Y, D, O, S, I>
impl<T: Real, Y: State<T>, D: CallBackData, const O: usize, const S: usize, const I: usize> ExplicitRungeKutta<Delay, Fixed, T, Y, D, O, S, I>
Source§impl<E, T: Real, Y: State<T>, D: CallBackData> ExplicitRungeKutta<E, Fixed, T, Y, D, 1, 1, 1>
impl<E, T: Real, Y: State<T>, D: CallBackData> ExplicitRungeKutta<E, Fixed, T, Y, D, 1, 1, 1>
Source§impl<E, T: Real, Y: State<T>, D: CallBackData> ExplicitRungeKutta<E, Fixed, T, Y, D, 2, 2, 2>
impl<E, T: Real, Y: State<T>, D: CallBackData> ExplicitRungeKutta<E, Fixed, T, Y, D, 2, 2, 2>
Source§impl<E, T: Real, Y: State<T>, D: CallBackData> ExplicitRungeKutta<E, Fixed, T, Y, D, 2, 2, 2>
impl<E, T: Real, Y: State<T>, D: CallBackData> ExplicitRungeKutta<E, Fixed, T, Y, D, 2, 2, 2>
Source§impl<E, T: Real, Y: State<T>, D: CallBackData> ExplicitRungeKutta<E, Fixed, T, Y, D, 2, 2, 2>
impl<E, T: Real, Y: State<T>, D: CallBackData> ExplicitRungeKutta<E, Fixed, T, Y, D, 2, 2, 2>
Source§impl<E, T: Real, Y: State<T>, D: CallBackData> ExplicitRungeKutta<E, Fixed, T, Y, D, 4, 4, 4>
impl<E, T: Real, Y: State<T>, D: CallBackData> ExplicitRungeKutta<E, Fixed, T, Y, D, 4, 4, 4>
Source§impl<E, T: Real, Y: State<T>, D: CallBackData> ExplicitRungeKutta<E, Fixed, T, Y, D, 4, 4, 4>
impl<E, T: Real, Y: State<T>, D: CallBackData> ExplicitRungeKutta<E, Fixed, T, Y, D, 4, 4, 4>
Sourcepub fn three_eighths(h0: T) -> Self
pub fn three_eighths(h0: T) -> Self
Creates the three-eighths rule 4th order Runge-Kutta method.
Source§impl<E, F, T: Real, Y: State<T>, D: CallBackData, const O: usize, const S: usize, const I: usize> ExplicitRungeKutta<E, F, T, Y, D, O, S, I>
impl<E, F, T: Real, Y: State<T>, D: CallBackData, const O: usize, const S: usize, const I: usize> ExplicitRungeKutta<E, F, T, Y, D, O, S, I>
Sourcepub fn max_rejects(self, max_rejects: usize) -> Self
pub fn max_rejects(self, max_rejects: usize) -> Self
Set the maximum number of consecutive rejected steps before declaring stiffness
Sourcepub fn safety_factor(self, safety_factor: T) -> Self
pub fn safety_factor(self, safety_factor: T) -> Self
Set the safety factor for step size control (default: 0.9)
Sourcepub fn min_scale(self, min_scale: T) -> Self
pub fn min_scale(self, min_scale: T) -> Self
Set the minimum scale factor for step size changes (default: 0.2)
Sourcepub fn max_scale(self, max_scale: T) -> Self
pub fn max_scale(self, max_scale: T) -> Self
Set the maximum scale factor for step size changes (default: 10.0)
Sourcepub fn dense_stages(&self) -> usize
pub fn dense_stages(&self) -> usize
Get the number of terms in the dense output interpolation polynomial
Trait Implementations§
Source§impl<E, F, T: Real, Y: State<T>, D: CallBackData, const O: usize, const S: usize, const I: usize> Default for ExplicitRungeKutta<E, F, T, Y, D, O, S, I>
impl<E, F, T: Real, Y: State<T>, D: CallBackData, const O: usize, const S: usize, const I: usize> Default for ExplicitRungeKutta<E, F, T, Y, D, O, S, I>
Source§impl<const L: usize, T: Real, Y: State<T>, H: Fn(T) -> Y, D: CallBackData, const O: usize, const S: usize, const I: usize> DelayNumericalMethod<L, T, Y, H, D> for ExplicitRungeKutta<Delay, Adaptive, T, Y, D, O, S, I>
impl<const L: usize, T: Real, Y: State<T>, H: Fn(T) -> Y, D: CallBackData, const O: usize, const S: usize, const I: usize> DelayNumericalMethod<L, T, Y, H, D> for ExplicitRungeKutta<Delay, Adaptive, T, Y, D, O, S, I>
Source§fn init<F>(
&mut self,
dde: &F,
t0: T,
tf: T,
y0: &Y,
phi: &H,
) -> Result<Evals, Error<T, Y>>where
F: DDE<L, T, Y, D>,
fn init<F>(
&mut self,
dde: &F,
t0: T,
tf: T,
y0: &Y,
phi: &H,
) -> Result<Evals, Error<T, Y>>where
F: DDE<L, T, Y, D>,
Initialize the solver before integration Read more
Source§fn step<F>(&mut self, dde: &F, phi: &H) -> Result<Evals, Error<T, Y>>where
F: DDE<L, T, Y, D>,
fn step<F>(&mut self, dde: &F, phi: &H) -> Result<Evals, Error<T, Y>>where
F: DDE<L, T, Y, D>,
Advance the solution by one step Read more
Source§fn set_status(&mut self, status: Status<T, Y, D>)
fn set_status(&mut self, status: Status<T, Y, D>)
Set solver status
Source§impl<const L: usize, T: Real, Y: State<T>, H: Fn(T) -> Y, D: CallBackData, const O: usize, const S: usize, const I: usize> DelayNumericalMethod<L, T, Y, H, D> for ExplicitRungeKutta<Delay, DormandPrince, T, Y, D, O, S, I>
impl<const L: usize, T: Real, Y: State<T>, H: Fn(T) -> Y, D: CallBackData, const O: usize, const S: usize, const I: usize> DelayNumericalMethod<L, T, Y, H, D> for ExplicitRungeKutta<Delay, DormandPrince, T, Y, D, O, S, I>
Source§fn init<F>(
&mut self,
dde: &F,
t0: T,
tf: T,
y0: &Y,
phi: &H,
) -> Result<Evals, Error<T, Y>>where
F: DDE<L, T, Y, D>,
fn init<F>(
&mut self,
dde: &F,
t0: T,
tf: T,
y0: &Y,
phi: &H,
) -> Result<Evals, Error<T, Y>>where
F: DDE<L, T, Y, D>,
Initialize the solver before integration Read more
Source§fn step<F>(&mut self, dde: &F, phi: &H) -> Result<Evals, Error<T, Y>>where
F: DDE<L, T, Y, D>,
fn step<F>(&mut self, dde: &F, phi: &H) -> Result<Evals, Error<T, Y>>where
F: DDE<L, T, Y, D>,
Advance the solution by one step Read more
Source§fn set_status(&mut self, status: Status<T, Y, D>)
fn set_status(&mut self, status: Status<T, Y, D>)
Set solver status
Source§impl<const L: usize, T: Real, Y: State<T>, H: Fn(T) -> Y, D: CallBackData, const O: usize, const S: usize, const I: usize> DelayNumericalMethod<L, T, Y, H, D> for ExplicitRungeKutta<Delay, Fixed, T, Y, D, O, S, I>
impl<const L: usize, T: Real, Y: State<T>, H: Fn(T) -> Y, D: CallBackData, const O: usize, const S: usize, const I: usize> DelayNumericalMethod<L, T, Y, H, D> for ExplicitRungeKutta<Delay, Fixed, T, Y, D, O, S, I>
Source§fn init<F>(
&mut self,
dde: &F,
t0: T,
tf: T,
y0: &Y,
phi: &H,
) -> Result<Evals, Error<T, Y>>where
F: DDE<L, T, Y, D>,
fn init<F>(
&mut self,
dde: &F,
t0: T,
tf: T,
y0: &Y,
phi: &H,
) -> Result<Evals, Error<T, Y>>where
F: DDE<L, T, Y, D>,
Initialize the solver before integration Read more
Source§fn step<F>(&mut self, dde: &F, phi: &H) -> Result<Evals, Error<T, Y>>where
F: DDE<L, T, Y, D>,
fn step<F>(&mut self, dde: &F, phi: &H) -> Result<Evals, Error<T, Y>>where
F: DDE<L, T, Y, D>,
Advance the solution by one step Read more
Source§fn set_status(&mut self, status: Status<T, Y, D>)
fn set_status(&mut self, status: Status<T, Y, D>)
Set solver status
Source§impl<T: Real, Y: State<T>, D: CallBackData, const O: usize, const S: usize, const I: usize> Interpolation<T, Y> for ExplicitRungeKutta<Delay, Adaptive, T, Y, D, O, S, I>
impl<T: Real, Y: State<T>, D: CallBackData, const O: usize, const S: usize, const I: usize> Interpolation<T, Y> for ExplicitRungeKutta<Delay, Adaptive, T, Y, D, O, S, I>
Source§fn interpolate(&mut self, t_interp: T) -> Result<Y, Error<T, Y>>
fn interpolate(&mut self, t_interp: T) -> Result<Y, Error<T, Y>>
Interpolates the solution at a given time t_interp.
Source§impl<T: Real, Y: State<T>, D: CallBackData, const O: usize, const S: usize, const I: usize> Interpolation<T, Y> for ExplicitRungeKutta<Delay, DormandPrince, T, Y, D, O, S, I>
impl<T: Real, Y: State<T>, D: CallBackData, const O: usize, const S: usize, const I: usize> Interpolation<T, Y> for ExplicitRungeKutta<Delay, DormandPrince, T, Y, D, O, S, I>
Source§impl<T: Real, Y: State<T>, D: CallBackData, const O: usize, const S: usize, const I: usize> Interpolation<T, Y> for ExplicitRungeKutta<Delay, Fixed, T, Y, D, O, S, I>
impl<T: Real, Y: State<T>, D: CallBackData, const O: usize, const S: usize, const I: usize> Interpolation<T, Y> for ExplicitRungeKutta<Delay, Fixed, T, Y, D, O, S, I>
Source§fn interpolate(&mut self, t_interp: T) -> Result<Y, Error<T, Y>>
fn interpolate(&mut self, t_interp: T) -> Result<Y, Error<T, Y>>
Interpolates the solution at time t_interp within the last accepted step.
Source§impl<T: Real, Y: State<T>, D: CallBackData, const O: usize, const S: usize, const I: usize> Interpolation<T, Y> for ExplicitRungeKutta<Ordinary, Adaptive, T, Y, D, O, S, I>
impl<T: Real, Y: State<T>, D: CallBackData, const O: usize, const S: usize, const I: usize> Interpolation<T, Y> for ExplicitRungeKutta<Ordinary, Adaptive, T, Y, D, O, S, I>
Source§impl<T: Real, Y: State<T>, D: CallBackData, const O: usize, const S: usize, const I: usize> Interpolation<T, Y> for ExplicitRungeKutta<Ordinary, DormandPrince, T, Y, D, O, S, I>
impl<T: Real, Y: State<T>, D: CallBackData, const O: usize, const S: usize, const I: usize> Interpolation<T, Y> for ExplicitRungeKutta<Ordinary, DormandPrince, T, Y, D, O, S, I>
Source§impl<T: Real, Y: State<T>, D: CallBackData, const O: usize, const S: usize, const I: usize> Interpolation<T, Y> for ExplicitRungeKutta<Ordinary, Fixed, T, Y, D, O, S, I>
impl<T: Real, Y: State<T>, D: CallBackData, const O: usize, const S: usize, const I: usize> Interpolation<T, Y> for ExplicitRungeKutta<Ordinary, Fixed, T, Y, D, O, S, I>
Source§impl<T: Real, Y: State<T>, D: CallBackData, const O: usize, const S: usize, const I: usize> Interpolation<T, Y> for ExplicitRungeKutta<Stochastic, Fixed, T, Y, D, O, S, I>
impl<T: Real, Y: State<T>, D: CallBackData, const O: usize, const S: usize, const I: usize> Interpolation<T, Y> for ExplicitRungeKutta<Stochastic, Fixed, T, Y, D, O, S, I>
Source§impl<T: Real, Y: State<T>, D: CallBackData, const O: usize, const S: usize, const I: usize> OrdinaryNumericalMethod<T, Y, D> for ExplicitRungeKutta<Ordinary, Adaptive, T, Y, D, O, S, I>
impl<T: Real, Y: State<T>, D: CallBackData, const O: usize, const S: usize, const I: usize> OrdinaryNumericalMethod<T, Y, D> for ExplicitRungeKutta<Ordinary, Adaptive, T, Y, D, O, S, I>
Source§fn init<F>(
&mut self,
ode: &F,
t0: T,
tf: T,
y0: &Y,
) -> Result<Evals, Error<T, Y>>where
F: ODE<T, Y, D>,
fn init<F>(
&mut self,
ode: &F,
t0: T,
tf: T,
y0: &Y,
) -> Result<Evals, Error<T, Y>>where
F: ODE<T, Y, D>,
Initialize the solver before integration Read more
Source§fn step<F>(&mut self, ode: &F) -> Result<Evals, Error<T, Y>>where
F: ODE<T, Y, D>,
fn step<F>(&mut self, ode: &F) -> Result<Evals, Error<T, Y>>where
F: ODE<T, Y, D>,
Advance the solution by one step Read more
Source§fn set_status(&mut self, status: Status<T, Y, D>)
fn set_status(&mut self, status: Status<T, Y, D>)
Set solver status
Source§impl<T: Real, Y: State<T>, D: CallBackData, const O: usize, const S: usize, const I: usize> OrdinaryNumericalMethod<T, Y, D> for ExplicitRungeKutta<Ordinary, DormandPrince, T, Y, D, O, S, I>
impl<T: Real, Y: State<T>, D: CallBackData, const O: usize, const S: usize, const I: usize> OrdinaryNumericalMethod<T, Y, D> for ExplicitRungeKutta<Ordinary, DormandPrince, T, Y, D, O, S, I>
Source§fn init<F>(
&mut self,
ode: &F,
t0: T,
tf: T,
y0: &Y,
) -> Result<Evals, Error<T, Y>>where
F: ODE<T, Y, D>,
fn init<F>(
&mut self,
ode: &F,
t0: T,
tf: T,
y0: &Y,
) -> Result<Evals, Error<T, Y>>where
F: ODE<T, Y, D>,
Initialize the solver before integration Read more
Source§fn step<F>(&mut self, ode: &F) -> Result<Evals, Error<T, Y>>where
F: ODE<T, Y, D>,
fn step<F>(&mut self, ode: &F) -> Result<Evals, Error<T, Y>>where
F: ODE<T, Y, D>,
Advance the solution by one step Read more
Source§fn set_status(&mut self, status: Status<T, Y, D>)
fn set_status(&mut self, status: Status<T, Y, D>)
Set solver status
Source§impl<T: Real, Y: State<T>, D: CallBackData, const O: usize, const S: usize, const I: usize> OrdinaryNumericalMethod<T, Y, D> for ExplicitRungeKutta<Ordinary, Fixed, T, Y, D, O, S, I>
impl<T: Real, Y: State<T>, D: CallBackData, const O: usize, const S: usize, const I: usize> OrdinaryNumericalMethod<T, Y, D> for ExplicitRungeKutta<Ordinary, Fixed, T, Y, D, O, S, I>
Source§fn init<F>(
&mut self,
ode: &F,
t0: T,
tf: T,
y0: &Y,
) -> Result<Evals, Error<T, Y>>where
F: ODE<T, Y, D>,
fn init<F>(
&mut self,
ode: &F,
t0: T,
tf: T,
y0: &Y,
) -> Result<Evals, Error<T, Y>>where
F: ODE<T, Y, D>,
Initialize the solver before integration Read more
Source§fn step<F>(&mut self, ode: &F) -> Result<Evals, Error<T, Y>>where
F: ODE<T, Y, D>,
fn step<F>(&mut self, ode: &F) -> Result<Evals, Error<T, Y>>where
F: ODE<T, Y, D>,
Advance the solution by one step Read more
Source§fn set_status(&mut self, status: Status<T, Y, D>)
fn set_status(&mut self, status: Status<T, Y, D>)
Set solver status
Source§impl<T: Real, Y: State<T>, D: CallBackData, const O: usize, const S: usize, const I: usize> StochasticNumericalMethod<T, Y, D> for ExplicitRungeKutta<Stochastic, Fixed, T, Y, D, O, S, I>
impl<T: Real, Y: State<T>, D: CallBackData, const O: usize, const S: usize, const I: usize> StochasticNumericalMethod<T, Y, D> for ExplicitRungeKutta<Stochastic, Fixed, T, Y, D, O, S, I>
Source§fn init<F>(
&mut self,
sde: &mut F,
t0: T,
tf: T,
y0: &Y,
) -> Result<Evals, Error<T, Y>>where
F: SDE<T, Y, D>,
fn init<F>(
&mut self,
sde: &mut F,
t0: T,
tf: T,
y0: &Y,
) -> Result<Evals, Error<T, Y>>where
F: SDE<T, Y, D>,
Initialize the solver before integration Read more
Source§fn step<F>(&mut self, sde: &mut F) -> Result<Evals, Error<T, Y>>where
F: SDE<T, Y, D>,
fn step<F>(&mut self, sde: &mut F) -> Result<Evals, Error<T, Y>>where
F: SDE<T, Y, D>,
Advance the solution by one step Read more
Source§fn set_status(&mut self, status: Status<T, Y, D>)
fn set_status(&mut self, status: Status<T, Y, D>)
Set solver status
Auto Trait Implementations§
impl<E, F, T, Y, D, const O: usize, const S: usize, const I: usize> Freeze for ExplicitRungeKutta<E, F, T, Y, D, O, S, I>
impl<E, F, T, Y, D, const O: usize, const S: usize, const I: usize> RefUnwindSafe for ExplicitRungeKutta<E, F, T, Y, D, O, S, I>
impl<E, F, T, Y, D, const O: usize, const S: usize, const I: usize> Send for ExplicitRungeKutta<E, F, T, Y, D, O, S, I>
impl<E, F, T, Y, D, const O: usize, const S: usize, const I: usize> Sync for ExplicitRungeKutta<E, F, T, Y, D, O, S, I>
impl<E, F, T, Y, D, const O: usize, const S: usize, const I: usize> Unpin for ExplicitRungeKutta<E, F, T, Y, D, O, S, I>
impl<E, F, T, Y, D, const O: usize, const S: usize, const I: usize> UnwindSafe for ExplicitRungeKutta<E, F, T, Y, D, O, S, I>
Blanket Implementations§
Source§impl<T> BorrowMut<T> for Twhere
T: ?Sized,
impl<T> BorrowMut<T> for Twhere
T: ?Sized,
Source§fn borrow_mut(&mut self) -> &mut T
fn borrow_mut(&mut self) -> &mut T
Mutably borrows from an owned value. Read more
Source§impl<SS, SP> SupersetOf<SS> for SPwhere
SS: SubsetOf<SP>,
impl<SS, SP> SupersetOf<SS> for SPwhere
SS: SubsetOf<SP>,
Source§fn to_subset(&self) -> Option<SS>
fn to_subset(&self) -> Option<SS>
The inverse inclusion map: attempts to construct
self from the equivalent element of its
superset. Read moreSource§fn is_in_subset(&self) -> bool
fn is_in_subset(&self) -> bool
Checks if
self is actually part of its subset T (and can be converted to it).Source§fn to_subset_unchecked(&self) -> SS
fn to_subset_unchecked(&self) -> SS
Use with care! Same as
self.to_subset but without any property checks. Always succeeds.Source§fn from_subset(element: &SS) -> SP
fn from_subset(element: &SS) -> SP
The inclusion map: converts
self to the equivalent element of its superset.