API reference
BSeries.jl API
BSeries.BSeries
— ModuleBSeries
A collection of functionality around B-series in Julia. See
- Philippe Chartier, Ernst Hairer, Gilles Vilmart (2010) Algebraic Structures of B-series. Foundations of Computational Mathematics DOI: 10.1007/s10208-010-9065-1
API Documentation
The API of BSeries.jl is documented in the online documentation. Information on each function is available in their docstrings.
BSeries.jl re-exports everything from RootedTrees.jl. However, if you rely on functionality from that package, you should also include it explicitly in your project dependencies to track breaking changes, since the version numbers of RootedTrees.jl and BSeries.jl are not necessarily synchronized.
The main API of BSeries.jl consists of the following components.
- B-series behave like
AbstractDict
s mappingRootedTree
s to coefficients. - The B-series of time integration methods such as Runge-Kutta methods can be constructed by the function
bseries
. - Vector space operations (addition/subtraction and multiplication by scalars) are available.
- The algebraic structures of the composition law and the substitution law are implemented via
compose
andsubstitute
. - Backward error analysis can be performed via
modified_equation
s andmodifying_integrator
s.
Please consult the documentation or the docstrings for further information.
Please note that B-series analysis is most conveniently applied to the autonomous form of ordinary differential equations (ODEs). Thus, BSeries.jl and RootedTrees.jl usually assume that time integration methods give the same result, independent of whether an ODE is written in an autonomous or a non-autonomous form. For Runge-Kutta methods, this means that the usual row-sum assumption is used.
Referencing
If you use BSeries.jl for your research, please cite it using the bibtex entry
@article{ketcheson2023computing,
title={Computing with {B}-series},
author={Ketcheson, David I and Ranocha, Hendrik},
journal={ACM Transactions on Mathematical Software},
volume={49},
number={2},
year={2023},
month={06},
doi={10.1145/3573384},
eprint={2111.11680},
eprinttype={arXiv},
eprintclass={math.NA}
}
In addition, you can also refer to BSeries.jl directly as
@misc{ranocha2021bseries,
title={{BSeries.jl}: {C}omputing with {B}-series in {J}ulia},
author={Ranocha, Hendrik and Ketcheson, David I},
year={2021},
month={09},
howpublished={\url{https://github.com/ranocha/BSeries.jl}},
doi={10.5281/zenodo.5534602}
}
License and contributing
This project is licensed under the MIT license (see License). Since it is an open-source project, we are very happy to accept contributions from the community. Please refer to Contributing for more details.
BSeries.AverageVectorFieldMethod
— TypeAverageVectorFieldMethod([T=Rational{Int}])
Construct a representation of the average vector field (AVF) method using coefficients of type T
. You can pass it as argument to bseries
to construct the corresponding B-series.
Examples
We can generate this as follows.
julia> series = bseries(AverageVectorFieldMethod(), 3)
TruncatedBSeries{RootedTree{Int64, Vector{Int64}}, Rational{Int64}} with 5 entries:
RootedTree{Int64}: Int64[] => 1
RootedTree{Int64}: [1] => 1
RootedTree{Int64}: [1, 2] => 1//2
RootedTree{Int64}: [1, 2, 3] => 1//4
RootedTree{Int64}: [1, 2, 2] => 1//3
References
The B-series of the average vector field (AVF) method is given by $b(.) = 1$ and $b([t_1, ..., t_n]) = b(t_1)...b(t_n) / (n + 1)$, see
- Elena Celledoni, Robert I. McLachlan, David I. McLaren, Brynjulf Owren, G. Reinout W. Quispel, and William M. Wright. "Energy-preserving Runge-Kutta methods." ESAIM: Mathematical Modelling and Numerical Analysis 43, no. 4 (2009): 645-649. DOI: 10.1051/m2an/2009020
BSeries.ContinuousStageRungeKuttaMethod
— TypeContinuousStageRungeKuttaMethod(M)
A struct that describes a continuous stage Runge-Kutta (CSRK) method. It can be constructed by passing the parameter matrix M
as described by Miyatake and Butcher (2016). You can compute the B-series of the method by bseries
.
References
- Yuto Miyatake and John C. Butcher. "A characterization of energy-preserving methods and the construction of parallel integrators for Hamiltonian systems." SIAM Journal on Numerical Analysis 54, no. 3 (2016): DOI: 10.1137/15M1020861
BSeries.ExactSolution
— TypeExactSolution{V}()
Lazy representation of the B-series of the exact solution of an ordinary differential equation using coefficients of type at least as representative as V
.
BSeries.ExactSolution
— MethodExactSolution(series_integrator)
A representation of the B-series of the exact solution of an ODE using the same type of coefficients as the B-series series_integrator
.
BSeries.MultirateInfinitesimalSplitMethod
— TypeMultirateInfinitesimalSplitMethod(A, D, G, c)
References
- Knoth, Oswald, and Joerg Wensch. "Generalized split-explicit Runge-Kutta methods for the compressible Euler equations". Monthly Weather Review 142, no. 5 (2014): 2067-2081. DOI: 10.1175/MWR-D-13-00068.1
This code is considered to be experimental at the moment and can change any time.
BSeries.TruncatedBSeries
— TypeTruncatedBSeries
A struct that can describe B-series of both numerical integration methods (where the coefficient of the empty tree is unity) and right-hand sides of an ordinary differential equation and perturbations thereof (where the coefficient of the empty tree is zero) up to a prescribed order
.
Generally, this kind of struct
should be constructed via bseries
or one of the other functions returning a B-series, e.g., modified_equation
or modifying_integrator
.
BSeries.bseries
— Functionbseries(f::Function, order, iterator_type=RootedTreeIterator)
Return a truncated B-series up to the specified order
with coefficients determined by f
. The type of rooted trees is determined by the iterator_type
, which can be RootedTreeIterator
or BicoloredRootedTreeIterator
. Calling f(t, series)
needs to return the coefficient of the rooted tree t
of the desired series in a type-stable manner. For the empty tree, f
is called as f(t, nothing)
. Otherwise, the series constructed so far is passed as second argument, allowing one to access values of lower-order trees.
The coefficients of the B-series returned by this method need to be multiplied by a power of the time step divided by the symmetry
of the rooted tree and multiplied by the corresponding elementary differential of the input vector field $f$. See also evaluate
.
Examples
The B-series of the average vector field (AVF) method is given by $b(.) = 1$ and $b([t_1, ..., t_n]) = b(t_1)...b(t_n) / (n + 1)$, see
- Elena Celledoni, Robert I. McLachlan, David I. McLaren, Brynjulf Owren, G. Reinout W. Quispel, and William M. Wright. "Energy-preserving Runge-Kutta methods." ESAIM: Mathematical Modelling and Numerical Analysis 43, no. 4 (2009): 645-649. DOI: 10.1051/m2an/2009020
We can generate this as follows.
julia> series = bseries(3) do t, series
if order(t) in (0, 1)
return 1 // 1
else
v = 1 // 1
n = 0
for subtree in SubtreeIterator(t)
v *= series[subtree]
n += 1
end
return v / (n + 1)
end
end
TruncatedBSeries{RootedTree{Int64, Vector{Int64}}, Rational{Int64}} with 5 entries:
RootedTree{Int64}: Int64[] => 1
RootedTree{Int64}: [1] => 1
RootedTree{Int64}: [1, 2] => 1//2
RootedTree{Int64}: [1, 2, 3] => 1//4
RootedTree{Int64}: [1, 2, 2] => 1//3
BSeries.bseries
— Methodbseries(ark::AdditiveRungeKuttaMethod, order)
Compute the B-series of the additive Runge-Kutta method ark
up to a prescribed integer order
.
The coefficients of the B-series returned by this method need to be multiplied by a power of the time step divided by the symmetry
of the colored rooted tree and multiplied by the corresponding elementary differential of the input vector fields $f^\nu$. See also evaluate
.
BSeries.bseries
— Methodbseries(csrk::ContinuousStageRungeKuttaMethod, order)
Compute the B-series of the ContinuousStageRungeKuttaMethod
csrk
up to the prescribed integer order
as described by Miyatake & Butcher (2016).
The coefficients of the B-series returned by this method need to be multiplied by a power of the time step divided by the symmetry
of the rooted tree and multiplied by the corresponding elementary differential of the input vector field $f$. See also evaluate
.
Examples
The AverageVectorFieldMethod
is given by the parameter matrix with single entry one.
julia> M = fill(1//1, 1, 1)
1×1 Matrix{Rational{Int64}}:
1
julia> series = bseries(ContinuousStageRungeKuttaMethod(M), 4)
TruncatedBSeries{RootedTree{Int64, Vector{Int64}}, Rational{Int64}} with 9 entries:
RootedTree{Int64}: Int64[] => 1
RootedTree{Int64}: [1] => 1
RootedTree{Int64}: [1, 2] => 1//2
RootedTree{Int64}: [1, 2, 3] => 1//4
RootedTree{Int64}: [1, 2, 2] => 1//3
RootedTree{Int64}: [1, 2, 3, 4] => 1//8
RootedTree{Int64}: [1, 2, 3, 3] => 1//6
RootedTree{Int64}: [1, 2, 3, 2] => 1//6
RootedTree{Int64}: [1, 2, 2, 2] => 1//4
julia> series - bseries(AverageVectorFieldMethod(), order(series))
TruncatedBSeries{RootedTree{Int64, Vector{Int64}}, Rational{Int64}} with 9 entries:
RootedTree{Int64}: Int64[] => 0
RootedTree{Int64}: [1] => 0
RootedTree{Int64}: [1, 2] => 0
RootedTree{Int64}: [1, 2, 3] => 0
RootedTree{Int64}: [1, 2, 2] => 0
RootedTree{Int64}: [1, 2, 3, 4] => 0
RootedTree{Int64}: [1, 2, 3, 3] => 0
RootedTree{Int64}: [1, 2, 3, 2] => 0
RootedTree{Int64}: [1, 2, 2, 2] => 0
References
- Yuto Miyatake and John C. Butcher. "A characterization of energy-preserving methods and the construction of parallel integrators for Hamiltonian systems." SIAM Journal on Numerical Analysis 54, no. 3 (2016): DOI: 10.1137/15M1020861
BSeries.bseries
— Methodbseries(mis::MultirateInfinitesimalSplitMethod, order)
Compute the B-series of the multirate infinitesimal split method mis
up to a prescribed integer order
.
The coefficients of the B-series returned by this method need to be multiplied by a power of the time step divided by the symmetry
of the colored rooted tree and multiplied by the corresponding elementary differential of the input vector fields $f^\nu$. See also evaluate
.
BSeries.bseries
— Methodbseries(ros::RosenbrockMethod, order)
Compute the B-series of the Rosenbrock method ros
up to a prescribed integer order
.
The coefficients of the B-series returned by this method need to be multiplied by a power of the time step divided by the symmetry
of the rooted tree and multiplied by the corresponding elementary differential of the input vector field $f$. See also evaluate
.
BSeries.bseries
— Methodbseries(rk::RungeKuttaMethod, order)
bseries(A::AbstractMatrix, b::AbstractVector, c::AbstractVector, order)
Compute the B-series of the Runge-Kutta method rk
with Butcher coefficients A, b, c
up to a prescribed integer order
.
The coefficients of the B-series returned by this method need to be multiplied by a power of the time step divided by the symmetry
of the rooted tree and multiplied by the corresponding elementary differential of the input vector field $f$. See also evaluate
.
BSeries.bseries
— Methodbseries(avf::AverageVectorFieldMethod, order)
Compute the B-series of the AverageVectorFieldMethod
up to a prescribed integer order
.
The coefficients of the B-series returned by bseries
need to be multiplied by a power of the time step divided by the symmetry
of the rooted tree and multiplied by the corresponding elementary differential of the input vector field $f$. See also evaluate
.
BSeries.compose
— Methodcompose(b, a, t::RootedTree)
Compute the coefficient corresponding to the tree t
of the B-series that is formed by composing the B-series a
with the B-series b
. It is assumed that the B-series b
has the coefficient unity of the empty tree.
References
Section 3.1 of
- Philippe Chartier, Ernst Hairer, Gilles Vilmart (2010) Algebraic Structures of B-series. Foundations of Computational Mathematics DOI: 10.1007/s10208-010-9065-1
BSeries.compose
— Methodcompose(b, a; normalize_stepsize=false)
Compose the B-series a
with the B-series b
. It is assumed that the B-series b
has the coefficient unity of the empty tree.
In the notation of Chartier, Hairer and Vilmart (2010), we have compose(b, a) = b ⋅ a
. Note that this means that method b
is applied first, followed by method a
.
If normalize_stepsize = true
, the coefficients of the returned B-series are divided by 2^order(t)
for each rooted tree t
. This normalizes the step size so that the resulting numerical integrator B-series uses the same step size as the input series (instead of a doubled step size).
References
Section 3.1 of
- Philippe Chartier, Ernst Hairer, Gilles Vilmart (2010) Algebraic Structures of B-series. Foundations of Computational Mathematics DOI: 10.1007/s10208-010-9065-1
BSeries.compose
— Methodcompose(b1, b2, bs...; normalize_stepsize=false)
Compose the B-series b1
, b2
, bs...
. It is assumed that all B-series have the coefficient unity of the empty tree.
In the notation of Chartier, Hairer and Vilmart (2010), we have compose(b1, b2, b3) = b1 ⋅ b2 ⋅ b3
. Note that this product is associative and has to be read from left to right, i.e., method b1
is applied first, followed by b2
, bs...
.
If normalize_stepsize = true
, the coefficients of the returned B-series are divided by n^order(t)
for each rooted tree t
, where n
is the total number of composed B-series. This normalizes the step size so that the resulting numerical integrator B-series uses the same step size as the input series (instead of an n
-fold step size).
References
Section 3.1 of
- Philippe Chartier, Ernst Hairer, Gilles Vilmart (2010) Algebraic Structures of B-series. Foundations of Computational Mathematics DOI: 10.1007/s10208-010-9065-1
BSeries.compute_derivative
— Functioncompute_derivative(expression, variable)
Internal function specialized on symbolic variables and expressions from
if these packages are loaded (via Requires.jl or weak dependencies on Julia v1.9 and newer).
BSeries.elementary_differentials
— Methodelementary_differentials(f::AbstractVector, u, order)
Compute all elementary differentials of the vector field f
with independent variables u
up to the given order
. The return value can be indexed by rooted trees to obtain the corresponding elementary differential.
BSeries.elementary_differentials
— Methodelementary_differentials(fs::NTuple{2, AbstractVector}, u, order)
Compute all elementary differentials of the sum of the two vector fields f
with independent variables u
up to the given order
. The return value can be indexed by (bi-) colored rooted trees to obtain the corresponding elementary differential.
BSeries.energy_preserving_order
— Methodenergy_preserving_order(rk::RungeKuttaMethod, max_order)
This function checks up to which order a Runge-Kutta method rk
is energy-preserving for Hamiltonian problems. It requires a max_order
so that it does not run forever if the order up to which the method is energy-preserving is too big or infinite.
See also is_energy_preserving
.
BSeries.evaluate
— Functionevaluate(f, u, dt, series, reduce_order_by=0)
Evaluate the B-series series
specialized to the ordinary differential equation $u'(t) = f(u(t))$ with vector field f
and dependent variables u
for a time step size dt
.
Here, u
is assumed to be a vector of symbolic variables and f
is assumed to be a vector of expressions in these variables for plain B-series. For B-series with colored trees, f
must be a tuple of vectors of expressions in the variables u
. Currently, symbolic variables from
are supported.
The powers of dt
can be controlled by reduce_order_by
to make them different from the usual order(t)
for a rooted tree t
. This can be useful in the context of modified_equation
s or modifying_integrator
s, where the B-series coefficients are those of $h fₕ$, i.e., they contain an additional power of dt
. In this case, the B-series of the vector field can be obtained using reduce_order_by = 1
.
References
Section 3.2 of
- Philippe Chartier, Ernst Hairer, Gilles Vilmart (2010) Algebraic Structures of B-series. Foundations of Computational Mathematics DOI: 10.1007/s10208-010-9065-1
BSeries.is_energy_preserving
— Methodis_energy_preserving(series_integrator)::Bool
This function checks whether the B-series series_integrator
of a time integration method is energy-preserving for Hamiltonian systems - up to the order
of series_integrator
.
References
This code is based on the Theorem 2 of
- Elena Celledoni, Robert I. McLachlan, Brynjulf Owren, and G. R. W. Quispel. "Energy-preserving integrators and the structure of B-series." Foundations of Computational Mathematics 10 (2010): 673-693. DOI: 10.1007/s10208-010-9073-1
BSeries.is_energy_preserving
— Methodis_energy_preserving(rk::RungeKuttaMethod, order)::Bool
This function checks whether the Runge-Kutta method rk
is energy-preserving for Hamiltonian systems up to a given order
.
BSeries.is_symplectic
— Methodis_symplectic(series_integrator; kwargs...)::Bool
This function checks whether the B-series series_integrator
of a time integration method is symplectic (conserves quadratic invariants) - up to the order
of series_integrator
.
By default, the comparison of the coefficients entering the conditions is performed using isequal
. If keyword arguments such as absolute/relative tolerances atol
/rtol
are given or floating point numbers are used, the comparison is performed using isapprox
and the keyword arguments kwargs...
are forwarded.
See also order_of_symplecticity
.
BSeries.modified_equation
— Methodmodified_equation(f, u, dt, series_integrator)
Compute the B-series of the modified_equation
of the time integration method with B-series series_integrator
with respect to the ordinary differential equation $u'(t) = f(u(t))$ with vector field f
and dependent variables u
for a time step size dt
.
Here, u
is assumed to be a vector of symbolic variables and f
is assumed to be a vector of expressions in these variables for plain B-series. For B-series with colored trees, f
must be a tuple of vectors of expressions in the variables u
. Currently, symbolic variables from
are supported.
BSeries.modified_equation
— Methodmodified_equation(f, u, dt, rk::RungeKuttaMethod, order)
modified_equation(f, u, dt,
A::AbstractMatrix, b::AbstractVector, c::AbstractVector,
order)
Compute the B-series of the modified_equation
of the Runge-Kutta method rk
with Butcher coefficients A, b, c
up to the prescribed order
with respect to the ordinary differential equation $u'(t) = f(u(t))$ with vector field f
and dependent variables u
for a time step size dt
.
Here, u
is assumed to be a vector of symbolic variables and f
is assumed to be a vector of expressions in these variables for plain B-series. For B-series with colored trees, f
must be a tuple of vectors of expressions in the variables u
. Currently, symbolic variables from
are supported.
BSeries.modified_equation
— Methodmodified_equation(series_integrator)
Compute the B-series of the modified equation of the time integration method with B-series series_integrator
.
Given an ordinary differential equation (ODE) $u'(t) = f(u(t))$ and a Runge-Kutta method, the idea is to interpret the numerical solution with given time step size as exact solution of a modified ODE $u'(t) = fₕ(u(t))$. This method returns the B-series of $h fₕ$.
The coefficients of the B-series returned by this method need to be multiplied by a power of the time step divided by the symmetry
of the rooted tree and multiplied by the corresponding elementary differential of the input vector field $f$. See also evaluate
.
References
Section 3.2 of
- Philippe Chartier, Ernst Hairer, Gilles Vilmart (2010) Algebraic Structures of B-series. Foundations of Computational Mathematics DOI: 10.1007/s10208-010-9065-1
BSeries.modified_equation
— Methodmodified_equation(rk::RungeKuttaMethod, order)
modified_equation(A::AbstractMatrix, b::AbstractVector, c::AbstractVector,
order)
Compute the B-series of the modified_equation
of the Runge-Kutta method rk
with Butcher coefficients A, b, c
up to the prescribed order
.
The coefficients of the B-series returned by this method need to be multiplied by a power of the time step divided by the symmetry
of the rooted tree and multiplied by the corresponding elementary differential of the input vector field $f$. See also evaluate
.
BSeries.modifying_integrator
— Methodmodifying_integrator(f, u, dt, series_integrator)
Compute the B-series of the modifying_integrator
equation of the time integration method with B-series series_integrator
with respect to the ordinary differential equation $u'(t) = f(u(t))$ with vector field f
and dependent variables u
for a time step size dt
.
Here, u
is assumed to be a vector of symbolic variables and f
is assumed to be a vector of expressions in these variables for plain B-series. For B-series with colored trees, f
must be a tuple of vectors of expressions in the variables u
. Currently, symbolic variables from
are supported.
BSeries.modifying_integrator
— Methodmodifying_integrator(f, u, dt, rk::RungeKuttaMethod, order)
modifying_integrator(f, u, dt,
A::AbstractMatrix, b::AbstractVector, c::AbstractVector,
order)
Compute the B-series of the modifying_integrator
equation of the Runge-Kutta method with Butcher coefficients A, b, c
up to the prescribed order
with respect to the ordinary differential equation $u'(t) = f(u(t))$ with vector field f
and dependent variables u
for a time step size dt
.
Here, u
is assumed to be a vector of symbolic variables and f
is assumed to be a vector of expressions in these variables for plain B-series. For B-series with colored trees, f
must be a tuple of vectors of expressions in the variables u
. Currently, symbolic variables from
are supported.
BSeries.modifying_integrator
— Methodmodifying_integrator(series_integrator)
Compute the B-series of a "modifying integrator" equation of the time integration method with B-series series_integrator
.
Given an ordinary differential equation (ODE) $u'(t) = f(u(t))$ and a Runge-Kutta method, the idea is to find a modified ODE $u'(t) = fₕ(u(t))$ such that the numerical solution with given time step size is the exact solution of the original ODE. This method returns the B-series of $h fₕ$.
The coefficients of the B-series returned by this method need to be multiplied by a power of the time step divided by the symmetry
of the rooted tree and multiplied by the corresponding elementary differential of the input vector field $f$. See also evaluate
.
References
Section 3.2 of
- Philippe Chartier, Ernst Hairer, Gilles Vilmart (2010) Algebraic Structures of B-series. Foundations of Computational Mathematics DOI: 10.1007/s10208-010-9065-1
BSeries.modifying_integrator
— Methodmodifying_integrator(rk::RungeKuttaMethod, order)
modifying_integrator(A::AbstractMatrix, b::AbstractVector, c::AbstractVector,
order)
Compute the B-series of the modifying_integrator
equation of the Runge-Kutta method with Butcher coefficients A, b, c
up to the prescribed order
.
The coefficients of the B-series returned by this method need to be multiplied by a power of the time step divided by the symmetry
of the rooted tree and multiplied by the corresponding elementary differential of the input vector field $f$. See also evaluate
.
BSeries.order_of_accuracy
— Methodorder_of_accuracy(series; kwargs...)
Determine the order of accuracy of the B-series series
. By default, the comparison with the coefficients of the exact solution is performed using isequal
. If keyword arguments such as absolute/relative tolerances atol
/rtol
are given or floating point numbers are used, the comparison is performed using isapprox
and the keyword arguments kwargs...
are forwarded.
See also order
, ExactSolution
, order_of_symplecticity
.
BSeries.order_of_symplecticity
— Methodorder_of_symplecticity(series_integrator; kwargs...)
Determine the order of symplecticity of the B-series series_integrator
, i.e., the order up to which quadratic invariants are conserved. By default, the comparison of the coefficients entering the conditions is performed using isequal
. If keyword arguments such as absolute/relative tolerances atol
/rtol
are given or floating point numbers are used, the comparison is performed using isapprox
and the keyword arguments kwargs...
are forwarded.
See also is_symplectic
, order
, order_of_accuracy
.
BSeries.renormalize!
— Methodrenormalize!(series)
This function modifies a B-series by dividing each coefficient by the symmetry
of the corresponding tree.
This breaks assumptions made on the representation of a B-series. The modified B-series should not be passed to any other function assuming the default normalization.
BSeries.satisfied_for_trees_up_to_order
— Functionsatisfied_for_trees_of_order(condition, series, order,
iterator = RootedTreeIterator)
Checks whether a given condition
is satisfied for all pairs of trees t1
and t2
with given order == order(t1) + order(t2)
for a given series
. Which trees are considered is determined by the iterator
.
The condition
is called as condition(series, t1, t2)
and should return true
if the condition is satisfied and false
otherwise.
BSeries.substitute
— Methodsubstitute(b, a, t::AbstractRootedTree)
Compute the coefficient corresponding to the tree t
of the B-series that is formed by substituting the B-series b
into the B-series a
. It is assumed that the B-series b
has the coefficient zero of the empty tree.
References
Section 3.2 of
- Philippe Chartier, Ernst Hairer, Gilles Vilmart (2010) Algebraic Structures of B-series. Foundations of Computational Mathematics DOI: 10.1007/s10208-010-9065-1
BSeries.substitute
— Methodsubstitute(b, a)
Substitute the B-series b
into the B-series a
. It is assumed that the B-series b
has the coefficient zero of the empty tree.
In the notation of Chartier, Hairer and Vilmart (2010), we have substitute(b, a) = b ★ a
.
References
Section 3.2 of
- Philippe Chartier, Ernst Hairer, Gilles Vilmart (2010) Algebraic Structures of B-series. Foundations of Computational Mathematics DOI: 10.1007/s10208-010-9065-1
RootedTrees.elementary_weight
— Methodelementary_weight(t::RootedTree, csrk::ContinuousStageRungeKuttaMethod)
Compute the elementary weight Φ(t
) of the ContinuousStageRungeKuttaMethod
csrk
for a rooted tree t
`.
References
- Yuto Miyatake and John C. Butcher. "A characterization of energy-preserving methods and the construction of parallel integrators for Hamiltonian systems." SIAM Journal on Numerical Analysis 54, no. 3 (2016): DOI: 10.1137/15M1020861
RootedTrees.order
— Methodorder(series::TruncatedBSeries)
The maximal order
of a rooted tree with non-vanishing coefficient in the truncated B-series series
.
See also order_of_accuracy
.
RootedTrees.jl API
RootedTrees.RootedTrees
— ModuleRootedTrees
A collection of functionality around rooted trees to generate order conditions for Runge-Kutta methods in Julia. This package also provides basic functionality for BSeries.jl.
API Documentation
The API of RootedTrees.jl is documented in the following. Additional information on each function is available in their docstrings and in the online documentation.
Construction
RootedTree
s are represented using level sequences, i.e., AbstractVector
s containing the distances of the nodes from the root, see
- Beyer, Terry, and Sandra Mitchell Hedetniemi. "Constant time generation of rooted trees". SIAM Journal on Computing 9.4 (1980): 706-712. DOI: 10.1137/0209055
RootedTree
s can be constructed from their level sequence using
julia> t = rootedtree([1, 2, 3, 2])
RootedTree{Int64}: [1, 2, 3, 2]
In the notation of Butcher (Numerical Methods for ODEs, 2016), this tree can be written as [[τ] τ]
or (τ ∘ τ) ∘ (τ ∘ τ)
, where ∘
is the non-associative Butcher product of RootedTree
s, which is also implemented.
To get the representation of a RootedTree
introduced by Butcher, use butcher_representation
:
julia> t = rootedtree([1, 2, 3, 4, 3, 3, 2, 2, 2, 2, 2])
RootedTree{Int64}: [1, 2, 3, 4, 3, 3, 2, 2, 2, 2, 2]
julia> butcher_representation(t)
"[[[τ]τ²]τ⁵]"
There are also some simple plot recipes for Plots.jl. Thus, you can visualize a rooted tree t
using plot(t)
when using Plots
.
Additionally, there is an un-exported function RootedTrees.latexify
that can generate LaTeX code for a rooted tree t
based on the LaTeX package forest. The relevant code that needs to be included in the preamble can be obtained from the docstring of RootedTrees.latexify
(type ?
and RootedTrees.latexify
in the Julia REPL). The same format is used when you are using Latexify
and their function latexify
, see Latexify.jl.
Iteration over RootedTree
s
A RootedTreeIterator(order::Integer)
can be used to iterate efficiently over all RootedTree
s of a given order
.
Be careful that the iterator is stateful for efficiency reasons, so you might need to use copy
appropriately, e.g.,
julia> map(identity, RootedTreeIterator(4))
4-element Array{RootedTrees.RootedTree{Int64,Array{Int64,1}},1}:
RootedTree{Int64}: [1, 2, 2, 2]
RootedTree{Int64}: [1, 2, 2, 2]
RootedTree{Int64}: [1, 2, 2, 2]
RootedTree{Int64}: [1, 2, 2, 2]
julia> map(copy, RootedTreeIterator(4))
4-element Array{RootedTrees.RootedTree{Int64,Array{Int64,1}},1}:
RootedTree{Int64}: [1, 2, 3, 4]
RootedTree{Int64}: [1, 2, 3, 3]
RootedTree{Int64}: [1, 2, 3, 2]
RootedTree{Int64}: [1, 2, 2, 2]
Functions on Trees
The usual functions on RootedTree
s are implemented, cf. Butcher (Numerical Methods for ODEs, 2016).
order(t::RootedTree)
: The order of aRootedTree
, i.e., the length of its level sequence.σ(t::RootedTree)
orsymmetry(t)
: The symmetryσ
of a rooted tree, i.e., the order of the group of automorphisms on a particular labelling (of the vertices) oft
.γ(t::RootedTree)
ordensity(t)
: The densityγ(t)
of a rooted tree, i.e., the product over all vertices oft
of the order of the subtree rooted at that vertex.α(t::RootedTree)
: The number of monotonic labelings oft
not equivalent under the symmetry group.β(t::RootedTree)
: The total number of labelings oft
not equivalent under the symmetry group.
Additionally, functions on trees connected to Runge-Kutta methods are implemented.
elementary_weight(t, A, b, c)
: Compute the elementary weight Φ(t
) oft::RootedTree
for the Butcher coefficientsA, b, c
of a Runge-Kutta method.derivative_weight(t, A, b, c)
: Compute the derivative weight (ΦᵢD)(t
) oft
for the Butcher coefficientsA, b, c
of a Runge-Kutta method.residual_order_condition(t, A, b, c)
: The residual of the order condition(Φ(t) - 1/γ(t)) / σ(t)
with elementary weightΦ(t)
, densityγ(t)
, and symmetryσ(t)
of the rooted treet
for the Runge-Kutta method with Butcher coefficientsA, b, c
.
Brief Changelog
- v2.16: The LaTeX printing of rooted trees changed to allow representing colored rooted trees. Please update your LaTeX preamble as described in the docstring of
RootedTrees.latexify
. - v2.0: Rooted trees are considered up to isomorphisms introduced by shifting each coefficient of their level sequence by the same number.
Referencing
If you use RootedTrees.jl for your research, please cite the paper
@article{ketcheson2023computing,
title={Computing with {B}-series},
author={Ketcheson, David I and Ranocha, Hendrik},
journal={ACM Transactions on Mathematical Software},
volume={49},
number={2},
year={2023},
month={06},
doi={10.1145/3573384},
eprint={2111.11680},
eprinttype={arXiv},
eprintclass={math.NA}
}
In addition, you can also refer to RootedTrees.jl directly as
@misc{ranocha2019rootedtrees,
title={{RootedTrees.jl}: {A} collection of functionality around rooted trees
to generate order conditions for {R}unge-{K}utta methods in {J}ulia
for differential equations and scientific machine learning ({SciM}L)},
author={Ranocha, Hendrik and contributors},
year={2019},
month={05},
howpublished={\url{https://github.com/SciML/RootedTrees.jl}},
doi={10.5281/zenodo.5534590}
}
RootedTrees.AdditiveRungeKuttaMethod
— TypeAdditiveRungeKuttaMethod(rks)
AdditiveRungeKuttaMethod(As, bs, cs=map(A -> vec(sum(A, dims=2)), As))
Represent an additive Runge-Kutta method with collections of Butcher coefficients As
, bs
, and cs
. Alternatively, you can pass a collection of RungeKuttaMethod
s to the constructor. If the cs
are not provided, the usual "row sum" requirement of consistency with autonomous problems is applied.
An additive Runge-Kutta method applied to the ODE problem
\[ u'(t) = \sum_\nu f^\nu(t, u(t))\]
has the form
\[\begin{aligned} y^i &= u^n + \Delta t \sum_\nu \sum_j a^\nu_{i,j} f^\nu(y^i), \\ u^{n+1} &= u^n + \Delta t \sum_\nu \sum_i b^\nu_{i} f^\nu(y^i). \end{aligned}\]
In particular, additive Runge-Kutta methods are a superset of partitioned RK methods, which are applied to partitioned problems of the form
\[ (u^1)'(t) = f^1(t, u^1, u^2), \quad (u^2)'(t) = f^2(t, u^1, u^2).\]
References
- A. L. Araujo, A. Murua, and J. M. Sanz-Serna. "Symplectic Methods Based on Decompositions". SIAM Journal on Numerical Analysis 34.5 (1997): 1926-1947. DOI: 10.1137/S0036142995292128
RootedTrees.BicoloredRootedTree
— TypeBicoloredRootedTree{T<:Integer}
Representation of bicolored rooted trees.
See also ColoredRootedTree
, RootedTree
, rootedtree
.
RootedTrees.BicoloredRootedTreeIterator
— TypeBicoloredRootedTreeIterator(order::Integer)
Iterator over all bi-colored rooted trees of given order
. The returned trees are views to an internal tree modified during the iteration. If the returned trees shall be stored or modified during the iteration, a copy
has to be made.
RootedTrees.ColoredRootedTree
— TypeColoredRootedTree(level_sequence, color_sequence, is_canonical::Bool=false)
Represents a colored rooted tree using its level sequence. The single-colored version is RootedTree
.
See also BicoloredRootedTree
, rootedtree
.
This is a low-overhead and unsafe constructor. Please consider calling rootedtree
instead.
References
- Terry Beyer and Sandra Mitchell Hedetniemi. "Constant time generation of rooted trees". SIAM Journal on Computing 9.4 (1980): 706-712. DOI: 10.1137/0209055
- A. L. Araujo, A. Murua, and J. M. Sanz-Serna. "Symplectic Methods Based on Decompositions". SIAM Journal on Numerical Analysis 34.5 (1997): 1926–1947. DOI: 10.1137/S0036142995292128
RootedTrees.PartitionForestIterator
— TypePartitionForestIterator(t::AbstractRootedTree, edge_set)
Lazy iterator representation of the partition_forest
of the rooted tree t
. Similar to RootedTreeIterator
, you should copy
the iterates if you want to store or modify them during the iteration since they may be views to internal caches.
See also partition_forest
, partition_skeleton
, and PartitionIterator
.
References
Section 2.3 of
- Philippe Chartier, Ernst Hairer, Gilles Vilmart (2010) Algebraic Structures of B-series. Foundations of Computational Mathematics DOI: 10.1007/s10208-010-9065-1
RootedTrees.PartitionIterator
— TypePartitionIterator(t::AbstractRootedTree)
Iterator over all partition forests and skeletons of the rooted tree t
. This is basically a pure iterator version of all_partitions
. In particular, the partition forest may only be realized as an iterator. Similar to RootedTreeIterator
, you should copy
the iterates if you want to store or modify them during the iteration since they may be views to internal caches.
See also partition_forest
, partition_skeleton
, and PartitionForestIterator
.
References
Section 2.3 of
- Philippe Chartier, Ernst Hairer, Gilles Vilmart (2010) Algebraic Structures of B-series. Foundations of Computational Mathematics DOI: 10.1007/s10208-010-9065-1
RootedTrees.RootedTree
— TypeRootedTree(level_sequence, is_canonical::Bool=false)
Represents a rooted tree using its level sequence.
This is a low-overhead and unsafe constructor. Please consider calling rootedtree
instead.
References
- Terry Beyer and Sandra Mitchell Hedetniemi. "Constant time generation of rooted trees". SIAM Journal on Computing 9.4 (1980): 706-712. DOI: 10.1137/0209055
RootedTrees.RootedTreeIterator
— TypeRootedTreeIterator(order::Integer)
Iterator over all rooted trees of given order
. The returned trees are views to an internal tree modified during the iteration. If the returned trees shall be stored or modified during the iteration, a copy
has to be made.
RootedTrees.RosenbrockMethod
— TypeRosenbrockMethod(γ, A, b, c=vec(sum(A, dims=2)))
Represent a Rosenbrock (or Rosenbrock-Wanner, ROW) method with coefficients γ
, A
, b
, and c
. If c
is not provided, the usual "row sum" requirement of consistency with autonomous problems is applied.
Reference
- Ernst Hairer, Gerhard Wanner. Solving ordinary differential equations II: Stiff and differential-algebraic problems. Springer, 2010. Section IV.7
RootedTrees.RungeKuttaMethod
— TypeRungeKuttaMethod(A, b, c=vec(sum(A, dims=2)))
Represent a Runge-Kutta method with Butcher coefficients A
, b
, and c
. If c
is not provided, the usual "row sum" requirement of consistency with autonomous problems is applied.
RootedTrees.SplittingIterator
— TypeSplittingIterator(t::RootedTree)
Iterator over all splitting forests and subtrees of the rooted tree t
. This is basically an iterator version of all_splittings
.
See also partition_forest
and partition_skeleton
.
References
Section 2.2 of
- Philippe Chartier, Ernst Hairer, Gilles Vilmart (2010) Algebraic Structures of B-series. Foundations of Computational Mathematics DOI: 10.1007/s10208-010-9065-1
RootedTrees.SubtreeIterator
— TypeSubtreeIterator(t::AbstractRootedTree)
Lazy iterator representation of the subtrees
of the rooted tree t
. Similar to RootedTreeIterator
, you should copy
the iterates if you want to store or modify them during the iteration since they may be views to internal caches.
Base.:==
— Method==(t1::ColoredRootedTree, t2::ColoredRootedTree)
Compares two rooted trees based on their level (first) and color (second) sequences while considering equivalence classes given by different root indices.
Base.:==
— Method==(t1::RootedTree, t2::RootedTree)
Compares two rooted trees based on their level sequences while considering equivalence classes given by different root indices.
Examples
julia> t1 = rootedtree([1, 2, 3]);
julia> t2 = rootedtree([2, 3, 4]);
julia> t3 = rootedtree([1, 2, 2]);
julia> t1 == t2
true
julia> t1 == t3
false
Base.:∘
— Methodt1 ∘ t2
The non-associative Butcher product of rooted trees. It is formed by adding an edge from the root of t1
to the root of t2
.
See also butcher_product!
.
Reference: Section 301 of
- Butcher, John Charles. Numerical methods for ordinary differential equations. John Wiley & Sons, 2016.
Base.isless
— Methodisless(t1::ColoredRootedTree, t2::ColoredRootedTree)
Compares two colored rooted trees using a lexicographical comparison of their level (first) and color (second) sequences while considering equivalence classes given by different root indices.
Base.isless
— Methodisless(t1::RootedTree, t2::RootedTree)
Compares two rooted trees using a lexicographical comparison of their level sequences while considering equivalence classes given by different root indices.
RootedTrees.all_partitions
— Methodall_partitions(t::RootedTree)
Create all partition forests and skeletons of a rooted tree t
. This returns vectors of the return values of partition_forest
and partition_skeleton
when looping over all possible edge sets.
See also PartitionIterator
.
References
Section 2.3 of
- Philippe Chartier, Ernst Hairer, Gilles Vilmart (2010) Algebraic Structures of B-series. Foundations of Computational Mathematics DOI: 10.1007/s10208-010-9065-1
RootedTrees.all_splittings
— Methodall_splittings(t::RootedTree)
Create all splitting forests and subtrees associated to ordered subtrees of a rooted tree t
.
See also SplittingIterator
.
References
Section 2.2 of
- Philippe Chartier, Ernst Hairer, Gilles Vilmart (2010) Algebraic Structures of B-series. Foundations of Computational Mathematics DOI: 10.1007/s10208-010-9065-1
RootedTrees.butcher_product!
— Methodbutcher_product!(t, t1, t2)
Compute the non-associative Butcher product t = t1 ∘ t2
of rooted trees in-place. It is formed by adding an edge from the root of t1
to the root of t2
.
See also ∘
(available as \circ
plus TAB).
Reference: Section 301 of
- Butcher, John Charles. Numerical methods for ordinary differential equations. John Wiley & Sons, 2016.
RootedTrees.butcher_representation
— Functionbutcher_representation(t::RootedTree)
Returns the representation of t::RootedTree
introduced by Butcher as a string. Thus, the rooted tree consisting whose only vertex is the root itself is represented as τ
. The representation of other trees is defined recursively; if t₁, t₂, ... tₙ
are the subtrees
of the rooted tree t
, it is represented as t = [t₁ t₂ ... tₙ]
. If multiple subtrees are the same, their number of occurrences is written as a power.
Examples
julia> rootedtree([1, 2, 3, 2]) |> butcher_representation
"[[τ]τ]"
julia> rootedtree([1, 2, 3, 3, 2]) |> butcher_representation
"[[τ²]τ]"
References
Section 300 of
- Butcher, John Charles. Numerical methods for ordinary differential equations. John Wiley & Sons, 2008.
RootedTrees.canonical_representation!
— Functioncanonical_representation!(t::AbstractRootedTree)
Change the representation of the rooted tree t
to the canonical one, i.e., the one with lexicographically biggest level sequence.
See also canonical_representation
.
RootedTrees.canonical_representation
— Methodcanonical_representation(t::AbstractRootedTree)
Returns a new tree using the canonical representation of the rooted tree t
, i.e., the one with lexicographically biggest level sequence.
See also canonical_representation!
.
RootedTrees.check_canonical
— Methodcheck_canonical(t::AbstractRootedTree)
Check whether t
is in canonical representation.
This function is considered to be an internal implementation detail and will not necessarily be stable.
RootedTrees.count_trees
— Methodcount_trees(order)
Counts all rooted trees of order
.
RootedTrees.density
— Methodγ(t::AbstractRootedTree)
density(t::AbstractRootedTree)
The density γ(t)
of a rooted tree, i.e., the product over all vertices of t
of the order of the subtree rooted at that vertex.
Reference: Section 301 of
- Butcher, John Charles. Numerical methods for ordinary differential equations. John Wiley & Sons, 2008.
RootedTrees.derivative_weight
— Methodderivative_weight(t::ColoredRootedTree, ark::AdditiveRungeKuttaMethod)
Compute the derivative weight (ΦᵢD)(t
) of the AdditiveRungeKuttaMethod
ark
for the colored rooted tree t
.
References
- A. L. Araujo, A. Murua, and J. M. Sanz-Serna. "Symplectic Methods Based on Decompositions". SIAM Journal on Numerical Analysis 34.5 (1997): 1926–1947. DOI: 10.1137/S0036142995292128
- Butcher, John Charles. Numerical methods for ordinary differential equations. John Wiley & Sons, 2008. Section 312
RootedTrees.derivative_weight
— Methodderivative_weight(t::RootedTree, ros::RosenbrockMethod)
Compute the derivative weight (ΦᵢD)(t
) of the RosenbrockMethod
ros
for the rooted tree t
.
RootedTrees.derivative_weight
— Methodderivative_weight(t::RootedTree, rk::RungeKuttaMethod)
Compute the derivative weight (ΦᵢD)(t
) of the RungeKuttaMethod
rk
with Butcher coefficients A, b, c
for the rooted tree t
.
Reference: Section 312 of
- Butcher, John Charles. Numerical methods for ordinary differential equations. John Wiley & Sons, 2008.
RootedTrees.elementary_differential_latexstring
— Methodelementary_differential_latexstring(t::RootedTree)
Returns the elementary differential as a LaTeXString
from the package LaTeXStrings.jl.
RootedTrees.elementary_weight
— Methodelementary_weight(t::ColoredRootedTree, ark::AdditiveRungeKuttaMethod)
Compute the elementary weight Φ(t
) of the AdditiveRungeKuttaMethod
ark
for a colored rooted tree t
.
References
- A. L. Araujo, A. Murua, and J. M. Sanz-Serna. "Symplectic Methods Based on Decompositions". SIAM Journal on Numerical Analysis 34.5 (1997): 1926–1947. DOI: 10.1137/S0036142995292128
- Butcher, John Charles. Numerical methods for ordinary differential equations. John Wiley & Sons, 2008. Section 312
RootedTrees.elementary_weight
— Methodelementary_weight(t::RootedTree, ros::RosenbrockMethod)
Compute the elementary weight Φ(t
) of the RosenbrockMethod
ros
for a rooted tree t
.
RootedTrees.elementary_weight
— Methodelementary_weight(t::RootedTree, rk::RungeKuttaMethod)
elementary_weight(t::RootedTree, A::AbstractMatrix, b::AbstractVector, c::AbstractVector)
Compute the elementary weight Φ(t
) of the RungeKuttaMethod
rk
with Butcher coefficients A, b, c
for a rooted tree t
`.
Reference: Section 312 of
- Butcher, John Charles. Numerical methods for ordinary differential equations. John Wiley & Sons, 2008.
RootedTrees.elementary_weight_latexstring
— Methodelementary_weight_latexstring(t::RootedTree)
Returns the elementary_weight as a LaTeXString
from the package LaTeXStrings.jl.
RootedTrees.latexify
— Methodlatexify(t::Union{RootedTree, BicoloredRootedTree})
Return a LaTeX representation of the rooted tree t
. This makes use of the LaTeX package forest and assumes that you use the following LaTeX code in the preamble.
% Classical and colored Butcher trees based on
% https://tex.stackexchange.com/a/673436
\usepackage{forest}
\forestset{
whitenode/.style={draw, circle, minimum size=0.5ex, inner sep=0pt},
blacknode/.style={draw, fill=black, circle, minimum size=0.5ex, inner sep=0pt},
colornode/.style={draw, fill=#1, circle, minimum size=0.5ex, inner sep=0pt},
colornode/.default={red}
}
\newcommand{\blankforrootedtree}{\rule{0pt}{0pt}}
\NewDocumentCommand\rootedtree{o}{\begin{forest}
for tree={grow'=90, thick, edge=thick, l sep=0.5ex, l=0pt, s sep=0.5ex},
delay={
where content={}{
for children={no edge, before drawing tree={for tree={y-=5pt}}}
}
{
where content={o}{content={\blankforrootedtree}, whitenode}{
where content={.}{content={\blankforrootedtree}, blacknode}{}
}
}
}
[#1]
\end{forest}}
To change the style of latexify
to a human-readable Butcher-representation, you can use RootedTrees.set_latexify_style
.
Examples
julia> rootedtree([1, 2, 2]) |> RootedTrees.latexify |> println
\rootedtree[.[.][.]]
julia> rootedtree([1, 2, 3, 3, 2]) |> RootedTrees.latexify |> println
\rootedtree[.[.[.][.]][.]]
RootedTrees.normalize_root!
— Functionnormalize_root!(t::AbstractRootedTree, root=one(eltype(t.level_sequence)))
Normalize the level sequence of the rooted tree t
such that the root is set to root
.
RootedTrees.order
— Methodorder(t::AbstractRootedTree)
The order
of a rooted tree t
, i.e., the length of its level sequence.
RootedTrees.partition_forest
— Methodpartition_forest(t::RootedTree, edge_set)
Form the partition forest of the rooted tree t
where edges marked with false
in the edge_set
are removed. The ith value in the Boolean iterable edge_set
corresponds to the edge connecting node i+1
in the level sequence to its parent.
See also partition_skeleton
, PartitionIterator
, and PartitionForestIterator
.
References
Section 2.3 of
- Philippe Chartier, Ernst Hairer, Gilles Vilmart (2010) Algebraic Structures of B-series. Foundations of Computational Mathematics DOI: 10.1007/s10208-010-9065-1
RootedTrees.partition_skeleton
— Methodpartition_skeleton(t::AbstractRootedTree, edge_set)
Form the partition skeleton of the rooted tree t
, i.e., the rooted tree obtained by contracting each tree of the partition forest to a single vertex and re-establishing the edges removed to obtain the partition forest.
See also partition_forest
and PartitionIterator
.
References
Section 2.3 (and Section 6.1 for colored trees) of
- Philippe Chartier, Ernst Hairer, Gilles Vilmart (2010) Algebraic Structures of B-series. Foundations of Computational Mathematics DOI: 10.1007/s10208-010-9065-1
RootedTrees.residual_order_condition
— Methodresidual_order_condition(t::ColoredRootedTree, ark::AdditiveRungeKuttaMethod)
The residual of the order condition (Φ(t) - 1/γ(t)) / σ(t)
with elementary_weight
Φ(t)
, density
γ(t)
, and symmetry
σ(t)
of the AdditiveRungeKuttaMethod
ark
for the colored rooted tree t
.
References
- A. L. Araujo, A. Murua, and J. M. Sanz-Serna. "Symplectic Methods Based on Decompositions". SIAM Journal on Numerical Analysis 34.5 (1997): 1926–1947. DOI: 10.1137/S0036142995292128
- Butcher, John Charles. Numerical methods for ordinary differential equations. John Wiley & Sons, 2008. Section 312
RootedTrees.residual_order_condition
— Methodresidual_order_condition(t::RootedTree, ros::RosenbrockMethod)
The residual of the order condition (Φ(t) - 1/γ(t)) / σ(t)
with elementary_weight
Φ(t)
, density
γ(t)
, and symmetry
σ(t)
of the RosenbrockMethod
ros
for the rooted tree t
.
Reference
- Ernst Hairer, Gerhard Wanner. Solving ordinary differential equations II: Stiff and differential-algebraic problems. Springer, 2010. Section IV.7
RootedTrees.residual_order_condition
— Methodresidual_order_condition(t::RootedTree, rk::RungeKuttaMethod)
The residual of the order condition (Φ(t) - 1/γ(t)) / σ(t)
with elementary_weight
Φ(t)
, density
γ(t)
, and symmetry
σ(t)
of the RungeKuttaMethod
rk
with Butcher coefficients A, b, c
for the rooted tree t
.
Reference: Section 315 of
- Butcher, John Charles. Numerical methods for ordinary differential equations. John Wiley & Sons, 2008.
RootedTrees.root_color
— Methodroot_color(t::ColoredRootedTree)
Return the color of the root of t
.
RootedTrees.rootedtree!
— Methodrootedtree!(level_sequence, color_sequence)
Construct a canonical ColoredRootedTree
object from a level_sequence
and a color_sequence
which may be modified in this process. See also rootedtree
.
References
- Terry Beyer and Sandra Mitchell Hedetniemi. "Constant time generation of rooted trees". SIAM Journal on Computing 9.4 (1980): 706-712. DOI: 10.1137/0209055
RootedTrees.rootedtree!
— Methodrootedtree!(level_sequence)
Construct a canonical RootedTree
object from a level_sequence
which may be modified in this process. See also rootedtree
.
This may modify the level_sequence
and further modifications of the level_sequence
may invalidate the rooted tree returned by this function. Please consider calling rootedtree
instead.
References
- Terry Beyer and Sandra Mitchell Hedetniemi. "Constant time generation of rooted trees". SIAM Journal on Computing 9.4 (1980): 706-712. DOI: 10.1137/0209055
RootedTrees.rootedtree
— Methodrootedtree(level_sequence, color_sequence)
Construct a canonical ColoredRootedTree
object from a level_sequence
and a color_sequence
, i.e., a vector of integers representing the levels of each node of the tree and a vector of associated colors (e.g., Bool
s or Integers
).
References
- Terry Beyer and Sandra Mitchell Hedetniemi. "Constant time generation of rooted trees". SIAM Journal on Computing 9.4 (1980): 706-712. DOI: 10.1137/0209055
RootedTrees.rootedtree
— Methodrootedtree(level_sequence)
Construct a canonical RootedTree
object from a level_sequence
, i.e., a vector of integers representing the levels of each node of the tree.
References
- Terry Beyer and Sandra Mitchell Hedetniemi. "Constant time generation of rooted trees". SIAM Journal on Computing 9.4 (1980): 706-712. DOI: 10.1137/0209055
RootedTrees.set_latexify_style
— MethodRootedTrees.set_latexify_style(style::String)
Set the style of rooted trees when using latexify
. Possible options are
- "butcher": print the
butcher_representation
of rooted trees - "forest": use the LaTeX macro
\rootedtree
described in the docstring oflatexify
This system is based on Preferences.jl.
RootedTrees.set_printing_style
— MethodRootedTrees.set_printing_style(style::String)
Set the printing style of rooted trees. Possible options are
- "butcher": print the
butcher_representation
of rooted trees - "sequence": print the level sequence representation
This system is based on Preferences.jl.
RootedTrees.subtrees
— Methodsubtrees(t::ColoredRootedTree)
Returns a vector of all subtrees of t
.
RootedTrees.subtrees
— MethodRootedTrees.symmetry
— Methodσ(t::AbstractRootedTree)
symmetry(t::AbstractRootedTree)
The symmetry σ
of a rooted tree t
, i.e., the order of the group of automorphisms on a particular labelling (of the vertices) of t
.
Reference: Section 301 of
- Butcher, John Charles. Numerical methods for ordinary differential equations. John Wiley & Sons, 2008.
RootedTrees.unsafe_copyto!
— Methodunsafe_copyto!(t_dst::AbstractRootedTree, dst_offset,
t_src::AbstractRootedTree, src_offset, N)
Copy N
nodes from t_src
starting at offset src_offset
to t_dst
starting at offset dst_offset
. The types of the rooted trees must match. For example, you cannot copy a ColoredRootedTree
to a RootedTree
.
This is an unsafe operation since the rooted tree t_dst
will not necessarily be in canonical representation afterwards, even if the corresponding flag of t_dst
is set. Use with caution!
This function is considered to be an internal implementation detail and will not necessarily be stable.
RootedTrees.unsafe_deleteat!
— Methodunsafe_deleteat!(t::AbstractRootedTree, i)
Delete the node i
from the rooted tree t
. This is an unsafe operation since the rooted tree will not necessarily be in canonical representation afterwards, even if the corresponding flag of t
is set. Use with caution!
This function is considered to be an internal implementation detail and will not necessarily be stable.
RootedTrees.unsafe_resize!
— Methodunsafe_resize!(t::AbstractRootedTree, n::Integer)
Resize the rooted tree t
to n
nodes. This is an unsafe operation since the rooted tree will not necessarily be in canonical representation afterwards, even if the corresponding flag of t
is set. Use with caution!
This function is considered to be an internal implementation detail and will not necessarily be stable.
RootedTrees.α
— Methodα(t::AbstractRootedTree)
The number of monotonic labelings of t
not equivalent under the symmetry group.
Reference: Section 302 of
- Butcher, John Charles. Numerical methods for ordinary differential equations. John Wiley & Sons, 2008.
RootedTrees.β
— Methodβ(t::AbstractRootedTree)
The total number of labelings of t
not equivalent under the symmetry group.
Reference: Section 302 of
- Butcher, John Charles. Numerical methods for ordinary differential equations. John Wiley & Sons, 2008.