OrdinaryDiffEqQPRK

Quadruple-precision parallel Runge-Kutta (QPRK) methods are high-order explicit solvers specifically designed for ultra-high precision computations using quad-precision arithmetic (Float128). These methods combine parallel evaluation capabilities with coefficients optimized for extended precision arithmetic. Note: These methods are still under-benchmarked and need more research.

Key Properties

QPRK methods provide:

  • Ultra-high-order accuracy (9th order) for maximum precision
  • Quadruple-precision optimization specifically designed for Float128
  • Parallel function evaluations for computational efficiency
  • Extreme precision capabilities for very demanding applications
  • Optimized coefficients for extended precision arithmetic

When to Use QPRK Methods

These methods are recommended for:

  • Ultra-high precision requirements demanding Float128 arithmetic
  • Extremely low tolerances (< 1e-20) where standard precision fails
  • Scientific applications requiring maximum possible accuracy
  • Parallel computing environments with quad-precision support
  • Research applications exploring limits of numerical precision
  • Long-time integration where error accumulation must be minimized to extreme levels

Important Requirements

Precision Requirements

  • Must use Float128 or higher precision number types
  • All problem components should support extended precision
  • Tolerances should match the precision capabilities (< 1e-20)

Computational Considerations

  • Slower than standard precision methods due to extended precision arithmetic
  • Higher memory usage due to extended precision
  • Limited hardware support for quad-precision operations

Mathematical Background

QPRK methods use tableaus with coefficients computed in extended precision to maintain accuracy throughout the ultra-high precision computation. The parallel structure allows independent function evaluations to be computed simultaneously.

Solver Selection Guide

Available methods

  • QPRK98: Ninth-order method optimized for quad-precision arithmetic with parallel evaluation

Usage guidelines

  • Essential to use Float128 for the state vector and parameters
  • Consider MultiFloats.jl for higher precision number types
  • Set very low tolerances (e.g., 1e-25) to utilize full precision
  • Test against alternatives like Feagin methods with BigFloat

Performance Considerations

  • Slower than standard precision methods due to extended precision arithmetic
  • Memory intensive due to extended precision storage
  • Hardware dependent - some architectures lack efficient quad-precision support

Alternative High-Precision Methods

For ultra-high precision, also consider:

  • Feagin methods with BigFloat for arbitrary precision
  • Arbitrary precision extrapolation methods
  • Verner methods with BigFloat for slightly lower but efficient precision
  • Taylor series methods with automatic differentiation for extreme precision

Usage Example

using OrdinaryDiffEqQPRK
# Ensure using Float128 for ultra-high precision
u0 = Float128[1.0, 0.0]  
tspan = (Float128(0.0), Float128(10.0))
prob = ODEProblem(f, u0, tspan)
sol = solve(prob, QPRK98(), abstol=1e-25, reltol=1e-25)

Installation

To be able to access the solvers in OrdinaryDiffEqQPRK, you must first install them use the Julia package manager:

using Pkg
Pkg.add("OrdinaryDiffEqQPRK")

This will only install the solvers listed at the bottom of this page. If you want to explore other solvers for your problem, you will need to install some of the other libraries listed in the navigation bar on the left.

Example usage

using OrdinaryDiffEqQPRK

function lorenz!(du, u, p, t)
    du[1] = 10.0 * (u[2] - u[1])
    du[2] = u[1] * (28.0 - u[3]) - u[2]
    du[3] = u[1] * u[2] - (8 / 3) * u[3]
end
u0 = [1.0; 0.0; 0.0]
tspan = (0.0, 100.0)
prob = ODEProblem(lorenz!, u0, tspan)
sol = solve(prob, QPRK98())

Full list of solvers

OrdinaryDiffEqQPRK.QPRK98Type
QPRK98(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
         step_limiter! = OrdinaryDiffEq.trivial_limiter!,
         thread = OrdinaryDiffEq.False())

Explicit Runge-Kutta Method. Runge–Kutta pairs of orders 9(8) for use in quadruple precision computations

Keyword Arguments

  • stage_limiter!: function of the form limiter!(u, integrator, p, t)
  • step_limiter!: function of the form limiter!(u, integrator, p, t)
  • thread: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()) or use multiple threads (thread = OrdinaryDiffEq.True()) when Julia is started with multiple threads.

References

Kovalnogov VN, Fedorov RV, Karpukhina TV, Simos TE, Tsitouras C. Runge–Kutta pairs of orders 9 (8) for use in quadruple precision computations. Numerical Algorithms, 2023. doi: https://doi.org/10.1007/s11075-023-01632-8

source