Invariant Foliations

Here we lay out the API for invariant foliations

FoliationsManifoldsAutoencoders.FoliationIdentifyFunction
FoliationIdentify(dataIN, dataOUT, Tstep, 
                  embedscales, SysName, freq;
                  orders = (P=9,Q=1,U=5,W=3), 
                  iterations = (f=4000, l=30), kappa = 0.2)

This is a convenience method that integrates all steps of calculating a reduced order model.

  1. Estimates the linear dynamics about the origin, only using part of the data.
  2. Calculates an invariant foliation about the invariant subspace corresponding, which has the closest frequancy to freq in Hz.
  3. Performs a normal form transformation of the invariant foliation
  4. Calculates a locally accurate complementary invariant foliation, containing the invariant manifold that has the same dynamics as the foliation calculated in point 2.
  5. Extracts the invariant manifold from the locally accurate invariant foliation.
  6. Performs a correction of the instantaneous frequency and damping values using the invariant manifold calculated in point 5.

Input:

  • dataIN is a two dimensional array, each column is an $\boldsymbol{x}_k$ value,
  • dataOUT is a two dimensional array, each column is an $\boldsymbol{y}_k$
  • Tstep is the time step between $\boldsymbol{x}_k$ and $\boldsymbol{y}_k$
  • embedscales, denoted by $\boldsymbol{w}$ is a matrix applied to $\boldsymbol{x}_k$, $\boldsymbol{y}_k$, used to calculate the amplitude of a signal
  • freq is the frequency that the invariant foliation is calculated for.
  • orders is a named tuple, specifies the polynomial order of $\boldsymbol{P}$, $\boldsymbol{Q}$ and $\boldsymbol{U}$.
  • iterations is a named tuple, f is the maximum itearion when solving for a foliation, l is the maximum number of iterations when solving for a locally accurate foliation.
  • kappa is the size of the neighbourhood of a the invariant manifold considered, when calculating a locally accurate foliation.

Output: is a tuple with elements

  1. vector of instantaneous frequencies
  2. vector of instantaneous damping
  3. vector of instantaneous amplitudes
  4. 1
  5. 1
  6. 1
  7. 1
  8. 1
  9. amplitude of each data point using $\left\| \boldsymbol{w}\, \boldsymbol{U}\left( \boldsymbol{W} \left( \boldsymbol{x}_k\right)\right) \right\|$
source
FoliationsManifoldsAutoencoders.ISFPadeManifoldType
ISFPadeManifold(mdim, ndim, Porder, Qorder, Uorder, B=nothing, field::AbstractNumbers=ℝ)

Returns an ISFPadeManifold object, which provides the matrix manifold structure for an invariant foliation. The invariance equation, where these appear is

\[\boldsymbol{P}\left(\boldsymbol{U}\left(\boldsymbol{x}\right)\right) = \boldsymbol{Q}\left(\boldsymbol{U}\left( \boldsymbol{F}(\boldsymbol{x}_{k})\right)\right)\]

where $\boldsymbol{U}$ is a polynomial with HT tensor coefficients, $\boldsymbol{P}$ is a general dense polynomial and $\boldsymbol{Q}$ is a near identity polynomial, that is $D \boldsymbol{Q}(\boldsymbol{0}) = \boldsymbol{I}$. The purpose of polynomial $\boldsymbol{Q}$ is to have a Pade approximated nonlinear map, $\boldsymbol{S} = \boldsymbol{Q}^{-1} \circ \boldsymbol{P}$. This can balance polynomials on both sides of the invariance equation. In practice, we did not find much use for it yet.

The parameters are

  • mdim: co-dimension of the foliation
  • ndim: the dimesnion of the underlying phase space
  • Porder: order of polynomial $\boldsymbol{P}$
  • Qorder: order of polynomial $\boldsymbol{Q}$
  • Uorder: order of polynomial $\boldsymbol{U}$
  • B: the matrix $\boldsymbol{W}_1$, such that $\boldsymbol{U} (\boldsymbol{W}_1 \boldsymbol{z})$ is constraing to be linear.
  • fields: dummy, a standard parameter of Manifolds.jl
source
Base.zeroMethod
zero(M::ISFPadeManifold)

Creates a zero ISFPadeManifold data representation.

source
FoliationsManifoldsAutoencoders.GaussSouthwellLinearSetupFunction
GaussSouthwellLinearSetup(Misf, Xisf, dataINorig, dataOUTorig, Tstep, nearest; perbox=2000, retbox=4, nbox=10, exclude = false)

Sets up the invariant foliation with linear estimates. The linear dynamics is estimated by the matrix $\boldsymbol{A}$, the left invariant subspace is approximated by the orthogonal matrix $\boldsymbol{U}_1$, the right invariant subspace is approximated by the orthogonal matrix $\boldsymbol{W}_1$. The linearised dynamics is the matrix $\boldsymbol{S}_1$, such that $\boldsymbol{U}_1 \boldsymbol{A}=\boldsymbol{S}_1 \boldsymbol{U}_1$ and $\boldsymbol{A} \boldsymbol{W}_1=\boldsymbol{W}_1 \boldsymbol{S}_1$.

The routine then sets $D \boldsymbol{U} (0) = \boldsymbol{U}_1$ and $D \boldsymbol{P} (0) = \boldsymbol{S}_1$. It also sets the constrint that $\boldsymbol{U} (\boldsymbol{W}_1 \boldsymbol{z})$ is linear.

source
FoliationsManifoldsAutoencoders.GaussSouthwellOptimFunction
GaussSouthwellOptim(Misf, Xisf, dataIN, dataOUT, scale, Tstep, nearest; name = "", maxit=8000, gradstop = 1e-10)

Solves the optimisation problem

\[\arg\min_{\boldsymbol{S},\boldsymbol{U}}\sum_{k=1}^{N}\left\Vert \boldsymbol{x}_{k}\right\Vert ^{-2}\left\Vert \boldsymbol{P}\left(\boldsymbol{U}\left(\boldsymbol{x}_{k}\right)\right)-\boldsymbol{Q}\left(\boldsymbol{U}\left(\boldsymbol{y}_{k}\right)\right)\right\Vert ^{2}.\]

The method is block coordinate descent, where we optimise for matrix coefficients of the representation in a cyclic manner and also by choosing the coefficient whose gradient is the gratest for optimisation.

source