mercoledì 27 agosto 2014

Dense Output

Dense Output

Specifying specific output times for the solution, should not affect the internal time steps that the solver uses. The basic idea of the Dense Output concept is to provide the solution at a given time \(s \in [t, t+dt]\) with the same order of accuracy as the solutions computed at the internal time points by using suitable interpolation methods.
Up to now only linear interpolation was performed and this significantly lowered the accuracy if a higher order solver was used.
I then implemented a series of interpolation function:

  • linear_interpolation:
    x_out = linear_interpolation (t, x, t_out)
    Given the time span \(t=[t_0, t_1]\) and the function values \(x=[x_0, x_1]\), it returns the linear interpolation value \(x_{out}\) at the point \(t_{out}\).
  • quadratic_interpolation:
    x_out = quadratic_interpolation (t, x, der, t_out)
    Given the time span \(t=[t_0, t_1]\), the function values \(x=[x_0, x_1]\) and the derivative of the function at the point \(x_0\), it returns the quadratic interpolation value \(x_{out}\) at the point \(t_{out}\).
  • hermite_cubic_interpolation:
    x_out = hermite_cubic_interpolation (t, x, der, t_out)
    Given the time span \(t=[t_0, t_1]\), the function values \(x=[x_0, x_1]\) and the derivatives of the function at both points \(x_0\) and \(x_1\), it returns the 3rd order approximation \(x_{out}\) at the point \(t_{out}\) by performing Hermite interpolation.
  • hermite_quartic_interpolation:
    x_out = hermite_quartic_interpolation (t, x, der, t_out)
    Given the time span \(t=[t_0, t_1]\), the function values \(x=[x_0, x_{1/2}, x_1]\) (where \(x_{1/2}\) is the value of the function at the time \(t_0+dt/2\)) and the derivatives of the function at the extremes \(x0\) and \(x1\), it returns the 4th order approximation \(x_{out}\) at the point \(t_{out}\) by performing Hermite interpolation.
  • dorpri_interpolation:
    x_out = dorpri_interpolation (t, x, k, t_out)
    This interpolation method is specific for the Dormand-Prince Runge-Kutta scheme. Given the time span \(t=[t_0, t_1]\), the function value \(x=x_0\) and the vector \(k\) with the function evaluations required in the Dormand-Prince method, it returns the 4th order approximation \(x_{out}\) at the point \(t_{out}\). For more information on the method have a look at Hairer, Noersett, Wanner, Solving Ordinary Differential Equations I: Nonstiff Problems (pag. 191-193).
  • hermite_quintic_interpolation:
    x_out = hermite_quintic_interpolation (t, x, der, t_out)
    Given the time span \(t=[t_0, t_1]\), the function values \(x=[x_0, x_{1/2}, x_1]\) and the derivatives of the function at each point, it returns the 5th order approximation \(x_{out}\) at the point \(t_{out}\) by performing Hermite interpolation.
These methods are then used to perform the Dense Output according to the order of the solver chosen. This is the piece of code in integrate_adaptive.m that performs the interpolation:

% if next tspan value is caught, update counter
if( (z(end) == tspan(counter)) || (abs (z(end) - tspan(counter)) / ...
    (max (abs (z(end)), abs (tspan(counter)))) < 8*eps) )
  counter++;

% if there is an element in time vector at which the solution is required
% the program must compute this solution before going on with next steps
elseif( vdirection*z(end) > vdirection*tspan(counter) )
% initializing counter for the following cycle
  i = 2;
  while ( i <= length (z) )

    % if next tspan value is caught, update counter
    if( (counter <= k) && ...
        ( (z(i) == tspan(counter)) || (abs (z(i) - tspan(counter)) / ...
        (max (abs (z(i)), abs (tspan(counter)))) < 8*eps)) )
      counter++;
    endif
    
    % else, loop until there are requested values inside this subinterval
    while((counter <= k) && (vdirection*z(i) > vdirection*tspan(counter)))
      % choose interpolation scheme according to order of the solver
      switch order
      case 1
        u_interp = linear_interpolation ([z(i-1) z(i)], [u(:,i-1) u(:,i)], tspan(counter));
      case 2
        if (~isempty (k_vals))
          der = k_vals(1);
        else
          der = feval (func, z(i-1) , u(:,i-1), options.vfunarguments{:});
        endif
        u_interp = quadratic_interpolation ([z(i-1) z(i)], [u(:,i-1) u(:,i)], der, tspan(counter));
      case 3
        % only ode23 - use k_vals
        u_interp = hermite_cubic_interpolation ([z(i-1) z(i)], [u(:,i-1) u(:,i)], [k_vals(:,1) k_vals(:,end)], tspan(counter));
      case 4
        % if ode45 used without local extrapolation this function doesn't require a new function evaluation.
        u_interp = dorpri_interpolation ([z(i-1) z(i)], [u(:,i-1) u(:,i)], k_vals, tspan(counter));
      case 5
        % ode45 with Dormand-Prince scheme:
        % 4th order approximation of y in t+dt/2 as proposed by Shampine in Lawrence, Shampine, "Some Practical Runge-Kutta Formulas", 1986.
        u_half = u(:,i-1) + 1/2*dt*((6025192743/30085553152)*k_vals(:,1) + (51252292925/65400821598)*k_vals(:,3) - (2691868925/45128329728)*k_vals(:,4) + (187940372067/1594534317056)*k_vals(:,5) - (1776094331/19743644256)*k_vals(:,6) + (11237099/235043384)*k_vals(:,7));
        u_interp = hermite_quartic_interpolation ([z(i-1) z(i)], [u(:,i-1) u_half u(:,i)], [k_vals(:,1) k_vals(:,end)], tspan(counter));

        % it is also possible to do a new function evaluation and the quintic hermite interpolator
        %f_half = feval (func, t+1/2*dt, u_half, options.vfunarguments{:});
        %u_interp = hermite_quintic_interpolation ([z(i-1) z(i)], [u(:,i-1) u_half u(:,i)], [k_vals(:,1) f_half k_vals(:,end)], tspan(counter));
      otherwise
        warning ('high order interpolation not yet implemented: using cubic iterpolation instead');
        der(:,1) = feval (func, z(i-1) , u(:,i-1), options.vfunarguments{:});
        der(:,2) = feval (func, z(i) , u(:,i), options.vfunarguments{:});
        u_interp = hermite_cubic_interpolation ([z(i-1) z(i)], [u(:,i-1) u(:,i)], der, tspan(counter));
      end

      % add the interpolated value of the solution
      u = [u(:,1:i-1), u_interp, u(:,i:end)];
            
      % add the time requested
      z = [z(1:i-1);tspan(counter);z(i:end)];

      % update counters
      counter++;
      i++;
    endwhile

    % if new time requested is not out of this interval
    if ((counter <= k) && (vdirection*z(end) > vdirection*tspan(counter)))
      % update the counter
      i++;
    else
      % else, stop the cycle and go on with the next iteration
      i = length (z) + 1;
    endif

  endwhile

endif

It is important to notice that:

  • The 1st order approximation doesn't require any additional function evaluation.
  • The 2nd order approximation may require the evaluation of the function at the current time. This can be avoided if the stepper already returns that value.
  • The only 3rd order solver implemented is ode23. The 3rd order approximation exploits the Runge-Kutta \(k\) values to avoid further function evaluations.
  • There are no 4th order schemes as yet implemented. However if ones were to use ode45 without local extrapolation then the dorpri_interpolation function can be used to obtain a 4th order approximation without any additional function evaluation. For any other 4th order scheme the hermite_quartic_interpolation function can be used.
  • For the 5th order method ode45, Shampine proposes to obtain a 4th order approximation at the middle point and to use quartic interpolation. It is however possible to directly do quintic interpolation but this require an additional function evaluation without (according to Shampine) a significant improvement.
  • For the higher order solvers (ode78), a suitable interpolator has not yet been implemented.
Finally, since I wrote the interpolation functions in such a way that they are independent of the number of output points requested, a possible improvement would be to compute all the values of the dense output between \(t\) and \(t+dt\) all at once instead of one value at a time.

lunedì 18 agosto 2014

FSAL - new stepper implementation

As stated in the previous post, the implementation of the steppers as it was did not allow the possibility to exploit the FSAL (First Same As Last) property of the Bogacki-Shampine algorithm (ode23) and of the Dormand-Prince algorithm (ode45).
The input and output arguments of the steppers have then be modified. As an example here is the runge_kutta_23 stepper:

function varargout = runge_kutta_23 (f, t, x, dt, varargin)
  options = varargin{1};
  k = zeros (size (x, 1), 4);

  if (nargin == 5) % only the options are passed
    k(:,1) = feval (f, t , x, options.vfunarguments{:});
  elseif (nargin == 6) % both the options and the k values are passed
    k(:,1) = varargin{2}(:,end); % FSAL property
  endif
  k(:,2) = feval (f, t + (1/2)*dt, x + dt*(1/2)*k(:,1), options.vfunarguments{:});
  k(:,3) = feval (f, t + (3/4)*dt, x + dt*(3/4)*k(:,2), options.vfunarguments{:});

  %# computing new time and new values for the unkwnowns
  varargout{1} = t + dt; %t_next
  varargout{2} = x + dt.*((2/9)*k(:,1) + (1/3)*k(:,2) + (4/9)*k(:,3)); % return the 3rd order approximation x_next

  %# if the estimation of the error is required
  if (nargout >= 3)
    %# new solution to be compared with the previous one
    k(:,4) = feval (f, t + dt, varargout{2}, options.vfunarguments{:});
    varargout{3} = x + dt.*((7/24)*k(:,1) + (1/4)*k(:,2) + (1/3)*k(:,3) + (1/8)*k(:,4)); %x_est
    varargout{4} = k;
  endif

endfunction

And the call within the solver becomes:
[s, y, y_est, k_vals] = stepper (func, z(end), u(:,end), dt, options, k_vals);

where k_vals has to be initialized for the first iteration as f(t, x).
This implementation will reduce the number of function evaluation for each step.

Furthermore, after some tests in MATLAB, the return values for the solution and the estimate in the  runge_kutta_23 and runge_kutta_45 steppers have been swapped to automatically perform local extrapolation. The MATLAB functions are in fact of order 3 and 5 respectively.

martedì 12 agosto 2014

Status of the code: bugfixes and new issues


ODESET and ODEGET


  • odeset and odeget functions have been slightly modified to be compliant with MATLAB. Each MATLAB option is present and all the options are tested. The coding style has been adapted to the GNU-Octave standard.
  • ode_struct_value_check: this function has been introduced by Roberto in addition to odepkg_structue_check. The relation between the two functions has to be clarified: in particular it is necessary to understand if it is really necessary to have two different functions or one is sufficient.


CHANGES TO THE STEPPERS


  • The runge_kutta_78 stepper has been implemented.
  • Two 4th order steppers have been implemented: runge_kutta_45_dopri (Dormand-Prince coefficients) and runge_kutta_45_fehlberg (Fehlberg coefficients).


CHANGES TO THE SOLVERS


  • ode78 solver has been updated to the new structure. It now exploits the runge_kutta_78 stepper.
  • A series of tests has been added to each solver to check all the functionalities and the all options. This has made me possible to detect some bugs that have been corrected. In particular the adaptive timestep evaluation had some issues that lead to the use of too short timesteps. This has been corrected and now the algorithm proposed in [1] is used.
  • Furthermore the current implementation uses linear interpolation to evaluate the solution at the user specified times. This leads to a considerable loss in accuracy and is not consistent with MATLAB (which guarantees the same order of accuracy of the scheme used). In [1] different methods are proposed for the dense output: these will be used as a reference for the implementation of a better solution.
  • In the released version of odepkg some of the solvers perform local extrapolation, that is the higher-order estimate is chosen as the solution. With the new stepper structure, as it is now, this choice effects all the solvers. It have to be decided whether to perform it or not (MATLAB doesn't seem to do it, thus I suggest to avoid it).
  • MATLAB implementation of ode45 uses the Dormand-Prince (DP) coefficients. In the released version of odepkg there exits two solvers: ode45 that uses the Fehlberg coefficients and ode54 that uses the DP coefficients. To be consistent with MATLAB, ode45 now uses the DP method. This makes the runge_kutta_45_fehlberg stepper and the ode54 solver, as it is now, redundant. Either their elimination or a change of the solver might be considered. However one of the advantages of DP coefficients is the possibility to reuse the last function evaluation at a given step as the first evaluation of the subsequent one. This is not easily done with the stepper structure introduced by Roberto.


CHANGES TO THE OPTIONS


  • InitialStep option has been modified to be MATLAB compatible (it must be a positive scalar).
  • RelTol defalut value has been changed to 1e-3 instead of 1e-6 to be MATLAB compatible.
  • MaxStep option has been implemented.
  • NormControl option has been implemented.


TODO

In addition to the general plan, a couple of issues need to be addressed:

  • Clarify the relation between ode_struct_value_check and odepkg_structue_check.
  • Decide if local extrapolation has to be used or not. My opinion (and the current implementation) is to avoid it to be compliant to what MATLAB seems to be doing.
  • Solve the dense output problem in a way that guarantees the consistency with MATLAB.
  • Consider if it's possible to reduce the number of function evaluation for the Dormand-Prince stepper (ode45) and the Bogacki-Shampine stepper (ode23) exploiting the FSAL property (first same as last).
  • Decide if in the future releases of odepkg ode54 has to be removed or maybe changed to become a Fehlberg solver.



[1] E. Hairer, S.P. N{\o}rsett, G. Wanner, Solving Ordinary Differential Equations, 1993, Springer.

martedì 29 luglio 2014

SoCiS 2014 - New timeline


On the occasion of SoCiS 2014 I will take over the work that has been done last year on odepkg and continue it. The final goal is to release a new stable version of odepkg and to insert the most common solver in core-Octave in such a way that everything is MATLAB-compatible.

The following list explains the main points of the timeline for my SOCIS project:

  1. Check of the current status of the code, in particular with respect to the current release on SourceForge. The two repository will be merged so that every test present in the old version will be added to the new code. Verify if there are missing features and add them if necessary.
  2. Comparison of the performance between the old and the new structure. In particular we expect that the introduction of the Levenshtein algorithm for the string comparisons will be a critical issue. If necessary implement levenshtein.m and fuzzy_compare.m in C++.
  3. Verify that the functions odeset and odeget are MATLAB-compatible and compliant to Octave core. Add the two functions to the core.
  4. Move ode45, ode23 and ode23s to Octave core.
  5. Implement ode15s solver. This solver is still missing in odepkg but is highly suitable for stiff problems.
  6. Move ode15s to Octave core.
  7. New release of odepkg.

domenica 15 dicembre 2013

RATTLE algorithm

RATTLE

Constrained mechanical systems form an important class of differential equations on manifolds. For all the theory that I present here and I've used for the implementation I refer to "Geometric Numerical Integration" by Hairer, Lubich and Wanner.

Consider a mechanical system described by position coordinates \(q_1,\dots q_d\) and suppose that the motion is constrained to satisfy \(g(\mathbf q) = \mathbf 0\) where \(g:\mathbb R^d \rightarrow \mathbb R^m \), with \(m \lt d \). So that, the equations of motion governing the system become: $$ \left\{ \begin{array}{l} \dot{\mathbf q} = \frac{\partial H}{\partial \mathbf p} \\ \dot{\mathbf p} = -\frac{\partial H}{\partial \mathbf q} - G(q)^T\lambda \\ g(\mathbf q) = \mathbf 0 \end{array}\right. $$ where \(G(\mathbf q) = \dfrac{\partial g(\mathbf q)}{\partial \mathbf q} \).

Symplectic Euler method can be extended to constrained systems but we focus on SHAKE and RATTLE algorithms. SHAKE is a 2-steps algorithm so that, since I'm implementing only 1-step algorithms and the overall structure of solvers and integrators is made for 1-step solvers, I implemented just RATTLE algorithm.

odeRATTLE

The RATTLE algorithm implemented works with any general Hamiltonian \(H(\mathbf q,\mathbf p \) and is defined as follows: $$ \left\{\begin{array}{l} \mathbf p_{n+\frac{1}{2}} = \mathbf p_n -\frac{h}{2}\big(H_{\mathbf q}(\mathbf q_n,\mathbf p_{n+\frac{1}{2}}) + G(\mathbf q_n)^T \mathbf {\lambda}_n \big) \\ \mathbf q_{n+1} = \mathbf q_n +\frac{h}{2} \big( H_{\mathbf p}(\mathbf q_n,\mathbf p_{n+\frac{1}{2}}) + H_{\mathbf p}(\mathbf q_{n+1},\mathbf p_{n+\frac{1}{2}}) \big) \\ g(\mathbf q_{n+1}) = \mathbf 0 \\ \mathbf p_{n+1} = \mathbf p_{n+\frac{1}{2}} -\frac{h}{2}\big(H_{\mathbf q}(\mathbf q_{n+1},\mathbf p_{n+\frac{1}{2}}) + G(\mathbf q_{n+1})^T \mathbf{\mu}_n \big) \\ G(\mathbf q_{n+1}) H_{\mathbf p}(\mathbf q_{n+1},\mathbf p_{n+1}) = \mathbf 0 \end{array}\right. $$ where \( h=\Delta t=t_{k+1}-t_k\) and \(\mathbf{\mu}_n \),\( \mathbf{\mu}_n \) are Lagrangian multipliers nedded to impose the constraints.

It can be demonstrated that this numerical method is symmetric, symplectic and convergent of order 2.

The following code represent it's implementation:

function [t_next,x_next,err]=rattle(f,t,x,dt,options)
  H = odeget(options,'HamiltonianHessFcn',[],'fast');
  GG = odeget(options,'ConstraintHessFcn',[],'fast');
  if( ~isempty(H) && ~isempty(GG) )
    fsolve_opts = optimset('Jacobian','on');
  else
    fsolve_opts = optimset('Jacobian','off');
  end
  g = odeget(options,'ConstraintFcn',[],'fast');
  G = odeget(options,'ConstraintGradFcn',[],'fast');
  c_nb = odeget(options,'ConstraintsNb',[],'fast');

  dim = length(x)/2;
  q0 = x(1:dim);
  p0 = x(dim+1:end);

  RATTLE = @(y)constr_sys(y,dim,c_nb,f,H,g,G,GG,t,q0,p0,dt);
  y0 = [q0;p0;p0;zeros(2*c_nb,1)];
  y0 = fsolve(RATTLE,y0,fsolve_opts);
  t_next = t+dt;
  x_next = [y0(1:dim);y0(2*dim+1:3*dim)];

  if(nargout==3)
    dt = dt/2;
    RATTLE = @(y)constr_sys(y,dim,c_nb,f,H,g,G,GG,t,q0,p0,dt);
    y0 = [q0;p0;p0;zeros(2*c_nb,1)];
    y0 = fsolve(RATTLE,y0,fsolve_opts);

    q0 = y0(1:dim);
    p0 = y0(2*dim+1:3*dim);
    t = t+dt;

    RATTLE = @(y)constr_sys(y,dim,c_nb,f,H,g,G,GG,t,q0,p0,dt);
    y0 = [q0;p0;p0;zeros(2*c_nb,1)];
    y0 = fsolve(RATTLE,y0,fsolve_opts);

    x_est = [y0(1:dim);y0(2*dim+1:3*dim)];
    err = norm(x_est-x_next,2);
  end
end

function [F,J] = constr_sys(y,dim,c_nb,f,H,g,G,GG,t,q0,p0,dt)
  F = zeros(3*dim+2*c_nb,1);
  F(1:dim) = y(1:dim) - q0 - (dt/2).*(f(t,[q0; ...
        y(dim+1:2*dim)])(1:dim) + f(t+dt,y(1:2*dim))(1:dim));
  F(dim+1:2*dim) = y(dim+1:2*dim) - p0 - (dt/2).*(f(t,[q0; ...
   y(dim+1:2*dim)])(dim+1:end) - G(q0)'*y(3*dim+1:3*dim+c_nb));
  F(2*dim+1:3*dim) = y(2*dim+1:3*dim) - y(dim+1:2*dim) - ...
             (dt/2)*(f(t+dt,y(1:2*dim))(dim+1:end) - ...
             G(y(1:dim))'*y(3*dim+c_nb+1:end));
  F(3*dim+1:3*dim+c_nb) = g(y(1:dim));
  F(3*dim+c_nb+1:end) = G(y(1:dim))*(f(t+dt,[y(1:dim); ...
                                 y(2*dim+1:3*dim)])(1:dim));

  if( nargout==2 )
    J = zeros(3*dim+2*c_nb,3*dim+2*c_nb);
    J(1:dim,1:dim) = eye(dim) - ...
              (dt/2)*(H(t+dt,y(1:2*dim))(dim+1:end,1:dim));
    J(1:dim,dim+1:2*dim) = -(dt/2)*(H(t,[q0; ...
          y(dim+1:2*dim)])(dim+1:end,dim+1:end) + ...
            H(t+dt,y(1:2*dim))(dim+1:end,dim+1:end));
    J(dim+1:2*dim,dim+1:2*dim) = eye(dim) + ...
        (dt/2)*(H(t,[q0;y(dim+1:2*dim)])(1:dim,dim+1:end));
    J(dim+1:2*dim,3*dim+1:3*dim+c_nb) = (dt/2)*G(q0)';
    J(2*dim+1:3*dim,1:dim) = (dt/2)*(H(t+dt, ...
                            y(1:2*dim))(1:dim,1:dim));
    for k = 1:1:c_nb
      J(2*dim+1:3*dim,1:dim) = J(2*dim+1:3*dim,1:dim) - ...
            (dt/2)*(y(3*dim+c_nb+k)*(GG(y(1:dim))(:,:,k)));
    end
    J(2*dim+1:3*dim,dim+1:2*dim) = -eye(dim) + ...
               (dt/2)*(H(t+dt,y(1:2*dim))(1:dim,dim+1:end));
    J(2*dim+1:3*dim,2*dim+1:3*dim) = eye(dim) + ...
               (dt/2)*(H(t+dt,y(1:2*dim))(1:dim,dim+1:end));
    J(2*dim+1:3*dim,3*dim+c_nb+1:end) = (dt/2)*G(y(1:dim))';
    J(3*dim+1:3*dim+c_nb,1:dim) = G(y(1:dim));
    J(3*dim+c_nb+1:end,1:dim) = G(y(1:dim))* ...
      (H(t+dt,[y(1:dim);y(2*dim+1:3*dim)])(dim+1:end,1:dim));
    for k = 1:1:c_nb
      J(3*dim+c_nb+k,1:dim) = J(3*dim+c_nb+k,1:dim) + ...
          ((GG(y(1:dim))(:,:,k))*(f(t+dt,[y(1:dim); ...
                     y(2*dim+1:3*dim)])(1:dim)))';
    end
    J(3*dim+c_nb+1:end,2*dim+1:3*dim) = G(y(1:dim))* ...
        (H(t+dt,[y(1:dim);y(2*dim+1:3*dim)]) ...
                  (dim+1:end,dim+1:end));
  end
end
It works with any number of constraint, unless this is equal or greater to system dimension. As usual, all the source code is available at my public repository octave-odepkg.

lunedì 9 dicembre 2013

Spectral Variational Integrators

Variational integrators

The variational integrators are a class of numerical methods for mechanical systems which comes from the discrete formulation of Hamilton's principle of stationary action. They can be applied to ODEs, PDEs and to both conservative and forced systems. In absence of forcing terms these methods preserve momenta related to symmetries of the problem and don't dissipate energy. So that they exhibit long-term stability and good long-term behaviour.

Considering a configuration manifold V, the discrete Lagrangian is a function from V to the real numbers space, which represents an approximation of the action between two fixed configurations: $$ L_d(\mathbf q_k,\mathbf q_{k+1}) \approx \int_{t_k}^{t_{k+1}} L(\mathbf q,\dot{\mathbf q};t) dt \hspace{0.5cm}\text{with}\hspace{0.4cm}\mathbf q_{k},\mathbf q_{k+1}\hspace{0.3cm}\text{fixed.}$$ From here, applying the Hamilton's principle, we can get the Euler-Lagrange discrete equations: $$D_2L_d(\mathbf q_{k-1},\mathbf q_k) + D_1 L_d(\mathbf q_k,\mathbf q_{k+1}) = 0\ ,$$ and thanks to the discrete Legendre transforms we get the discrete Hamilton's equation of motion: $$\left\{\begin{array}{l} \mathbf p_k = -D_1 L_d(\mathbf q_k,\mathbf q_{k+1}) \\ \mathbf p_{k+1} = D_2 L_d(\mathbf q_k,\mathbf q_{k+1}) \end{array}\right.\ ,$$ so that we pass from a 2-step system of order N to a 1-step system of order 2N. This system gives the updating map: $$ (\mathbf q_k,\mathbf p_k)\rightarrow(\mathbf q_{k+1},\mathbf p_{k+1})\ .$$ For all the theory behind this I refer to "Discrete mechanics and variational integrators" by Marsden and West.

Spectral variational integrators

To create a spectral variational integrator I considered a discretization of the configuration manifold on a n-dimensional functional space generated by the orthogonal basis of Legendre polynomials. So that, after rescaling the problem from \([t_k,t_{k+1}]\) (with \(t_{k+1}-t_k=h\)) to \([-1,1]\), we get: $$ \begin{array}{l} \mathbf q_k(z) = \sum_{i=0}^{n-1}\mathbf q_k^i l_i(z)\\ \dot{\mathbf q}_k(z) = \frac{2}{h} \sum_{i=0}^{n-1}\mathbf q_k^i \dot l_i(z) \end{array} \ .$$ Then I approximate the action using the Gaussian quadrature rule using \(m\) internal nodes, so putting all together we have: $$ \int_{t_k}^{t_{k+1}} L(\mathbf q,\dot{\mathbf q};t)dt\hspace{0.5cm} \approx\hspace{0.5cm} \frac{h}{2}\sum_{j=0}^{m-1} \omega_j L\big( \sum_{i=0}^{n-1}\mathbf q_k^i l_i , \frac{2}{h} \sum_{i=0}^{n-1}\mathbf q_k^i \dot l_i \big) $$ Now imposing Hamilton's principle of stationary action under the constraints: $$ \mathbf q_k = \sum_{i=0}^{n-1}\mathbf q_k^i l_i(-1) \hspace{1.5cm} \mathbf q_{k+1} = \sum_{i=0}^{n-1}\mathbf q_k^i l_i(1)\ ,$$ we obtain the system: $$ \left\{ \begin{array}{l} \sum_{j=0}^{m-1}\omega_j \bigg[ \mathbf p_j \dot{l}_s(z_j) - \frac{h}{2} l_s(z_j) \frac{\partial H}{\partial \mathbf q} \bigg ( \sum_{i=0}^{n-1}\mathbf q_k^i l_i(z_j),\mathbf p_j \bigg ) \bigg ] + l_s(-1)\mathbf p_k - l_s(1)\mathbf p_{k+1} = 0 \hspace{1cm} \forall s=0,\dots,n-1\\ \frac{\partial H}{\partial \mathbf p}\bigg (\sum_{i=0}^{n-1}\mathbf q_k^i l_i(z_r),\mathbf p_r\bigg ) -\frac{2}{h}\ \sum_{i=0}^{n-1} \mathbf q_k^i \dot{l}_i(z_r)=0 \hspace{1cm} \forall r=0,\dots,m-1 \\ \sum_{i=0}^{n-1} \mathbf q_k^i l_i(-1) - \mathbf q_k = 0\\ \sum_{i=0}^{n-1} \mathbf q_k^i l_i(1) - \mathbf q_{k+1} = 0 \end{array}\right. $$

odeSPVI

Within odeSPVI I implemented a geometric integrator for Hamiltonian systems, like odeSE and odeSV, which uses spectral variational integrators with any order for polynomials of the basis and with any number of internal quadrature nodes, for both unidimensional and multidimensional problems.

[T,Y]=odeSPVI(@hamilton_equations,time,initial_conditions,options)
This solver just uses the stepper spectral_var_int.m which is the function where I implemented the resolution of the system displayed above; hamilton_equations must be a function_handle or an inline function where the user has to implement Hamilton's equations of motion: $$ \left\{\begin{array}{l} \dot{\mathbf q} = \frac{\partial H}{\partial \mathbf p}(\mathbf q,\mathbf p)\\ \dot{\mathbf p} = - \frac{\partial H}{\partial \mathbf q}(\mathbf q,\mathbf p) + \mathbf f(\mathbf q,\mathbf p)\end{array}\right.$$ where \(H(\mathbf q,\mathbf p)\) is the Hamiltonian of the system and \(\mathbf f(\mathbf q,\mathbf p)\) is an optional forcing term.

This below is the implementation of the stepper:

function [t_next,x_next,err]=spectral_var_int(f,t,x,dt,options)
  fsolve_opts = optimset('Jacobian','on');

  q_dofs = odeget(options,'Q_DoFs',[],'fast');
  p_dofs = odeget(options,'P_DoFs',[],'fast');
  nodes = odeget(options,'Nodes',[],'fast');
  weights = odeget(options,'Weights',[],'fast');
  leg = odeget(options,'Legendre',[],'fast');
  deriv = odeget(options,'Derivatives',[],'fast');
  extremes = odeget(options,'Extremes',[],'fast');

  dim = length(x)/2;
  N = q_dofs+p_dofs+2;
  q0 = x(1:dim)(:)';
  p0 = x(dim+1:end)(:)';

  SPVI = @(y)svi_system(y,q_dofs,p_dofs,weights,leg, ...
                deriv,extreme,dim,N,f,t,q0,p0,dt);
  y0 = [q0;zeros(q_dofs-1,dim);ones(p_dofs,1)*p0;q0;p0];
  y0 = reshape(y0,dim*N,1);
  z = fsolve(SPVI,y0,fsolve_opts);
  %z = inexact_newton(SVI,y0,new_opts) %new_opts must be set 
  z = reshape(z,N,dim);

  y = [leg*z(1:q_dofs,1:dim),z((q_dofs+1):(end-2),1:dim)];
  x_next = [y;z(end-1,:),z(end,:)]';
  t_next = [t+(dt/2)*(nodes+1), t+dt]';

  if(nargout==3)
    q_dofs_err = odeget(options,'Q_DoFs_err',[],'fast');
    p_dofs_err = odeget(options,'P_DoFs_err',[],'fast');
    nodes_err = odeget(options,'Nodes_err',[],'fast');
    weights_err = odeget(options,'Weights_err',[],'fast');
    leg_err = odeget(options,'Legendre_err',[],'fast');
    deriv_err = odeget(options,'Derivatives_err',[],'fast');
    extreme_err = odeget(options,'Extremes_err',[],'fast');

    N_err = q_dofs_err+p_dofs_err+2;
    q1 = x_next(1:dim,end)(:)';
    p1 = x_next(dim+1:end,end)(:)';

    SPVI = @(y)svi_system(y,q_dofs_err,p_dofs_err, ...
         weights_err,leg_err,deriv_err,extreme_err, ...
         dim,N_err,f,t,q0,p0,dt);
    p_interp = zeros(p_dofs_err,dim);

    p_lower_order = [p0;z(q_dofs+1:q_dofs+p_dofs,:);p1];
    for i=1:1:p_dofs_err
      p_interp(i,:) = .5*(p_lower_order(i,:) + ...
                            p_lower_order(i+1,:));
    end

    y0 = [z(1:q_dofs,:);zeros(1,dim);p_interp;q1;p1];
    y0 = reshape(y0,dim*N_err,1);
    z = fsolve(SPVI,y0,fsolve_opts);
    %z = inexact_newton(SVI,y0,new_opts) %new_opts must be set 
    z = reshape(z,N_err,dim);

    x_est = [z(end-1,:),z(end,:)]';
    err = norm(x_est-x_next(:,end),2);
  end
end
At lines 22 and 56 I put a comment line to show that it is possible to solve the implicit system also with inexact_newton (which is the one explained in the previous post on backward Euler), so it's not mandatory to use only fsolve. I will make it an option in order to let the user choose which solver to use.

Another important aspect to point out is that in case of an adaptive timestep the error is estimated using a new solution on the same timestep but with polynomials one order higher and one more internal node for the quadrature rule. Furthermore, for this last solution, the starting guess for fsolve is chosen in a non-trivial way: at line 53 we see that y0 has in the first q_dofs_err-1 rows the same modal values calculated before for the new solution x_next and a row of zeros just below. Then the starting nodal values for the momenta are set (lines 47-51) as a trivial average of new solution nodal values. This can seem wrong but I empirically stated that the number of iterations done by fsolve are the same of cases in which the reinitialization of nodal values is more sophisticated, so that is less computationally expensive to do a trivial average.

In the following code is shown the implementation of the system to be solved:

function [F,Jac] = svi_system(y,q_dofs,p_dofs,w,L,D,C, ...
          dim,N,f,t,q0,p0,dt)
  F = zeros(N*dim,1);
  V = zeros(p_dofs,dim*2);
  X = zeros(dim*2,1);
  W = reshape(y,N,dim);
  for i = 1:1:p_dofs
    X = [L(i,:)*W(1:q_dofs,:),W(i+q_dofs,:)]';
    V(i,:) = f(t,X);
  end
  for i = 1:1:dim
    F((1+N*(i-1)):(q_dofs+N*(i-1)),1) = (ones(1,p_dofs)* ...
     (((w.*y((q_dofs+1+N*(i-1)):(q_dofs+p_dofs+N*(i-1))))* ...
     ones(1,q_dofs)).*D + (((dt/2).*w.*V(:,i+dim))* ...
     ones(1,q_dofs)).*L) + (p0(i)*ones(1,q_dofs)).*C(1,:) - ...
      (y(N*i)*ones(1,q_dofs)).*C(2,:))';    
    F(1+N*(i-1)+q_dofs:N*(i-1)+q_dofs+p_dofs,1) = V(:,i) - ...
            (2.0/dt)*(D*y((1+N*(i-1)):(q_dofs+N*(i-1))));
    F(N*i-1) = C(2,:)*y((1+N*(i-1)):(q_dofs+N*(i-1)))-y(N*i-1);
    F(N*i) = C(1,:)*y((1+N*(i-1)):(q_dofs+N*(i-1))) - q0(i);
  end
  if(nargout>1)
    warning('off','Octave:broadcast');
    flag = 0;
    try
      [~,H]=f(0,zeros(2*dim,1));
    catch
      flag = 1;
      warning();
    end
    DV = zeros((dim*2)^2,p_dofs);
    if( flag==1 )
      for i = 1:1:p_dofs
        X = [L(i,:)*W(1:q_dofs,:),W(i+q_dofs,:)]';
        DV(:,i) = differential(f,t,X);
      end
    else
      for i = 1:1:p_dofs
        X = [L(i,:)*W(1:q_dofs,:),W(i+q_dofs,:)]';
        [~,DV(:,i)] = f(t,X);
      end
    end
    DV = DV';
    Jac=zeros(N*dim,N*dim);
    for u=1:1:dim
      [...]
    end
  end
end 
Here, at line 46, I don't show the implementation of the Jacobian of the implicit system because, with this visualization style, it may appear very chaotic. The important aspect is that the user can use the implemented Jacobian to speedup the computation. I stated that it really allows to speedup the computation when using fsolve as solver. Furthermore the code tries to see if the user has implemented, in the function defining the Hamilton's equations, also the Hessian of the Hamiltonian (which is needed int the computation of the Jacobian of the system). If the function passes as second parameter the Hessian, then the program uses that one and the computation is faster, otherwise it approximates the Hessian with the following function and the computation will be slower:
function Hamilt = differential(f,t,x)
  f_x = f(t,x);
  dim_f = length(f_x);
  dim_x = length(x);
  if( dim_f ~= dim_x )
    error('not implemented yet');
  end
  Hamilt = zeros(dim_f*dim_x,1);
  delta = sqrt(eps);
  for i = 1:1:dim_f
    for j = i:1:dim_x
      Hamilt(i+(j-1)*dim_f) = (1/delta)*(f(t,x+delta* ...
         [zeros(j-1,1);1;zeros(dim_x-j,1)])(i) - f_x(i));
      Hamilt(j+(i-1)*dim_x) = Hamilt(i+(j-1)*dim_f);
    end
  end
end

The Hessian of the Hamiltonian must be stored in a particular way (that I have to optimize yet, but the actual one works fine too) which is showed in the following example which is the definition of the Hamilton's equations for the armonic oscillator:

function [dy,ddy] = armonic_oscillator(t,y)
  dy = zeros(size(y));
  d = length(y)/2;
  q = y(1:d);
  p= y(d+1:end);

  dy(1:d) = p;
  dy(d+1:end) = -q;

  H = eye(2*d);
  ddy = transf(H);
end

function v = transf(H)
  [r,c] = size(H);
  v = zeros(r*c,1);
  for i=1:1:r
    for j=1:1:c
      v(i+(j-1)*r) = H(i,j);
    end
  end
end