API

Finding Tikhonov-Fenichel Parameter Values

TikhonovFenichelReductions.tfpvs_and_varietiesMethod
tfpvs_and_varieties(
    problem::ReductionProblem;
    preset
) -> Tuple{Vector{Vector{Bool}}, Vector{Vector{Variety}}}

Find all slow-fast separations π⁺ that are TFPVs by using the necessary conditions

  • the affine variety V(f0) contains an irreducible component Y of dimension s
  • the s-th coefficient of the characteristic polynomial of D₁f(x,π⁺) is non-zero for x∈Y

preset can be used to restrict the search for TFPVs to slow-fast separations of rates with predefined slow/fast rates, e.g. for preset=(α=0, β=1) only TFPVs with α ∈ O(ε) and β ∈ O(1) are considered.

Description

The irreducible components are obtained by computing a minimal primary decomposition. The Jacobian at a point in an irreducible component Y is constructed symbolically by computing normal forms with respect to a Gröbner basis G, s.t. G generates the vanishing ideal of Y. To obtain all general TFPVs and not just slow-fast separations, one can use the function tfpvs_groebner.

See also: tfpvs_groebner, print_results, print_tfpvs, print_varieties

source
TikhonovFenichelReductions.tfpvs_groebnerMethod
tfpvs_groebner(
    problem::ReductionProblem
) -> Vector{Nemo.QQMPolyRingElem}

Use the necessary conditions imposed by the RHS of the polynomial ODE System and its Jacobian to filter out TFPV candidates. All TFPVs must result in the vanishing of the returned set of polynomials.

Description

If π⁺ is a TFPV yielding a reduction onto an s-dimensional slow manifold, there exist a point x₀, such that

  • f(x₀,π⁺)=0
  • for any k>s the determinants of all k×k minors of D₁f(x₀,π⁺) vanish

These properties can be used to filter out possible TFPV candidates.

We are interested in partial solutions of the system of polynomials defined by the conditions above. In particular, we only consider conditions on the parameters, since there might be multiple slow manifolds and reductions. Thus, we eliminate the state variables from the ideal generated by the polynomial conditions above. This function computes a generating set for this elimination ideal and all TFPVs lie in its vanishing set.

See also: tfpvs_and_varieties

source
TikhonovFenichelReductions.ReductionProblemType

Type that defines a Tikhonov-Fenichel reduction problem, i.e. a polynomial ODE system, for which all slow-fast separations of rates yielding a reduction onto dimension s are considered.

The type QQMPolyRingElem is used in Oscar.jl to represent elements of a polynomial ring with coefficients in .

Fields

  • f::Vector{Nemo.QQMPolyRingElem}: RHS of ODE system as a vector of polynomials

  • x::Vector{Nemo.QQMPolyRingElem}: Vector of state variables

  • p::Vector{Nemo.QQMPolyRingElem}: Vector of all parameters

  • s::Int64: Dimension of reduced system

  • J::AbstractAlgebra.Generic.MatSpaceElem{Nemo.QQMPolyRingElem}: Jacobian of f w.r.t. x

  • _f::Function: RHS of ODE system as a Julia function with signature _f(x,p)

  • _F::AbstractAlgebra.Generic.FracField{Nemo.QQMPolyRingElem}: ℚ(p,x): Fraction field over ℚ[p,x]

  • _Fp::AbstractAlgebra.Generic.RationalFunctionField{Nemo.QQFieldElem, Nemo.QQMPolyRingElem}: ℚ(p): Field of rational functions in the parameters

  • _Rx::AbstractAlgebra.MPolyRing{AbstractAlgebra.Generic.RationalFunctionFieldElem{Nemo.QQFieldElem, Nemo.QQMPolyRingElem}}: ℚ(p)[x]: Polynomial ring in x over the field Fp

  • _f_Rx::Vector{AbstractAlgebra.Generic.MPoly{AbstractAlgebra.Generic.RationalFunctionFieldElem{Nemo.QQFieldElem, Nemo.QQMPolyRingElem}}}: RHS of ODE system in polynomial ring Rx=ℚ(p)[x]

  • _x_Rx::Vector{AbstractAlgebra.Generic.MPoly{AbstractAlgebra.Generic.RationalFunctionFieldElem{Nemo.QQFieldElem, Nemo.QQMPolyRingElem}}}: Vector of state variables in polynomial ring Rx=ℚ(p)[x]

source
TikhonovFenichelReductions.ReductionProblemMethod
ReductionProblem(
    f::Function,
    x::Array{T<:Union{String, Symbol}, 1},
    p::Array{T<:Union{String, Symbol}, 1},
    s::Int64
) -> ReductionProblem

Constructor for ReductionProblem Type.

Arguments

  • f(x,p)::Function: Julia function defining the RHS of the ODE system
  • x::Vector{String}: Vector of state variables
  • p::Vector{String}: Vector of all parameters
  • s::Int: Dimension of reduced system

Description

This function is used to set up the problem of finding Tikhonov-Fenichel Parameter Values for dimension s. The names of all variables and parameters in the system are parsed to appropriate types in Oscar.jl, so that the necessary conditions for the existence of a reduction onto an s-dimensional slow manifold can be evaluated.

See also: tfpvs_and_varieties, tfpvs_groebner

source
TikhonovFenichelReductions.VarietyType

Type that holds information on the slow manifold as an irreducible component of the variety V(f0).

  • ideal::Oscar.MPolyIdeal: associated ideal in Rx=ℚ(p)[x]

  • gens_R::Vector{Nemo.QQMPolyRingElem}: generators of associated ideal parsed to R=ℚ[p,x]

  • groebner_basis::Vector{AbstractAlgebra.Generic.MPoly{AbstractAlgebra.Generic.RationalFunctionFieldElem{Nemo.QQFieldElem, Nemo.QQMPolyRingElem}}}: groebner basis of associated ideal in Rx=ℚ(p)[x]

  • groebner_basis_R::Vector{Nemo.QQMPolyRingElem}: groebner basis of associated ideal parsed to R=ℚ[p,x]

  • T::AbstractAlgebra.Generic.MatSpaceElem{AbstractAlgebra.Generic.MPoly{AbstractAlgebra.Generic.RationalFunctionFieldElem{Nemo.QQFieldElem, Nemo.QQMPolyRingElem}}}: Matrix T, such that gens(I)*T=groebner_basis

  • dim::Int64: Krull dimension of associated ideal

source

Computing a Tikhonov-Fenichel Reduction

TikhonovFenichelReductions.compute_reduction!Method
compute_reduction!(reduction::Reduction) -> Bool

Compute the reduced system after the slow manifold, non-singular point and product decomposition have been set successfully.

The function returns true if the reduced system was computed successfully. The reduction in raw form, i.e. before substituting the variables x according to the slow manifold is set reduction.g_raw while the s-dimensional reduction on the slow manifold is given by reduction.g. A safe getter function for this is get_reduced_system(reduction::Reduction).

See also: set_manifold!, set_decomposition!, set_point!, compute_all_reductions, Reduction, print_reduced_system

source
TikhonovFenichelReductions.get_explicit_manifoldMethod
get_explicit_manifold(
    problem::ReductionProblem,
    variety::Variety
) -> Tuple{Any, Any}

Heuristic approach to get an explicit (i.e. parameterized) representation of the slow manifold from a variety. If succesfull, the function returns the (attempted) explicit manifold together with a boolean value indicating whether the manifold could be computed automatically.

source
TikhonovFenichelReductions.get_varietiesMethod
get_varieties(
    problem::ReductionProblem,
    sf_separation::Vector{Bool}
) -> Tuple{Bool, Vector{Oscar.MPolyIdeal{AbstractAlgebra.Generic.MPoly{AbstractAlgebra.Generic.RationalFunctionFieldElem{Nemo.QQFieldElem, Nemo.QQMPolyRingElem}}}}, Union{Vector{Any}, Vector{Int64}, Vector{Nemo.NegInf}}}

Compute the irreducible components of V(f0) and their dimensions for a given TFPV candidate. If there exist a reduction, the corresponding slow manifold must be contained in one of these components.

Arguments

  • problem: Reduction problem type holding information on system and dimension of reduction.
  • sf_separation: Boolean index indicating slow-fast separation of rates (0: small, 1: large).

Description

This function can be used if one wants to check whether a particular slow-fast separation of rates yields a reduction for any dimension. If the dimension of an irreducible component of V(f0) differs from what was defined with ReductionProblem, the constructor Reduction can be called with the additional argument s specifying the dimension.

See also: Reduction

source
TikhonovFenichelReductions.set_decomposition!Method
set_decomposition!(
    reduction::Reduction,
    P::Union{AbstractAlgebra.Generic.MatSpaceElem, VecOrMat},
    Psi::Union{AbstractAlgebra.Generic.MatSpaceElem, Vector{Nemo.QQMPolyRingElem}}
) -> Bool

Set product decomposition f0=P⋅Psi locally satisfying V(f0)=V(Psi), where P is a matrix of rational functions and Psi is a vector of polynomials.

See also: set_manifold! set_point! Reduction

source
TikhonovFenichelReductions.set_decomposition!Method
set_decomposition!(
    reduction::Reduction,
    Psi::Union{Variety, Vector{Nemo.QQMPolyRingElem}}
) -> Bool

Try to automatically compute matrix of rational functions P from given vector of polynomials Psi, such that f0=P⋅Psi and V(f0)=V(Psi) holds locally.

NOTE: This always works if the drop in dimension r=n-s=1, but may fail for r>1 (if the number of generators for the irreducible component of V(f0) is greater than r).

Description

Psi can be chosen from r algebraically independent entries of f0. Practically, one can use the generators of the ideals defining the irreducible components of V(f0) as entries for Psi (possibly rewriting the rational equations as polynomials by multiplying appropriately with parameters occurring in a denominator).

source
TikhonovFenichelReductions.ReductionType

Type that holds all information to compute a Tikhonov-Fenichel reduction for a given slow-fast separation of rates.

Fields

  • problem::ReductionProblem: information on input system

  • sf_separation::Vector{Bool}: slow-fast separation (0: slow, 1: fast)

  • _p::Vector{Nemo.QQMPolyRingElem}: Parameters of the system where small ones are set to 0

  • f0::Vector{Nemo.QQMPolyRingElem}: RHS of system as vector with elements of ring R

  • f1::Vector{Nemo.QQMPolyRingElem}: Slow part of the system

  • higher_order_terms::Vector{Nemo.QQMPolyRingElem}: Second order terms in ε

  • Df0::AbstractAlgebra.Generic.MatSpaceElem{Nemo.QQMPolyRingElem}: Jacobian of f0

  • Df0_at_x0::AbstractAlgebra.Generic.MatSpaceElem{AbstractAlgebra.Generic.FracFieldElem{Nemo.QQMPolyRingElem}}: Jacobian of f0 at the non-singular point x0

  • T::AbstractAlgebra.Generic.PolyRing{AbstractAlgebra.Generic.FracFieldElem{Nemo.QQMPolyRingElem}}: Polynomial ring Q(p,x)[λ] for characteristic polynomial of Df0

  • chi::AbstractAlgebra.Generic.Poly{AbstractAlgebra.Generic.FracFieldElem{Nemo.QQMPolyRingElem}}: Characteristic polynomial of Df0

  • M::Vector{AbstractAlgebra.Generic.FracFieldElem{Nemo.QQMPolyRingElem}}: Components of the system on slow manifold

  • x0::Vector{AbstractAlgebra.Generic.FracFieldElem{Nemo.QQMPolyRingElem}}: Non-singular point in the irreducible component of V(f0) containing the slow manifold

  • P::AbstractAlgebra.Generic.MatSpaceElem{AbstractAlgebra.Generic.FracFieldElem{Nemo.QQMPolyRingElem}}: Matrix with rational functions, such that f0=P⋅Psi

  • Psi::AbstractAlgebra.Generic.MatSpaceElem{AbstractAlgebra.Generic.FracFieldElem{Nemo.QQMPolyRingElem}}: Vector with polynomials, such that f0=P⋅Psi

  • DPsi::AbstractAlgebra.Generic.MatSpaceElem{AbstractAlgebra.Generic.FracFieldElem{Nemo.QQMPolyRingElem}}: Jacobian of Psi

  • success::Vector{Bool}: Indicates whether slow manifold M, non-singular point x0 and product decomposition f0=P⋅Psi have been set successfully to allow the computation of the reduced system

  • idx_components::Vector{Bool}: Boolean indices of components that determine the flow on the slow manifold

  • g_raw::Vector{AbstractAlgebra.Generic.FracFieldElem{Nemo.QQMPolyRingElem}}: Reduced system in general form (before substituting variables according to the slow manifold)

  • g::Vector{AbstractAlgebra.Generic.FracFieldElem{Nemo.QQMPolyRingElem}}: Reduced system on the slow manifold

  • reduction_cached::Vector{Bool}: Whether g and g_raw are already computed

source
TikhonovFenichelReductions.ReductionMethod
Reduction(
    problem::ReductionProblem,
    sf_separation::Vector{Bool},
    V::Union{Variety, Vector{Nemo.QQMPolyRingElem}},
    M::Vector{<:AbstractAlgebra.RingElem};
    s
) -> Reduction

Convenience function that constructs an object of type Reduction and calls set_manifold! and set_decomposition!.

Arguments

  • problem: Reduction problem type holding information on system and dimension of reduction.
  • sf_separation: Boolean index indicating slow-fast separation of rates (0: small, 1: large).
  • V: Generators of affine variety corresponding to the the slow manifold
  • M: Slow manifold in explicit form
  • s::Int: (optional) Dimension of slow manifold. Can be specified if a reduction corresponding to a TFPV for dimension different from problem.s should be considered (e.g. for manually computing a reduction for a given slow-fast separation that is not necessarily obtained via tfpvs_and_varieties).

See also: set_manifold! set_decomposition! tfpvs_and_varieties

source
TikhonovFenichelReductions.ReductionMethod
Reduction(
    problem::ReductionProblem,
    sf_separation::Vector{Bool};
    s
) -> Reduction

Constructor for Reduction Type.

Arguments

  • problem: Reduction problem type holding information on system and dimension of reduction.
  • sf_separation: Boolean index indicating slow-fast separation of rates (0: small, 1: large).
  • s::Int: (optional) Dimension of slow manifold. Can be specified if a reduction corresponding to a TFPV for dimension different from problem.s should be considered (e.g. for manually computing a reduction for a given slow-fast separation that is not necessarily obtained via tfpvs_and_varieties).

See also: set_manifold! set_decomposition!

source

Computing Multiple Tikhonov-Fenichel Reductions Simultaneously

TikhonovFenichelReductions.compute_all_reductionsMethod
compute_all_reductions(
    problem::ReductionProblem,
    sf_separations::Vector{Vector{Bool}},
    varieties::Vector{Vector{Variety}},
    M::Union{Vector{Vector{AbstractAlgebra.Generic.FracFieldElem{Nemo.QQMPolyRingElem}}}, Vector{Vector{Nemo.QQMPolyRingElem}}};
    print
) -> Tuple{Vector{Vector{Reduction}}, Vector{Vector{Tuple{Int64, Int64}}}}

Compute reductions for TFPVs sf_separations onto each of the explicitly given slow manifolds in M. Returns the reductions R in the same format as varieties and indices idx such that M[k] corresponds to varieties[i][j] and R[i][j] for each (i,j) in idx[k].

Description

All possible choices of slow manifolds can be obtained with unique_varieties. For each of these (that one is interested in), a parameterized form, i.e. the slow manifold in phase space defined by the remaining components of the system, must be provided. In most cases this can be obtained with get_explicit_manifold. The backwards step, i.e. finding all manifolds for which this is an explicit description, is then done automatically with find_varieties.

See also: unique_varieties, compute_all_reductions, compute_reduction!, Reduction, set_manifold!, set_decomposition!

source
TikhonovFenichelReductions.find_varietiesMethod
find_varieties(
    problem::ReductionProblem,
    varieties::Vector{Vector{Variety}},
    M::Union{Vector{AbstractAlgebra.Generic.FracFieldElem{Nemo.QQMPolyRingElem}}, Vector{Nemo.QQMPolyRingElem}}
) -> Vector{Tuple{Int64, Int64}}

Find all varieties as returned by tfpvs_and_varieties for a given explicit manifold, i.e. a parametric description as a vector with the same size as problem.x.

See also: compute_all_reductions

source

Output

TikhonovFenichelReductions.print_resultsFunction
print_results(
  [io::IO,]
  problem::ReductionProblem,
  sf_separations::Vector{Vector{Bool}},
  V::Vector{Vector{Variety}};
  idx::Union{AbstractVector{Int}, Vector{Bool}}=1:length(V)
)

Print slow-fast separations that are TFPVs and generators for the corresponding irreducible components of V(f0) together with their Krull dimension as stored in Variety object.

Arguments

  • problem: ReductionProblem
  • sf_separations: Boolean indices defining all TFPVs π⁺ that are slow-fast separations (0: slow, 1: fast).
  • V: Generators for the irreducible component of the affine varietiy V(f(⋅,π⁺)) for each slow-fast separation and their respective dimension.
  • idx: (optional) index vector to include only certain TFPVs (boolean or numeric)

See also: tfpvs_and_varieties, print_tfpvs, print_varieties

source
TikhonovFenichelReductions.print_tfpvsFunction
print_tfpvs(
  [io::IO,] 
  problem::ReductionProblem,
  sf_separations::Vector{Vector{Bool}}; 
  latex::Bool=false,
  idx::Union{AbstractVector{Int}, Vector{Bool}}=1:length(sf_separations)
)

Print slow-fast separations to terminal or use latex=true to print LaTeX instead. Optionally only print subset defined by (boolean or numeric) index set idx.

source
TikhonovFenichelReductions.print_varietiesFunction
print_varieties(
  [io::IO,]
  V::Vector{Vector{Variety}};
  latex::Bool=false,
  idx::Union{AbstractVector{Int}, Vector{Bool}}=1:length(V)
)

Print generators of ideals corresponding to the irreducible components of varieties V(f0) for TFPV candidates and their dimension as returned by tfpvs_and_varieties (these are wrapped in the type Variety). Use keyword argument latex=true to print LaTeX code instead.

See also: tfpvs_and_varieties, print_tfpvs, print_results

source
TikhonovFenichelReductions.print_reduced_systemFunction
print_reduced_system(
  [io::IO,] 
  reduction::Reduction; 
  rewrite::Bool=true,
  factor::Bool=false, 
  latex::Bool=false,
  local_coordinates::Bool=true
)

Print the reduced system (after compute_reduction! has been called successfully on the Reduction object). If the slow manifold was successfully specified, this returns the system on the slow manifold. If local_coordinates=true, the reduced system is only given in s-dimensions (i.e. the local coordinates on the slow manifold). If the reduction could be computed but the slow manifold was not set successfully, the reduced system is shown in the original phase space.

Arguments

  • reduction: Reduction holding the reduced system
  • rewrite: Whether the RHS should be decomposed into polynomial and rational part (default is true), see rewrite_rational
  • factor: Whether the polynomial parts should be factorized
  • latex: Whether to print latex string
  • local_coordinates: Whether to present reduced system in local coordinates on slow manifold (i.e. as s dimensional ODE system)

See also: rewrite_rational, print_tfpvs, print_varieties

source
TikhonovFenichelReductions.rewrite_rationalFunction
rewrite_rational(
    term::AbstractAlgebra.Generic.FracFieldElem{Nemo.QQMPolyRingElem};
    factor
) -> Tuple{Nemo.QQMPolyRingElem, Nemo.QQMPolyRingElem, Nemo.QQMPolyRingElem}

Decompose a rational function f = p/q into polynomial and rational part, i.e. return f = h + r/q, where p,q,h,r are polynomials.

Arguments

  • factor: Whether all polynomials should be factorized (default is false)

See also: print_reduced_system

source