From Hopf / PD / Branch point to periodic orbits

From Hopf point to periodic orbits

In order to compute the bifurcated branch of periodic solutions at a Hopf bifurcation point, you need to choose a method to compute periodic orbits among:

Once you have decided which method to use, you can call the following:

continuation(br::ContResult, ind_HOPF::Int, _contParams::ContinuationPar,
	prob::AbstractPeriodicOrbitProblem ;
	δp = nothing, ampfactor = 1, kwargs...)

We refer to continuation for more information about the arguments. Here, we just say a few words about how we can specify prob::AbstractPeriodicOrbitProblem.

See Branch switching (Hopf point) for the precise method definition

Algorithm

The algorithm proceeds as follows. The normal form of the Hopf bifurcation is first computed. Then a predictor for the bifurcated branch of periodic orbits is generated from the normal form. Finally, this predictor is used as a guess for the computation of periodic orbits.

Example

The simplest example is from the getting-started section which we repeat partially below. Several examples are provided in example ODE. In the case of PDE, you can have a look at Brusselator or 2d Ginzburg-Landau equation.

We compute a branch with a Hopf bifurcation:

using BifurcationKit, Parameters, Plots

function Fsl(X, p)
    @unpack r, μ, ν, c3 = p
    u, v = X
    ua = u^2 + v^2
    [
        r * u - ν * v - ua * (c3 * u - μ * v)
        r * v + ν * u - ua * (c3 * v + μ * u)
    ]
end

par_sl = (r = 0.1, μ = 0., ν = 1.0, c3 = 1.0)
u0 = zeros(2)
prob = BifurcationProblem(Fsl, u0, par_sl, (@lens _.r))
opts = ContinuationPar()
br = continuation(prob, PALC(), opts, bothside = true)
 ┌─ Curve type: EquilibriumCont
 ├─ Number of points: 25
 ├─ Type of vectors: Vector{Float64}
 ├─ Parameter r starts at -1.0, ends at 1.0
 ├─ Algo: PALC
 └─ Special points:

- #  1, endpoint at r ≈ -1.00000000,                                                                     step =   0
- #  2,     hopf at r ≈ -0.00595553 ∈ (-0.00595553, +0.00299379), |δp|=9e-03, [converged], δ = ( 2,  2), step =   8
- #  3, endpoint at r ≈ +1.00000000,                                                                     step =  24

We then compute the branch of periodic solutions using orthogonal collocation (for example):

br_po = continuation(br, 2, opts,
        PeriodicOrbitOCollProblem(20, 5)
        )
plot(br, br_po)
Example block output

From Period-doubling point to curve of periodic orbits

Case of Shooting and Collocation

We provide an automatic branching procedure in this case. In essence, all you have to do is to call

continuation(br::ContResult, ind_PD::Int, _contParams::ContinuationPar;
    prm = true, detailed = true,
    kwargs...)

The option prm = true enforces that the period-doubling normal form is computed using the Poincaré return map ; this is only necessary in case of use of the collocation method. Indeed, in the case of the collocation method, an automatic procedure based on the Iooss normal form has yet to be implemented.

Case of Trapezoid method

We do not provide (for now) the automatic branching procedure for these bifurcations of periodic orbits. As a consequence, the user is asked to provide the amplitude of the bifurcated solution.

The call is as follows. Please note that a deflation is included in this method to simplify branch switching.

continuation(br::AbstractBranchResult, ind_PD::Int, contParams::ContinuationPar;
	δp = 0.1, ampfactor = 1, 
	usedeflation = false,
	kwargs...)

An example of use is provided in Lur'e problem.

Algorithm

The algorithm proceeds as follows. The normal form of the Period-doubling bifurcation is first computed. Then a predictor for the bifurcated branch of periodic orbits is generated from the normal form. Finally, this predictor is used as a guess for the computation of periodic orbits.

From Branch point to curve of periodic orbits

We do not provide (for now) the automatic branching procedure for these bifurcations of periodic orbits. As a consequence, the user is asked to provide the amplitude of the bifurcated solution.

We provide the branching method for the following methods to compute periodic orbits: PeriodicOrbitTrapProblem, PeriodicOrbitOCollProblem, ShootingProblem and PoincareShootingProblem. The call is as follows. Please note that a deflation is included in this method to simplify branch switching.

An example of use is provided in Lur'e problem.

continuation(br::AbstractBranchResult, ind_PD::Int, contParams::ContinuationPar;
	δp = 0.1, ampfactor = 1, usedeflation = false, kwargs...)