Using the an implicit solver based on Ariadne with Trixi.jl

using Trixi
using Theseus
using CairoMakie
using LinearAlgebra
import Ariadne: JacobianOperator

Notes: Must disable both Polyester and LoopVectorization for Enzyme to be able to differentiate Trixi.jl

LocalPreferences.jl

[Trixi]
loop_vectorization = false
backend = "static"
@assert Trixi._PREFERENCE_THREADING !== :polyester
@assert !Trixi._PREFERENCE_LOOPVECTORIZATION

Load Trixi Example

trixi_include(@__MODULE__, joinpath(examples_dir(), "tree_2d_dgsem", "elixir_advection_basic.jl"), sol = nothing);

u = copy(ode.u0)
du = zero(ode.u0)
res = zero(ode.u0)

F! = Theseus.nonlinear_problem(Theseus.ImplicitEuler(), ode.f)
J = JacobianOperator(F!, res, u, (ode.u0, 1.0, du, ode.p, 0.0, (), 1))

out = zero(u)
v = zero(u)

# precompile
mul!(u, J, v)
F!(res, u, (ode.u0, 1.0, du, ode.p, 0.0, (), 1))

@time mul!(u, J, v)
@time F!(res, u, (ode.u0, 1.0, du, ode.p, 0.0, (), 1))
[ Info: You just called `trixi_include`. Julia may now compile the code, please be patient.
  0.000197 seconds (83 allocations: 8.906 KiB)
  0.000141 seconds (7 allocations: 352 bytes)

Cost of time(mul!) ≈ 2 * time(F!)

Solve using ODE interface

sol_trbdf2 = solve(
    ode, Theseus.TRBDF2();
    dt = 1.0, # solve needs some value here but it will be overwritten by the stepsize_callback
    ode_default_options()..., callback = callbacks,
    # verbose=1,
    krylov_algo = :gmres,
    # krylov_kwargs=(;verbose=1)
);

████████╗██████╗ ██╗██╗  ██╗██╗
╚══██╔══╝██╔══██╗██║╚██╗██╔╝██║
   ██║   ██████╔╝██║ ╚███╔╝ ██║
   ██║   ██╔══██╗██║ ██╔██╗ ██║
   ██║   ██║  ██║██║██╔╝ ██╗██║
   ╚═╝   ╚═╝  ╚═╝╚═╝╚═╝  ╚═╝╚═╝

┌──────────────────────────────────────────────────────────────────────────────────────────────────┐
│ SemidiscretizationHyperbolic                                                                     │
│ ════════════════════════════                                                                     │
│ #spatial dimensions: ………………………… 2                                                                │
│ mesh: ………………………………………………………………… TreeMesh{2, Trixi.SerialTree{2, Float64}} with length 341        │
│ equations: …………………………………………………… LinearScalarAdvectionEquation2D                                  │
│ initial condition: ……………………………… initial_condition_convergence_test                               │
│ boundary conditions: ………………………… Trixi.BoundaryConditionPeriodic                                  │
│ source terms: …………………………………………… nothing                                                          │
│ solver: …………………………………………………………… DG                                                               │
│ total #DOFs per field: …………………… 4096                                                             │
└──────────────────────────────────────────────────────────────────────────────────────────────────┘

┌──────────────────────────────────────────────────────────────────────────────────────────────────┐
│ TreeMesh{2, Trixi.SerialTree{2, Float64}}                                                        │
│ ═════════════════════════════════════════                                                        │
│ center: …………………………………………………………… [0.0, 0.0]                                                       │
│ length: …………………………………………………………… 2.0                                                              │
│ periodicity: ……………………………………………… (true, true)                                                     │
│ current #cells: ……………………………………… 341                                                              │
│ #leaf-cells: ……………………………………………… 256                                                              │
│ maximum #cells: ……………………………………… 30000                                                            │
└──────────────────────────────────────────────────────────────────────────────────────────────────┘

┌──────────────────────────────────────────────────────────────────────────────────────────────────┐
│ LinearScalarAdvectionEquation2D                                                                  │
│ ═══════════════════════════════                                                                  │
│ #variables: ………………………………………………… 1                                                                │
│ │ variable 1: …………………………………………… scalar                                                           │
└──────────────────────────────────────────────────────────────────────────────────────────────────┘

┌──────────────────────────────────────────────────────────────────────────────────────────────────┐
│ DG{Float64}                                                                                      │
│ ═══════════                                                                                      │
│ basis: ……………………………………………………………… LobattoLegendreBasis{Float64}(polydeg=3)                         │
│ mortar: …………………………………………………………… LobattoLegendreMortarL2{Float64}(polydeg=3)                      │
│ surface integral: ………………………………… SurfaceIntegralWeakForm                                          │
│ │ surface flux: ……………………………………… FluxLaxFriedrichs(max_abs_speed)                                 │
│ volume integral: …………………………………… VolumeIntegralWeakForm                                           │
└──────────────────────────────────────────────────────────────────────────────────────────────────┘

┌──────────────────────────────────────────────────────────────────────────────────────────────────┐
│ AnalysisCallback                                                                                 │
│ ════════════════                                                                                 │
│ interval: ……………………………………………………… 100                                                              │
│ analyzer: ……………………………………………………… LobattoLegendreAnalyzer{Float64}(polydeg=6)                      │
│ │ error 1: …………………………………………………… l2_error                                                         │
│ │ error 2: …………………………………………………… linf_error                                                       │
│ │ integral 1: …………………………………………… entropy_timederivative                                           │
│ save analysis to file: …………………… no                                                               │
└──────────────────────────────────────────────────────────────────────────────────────────────────┘

┌──────────────────────────────────────────────────────────────────────────────────────────────────┐
│ SaveSolutionCallback                                                                             │
│ ════════════════════                                                                             │
│ interval: ……………………………………………………… 100                                                              │
│ solution variables: …………………………… cons2prim                                                        │
│ save initial solution: …………………… yes                                                              │
│ save final solution: ………………………… yes                                                              │
│ output directory: ………………………………… /home/runner/work/Ariadne.jl/Ariadne.jl/docs/build/generated/out │
└──────────────────────────────────────────────────────────────────────────────────────────────────┘

┌──────────────────────────────────────────────────────────────────────────────────────────────────┐
│ StepsizeCallback                                                                                 │
│ ════════════════                                                                                 │
│ CFL Advective: ………………………………………… Returns{Float64}(1.6)                                            │
│ CFL Diffusive: ………………………………………… Returns{Float64}(0.0)                                            │
│ Interval: ……………………………………………………… 1                                                                │
└──────────────────────────────────────────────────────────────────────────────────────────────────┘

┌──────────────────────────────────────────────────────────────────────────────────────────────────┐
│ Time integration                                                                                 │
│ ════════════════                                                                                 │
│ Start time: ………………………………………………… 0.0                                                              │
│ Final time: ………………………………………………… 1.0                                                              │
│ time integrator: …………………………………… TRBDF2                                                           │
│ adaptive: ……………………………………………………… false                                                            │
└──────────────────────────────────────────────────────────────────────────────────────────────────┘
┌──────────────────────────────────────────────────────────────────────────────────────────────────┐
│ Environment information                                                                          │
│ ═══════════════════════                                                                          │
│ #threads: ……………………………………………………… 1                                                                │
│ threading backend: ……………………………… static                                                           │
│ LoopVectorization: ……………………………… disabled                                                         │
└──────────────────────────────────────────────────────────────────────────────────────────────────┘

────────────────────────────────────────────────────────────────────────────────────────────────────
 Simulation running 'LinearScalarAdvectionEquation2D' with DGSEM(polydeg=3)
────────────────────────────────────────────────────────────────────────────────────────────────────
 #timesteps:                  0                run time:       1.08795400e+00 s
 Δt:             1.00000000e+00                └── GC time:    0.00000000e+00 s (0.000%)
 sim. time:      0.00000000e+00 (0.000%)       time/DOF/rhs!:  2.56275635e-08 s
                                               PID:                   Inf s
 #DOFs per field:          4096                alloc'd memory:        840.995 MiB
 #elements:                 256

 Variable:       scalar
 L2 error:       3.94120825e-06
 Linf error:     1.23670425e-05
 ∑∂S/∂U ⋅ Uₜ :  -2.19124035e-16
────────────────────────────────────────────────────────────────────────────────────────────────────

┌ Warning: TODO forward zero-set of memorycopy used memset rather than runtime type
│ Caused by:
│ Stacktrace:
│   [1] copy
│     @ ./array.jl:350
│   [2] unaliascopy
│     @ ./abstractarray.jl:1516
│   [3] unalias
│     @ ./abstractarray.jl:1500
│   [4] broadcast_unalias
│     @ ./broadcast.jl:946
│   [5] preprocess
│     @ ./broadcast.jl:953
│   [6] preprocess_args
│     @ ./broadcast.jl:955
│   [7] preprocess
│     @ ./broadcast.jl:952
│   [8] override_bc_copyto!
│     @ ~/.julia/packages/Enzyme/C5rSn/src/compiler/interpreter.jl:818
│   [9] copyto!
│     @ ./broadcast.jl:925
│  [10] materialize!
│     @ ./broadcast.jl:883
│  [11] materialize!
│     @ ./broadcast.jl:880
│  [12] TRBDF2
│     @ ~/work/Ariadne.jl/Ariadne.jl/libs/Theseus/src/Theseus.jl:178
└ @ Enzyme.Compiler ~/.julia/packages/Enzyme/C5rSn/src/rules/llvmrules.jl:812

────────────────────────────────────────────────────────────────────────────────────────────────────
 Simulation running 'LinearScalarAdvectionEquation2D' with DGSEM(polydeg=3)
────────────────────────────────────────────────────────────────────────────────────────────────────
 #timesteps:                 18                run time:       3.65153117e+01 s
 Δt:             5.55555556e-02                └── GC time:    5.05735319e-01 s (1.385%)
 sim. time:      1.00000000e+00 (100.000%)     time/DOF/rhs!:  2.11327499e-08 s
                                               PID:            1.61308350e-05 s
 #DOFs per field:          4096                alloc'd memory:        556.904 MiB
 #elements:                 256

 Variable:       scalar
 L2 error:       1.71224434e-04
 Linf error:     2.52822142e-04
 ∑∂S/∂U ⋅ Uₜ :  -1.97551223e-09
────────────────────────────────────────────────────────────────────────────────────────────────────

────────────────────────────────────────────────────────────────────────────────────
             Trixi.jl                      Time                    Allocations
                                  ───────────────────────   ────────────────────────
        Tot / % measured:              36.8s /   4.3%           1.95GiB /   1.8%

Section                   ncalls     time    %tot     avg     alloc    %tot      avg
────────────────────────────────────────────────────────────────────────────────────
I/O                            3    1.32s   83.7%   440ms   16.9MiB   47.6%  5.63MiB
  ~I/O~                        3    1.14s   72.0%   379ms   12.6MiB   35.5%  4.20MiB
  save solution                2    178ms   11.2%  88.8ms   3.68MiB   10.4%  1.84MiB
  get element variables        2   6.13ms    0.4%  3.07ms    620KiB    1.7%   310KiB
  get node variables           2   1.00μs    0.0%   501ns     0.00B    0.0%    0.00B
  save mesh                    2    521ns    0.0%   260ns     0.00B    0.0%    0.00B
analyze solution               2    212ms   13.4%   106ms   18.6MiB   52.4%  9.31MiB
rhs!                         533   46.0ms    2.9%  86.4μs   4.78KiB    0.0%    9.19B
  volume integral            533   30.8ms    2.0%  57.8μs     0.00B    0.0%    0.00B
  interface flux             533   4.51ms    0.3%  8.46μs     0.00B    0.0%    0.00B
  surface integral           533   3.40ms    0.2%  6.38μs     0.00B    0.0%    0.00B
  prolong2interfaces         533   3.28ms    0.2%  6.15μs     0.00B    0.0%    0.00B
  Jacobian                   533   2.82ms    0.2%  5.29μs     0.00B    0.0%    0.00B
  reset ∂u/∂t                533    547μs    0.0%  1.03μs     0.00B    0.0%    0.00B
  ~rhs!~                     533    535μs    0.0%  1.00μs   4.78KiB    0.0%    9.19B
  prolong2mortars            533   44.8μs    0.0%  84.1ns     0.00B    0.0%    0.00B
  prolong2boundaries         533   34.4μs    0.0%  64.6ns     0.00B    0.0%    0.00B
  mortar flux                533   31.5μs    0.0%  59.0ns     0.00B    0.0%    0.00B
  source terms               533   15.7μs    0.0%  29.4ns     0.00B    0.0%    0.00B
  boundary flux              533   15.6μs    0.0%  29.3ns     0.00B    0.0%    0.00B
calculate dt                  19   15.5μs    0.0%   818ns   2.67KiB    0.0%     144B
────────────────────────────────────────────────────────────────────────────────────

Plot the solution

We have to manually convert the sol since Theseus has it's own lightweight solution type.

plot(Trixi.PlotData2DTriangulated(sol_trbdf2.u[end], sol_trbdf2.prob.p))

Solve using OrdinaryDiffEqSDIRK

import OrdinaryDiffEqSDIRK
import DifferentiationInterface: AutoFiniteDiff
sol_sdrik = solve(
    ode, OrdinaryDiffEqSDIRK.TRBDF2(autodiff = AutoFiniteDiff());
    dt = 1.0, # solve needs some value here but it will be overwritten by the stepsize_callback
    ode_default_options()..., callback = callbacks,
    adaptive = false
);

████████╗██████╗ ██╗██╗  ██╗██╗
╚══██╔══╝██╔══██╗██║╚██╗██╔╝██║
   ██║   ██████╔╝██║ ╚███╔╝ ██║
   ██║   ██╔══██╗██║ ██╔██╗ ██║
   ██║   ██║  ██║██║██╔╝ ██╗██║
   ╚═╝   ╚═╝  ╚═╝╚═╝╚═╝  ╚═╝╚═╝

┌──────────────────────────────────────────────────────────────────────────────────────────────────┐
│ SemidiscretizationHyperbolic                                                                     │
│ ════════════════════════════                                                                     │
│ #spatial dimensions: ………………………… 2                                                                │
│ mesh: ………………………………………………………………… TreeMesh{2, Trixi.SerialTree{2, Float64}} with length 341        │
│ equations: …………………………………………………… LinearScalarAdvectionEquation2D                                  │
│ initial condition: ……………………………… initial_condition_convergence_test                               │
│ boundary conditions: ………………………… Trixi.BoundaryConditionPeriodic                                  │
│ source terms: …………………………………………… nothing                                                          │
│ solver: …………………………………………………………… DG                                                               │
│ total #DOFs per field: …………………… 4096                                                             │
└──────────────────────────────────────────────────────────────────────────────────────────────────┘

┌──────────────────────────────────────────────────────────────────────────────────────────────────┐
│ TreeMesh{2, Trixi.SerialTree{2, Float64}}                                                        │
│ ═════════════════════════════════════════                                                        │
│ center: …………………………………………………………… [0.0, 0.0]                                                       │
│ length: …………………………………………………………… 2.0                                                              │
│ periodicity: ……………………………………………… (true, true)                                                     │
│ current #cells: ……………………………………… 341                                                              │
│ #leaf-cells: ……………………………………………… 256                                                              │
│ maximum #cells: ……………………………………… 30000                                                            │
└──────────────────────────────────────────────────────────────────────────────────────────────────┘

┌──────────────────────────────────────────────────────────────────────────────────────────────────┐
│ LinearScalarAdvectionEquation2D                                                                  │
│ ═══════════════════════════════                                                                  │
│ #variables: ………………………………………………… 1                                                                │
│ │ variable 1: …………………………………………… scalar                                                           │
└──────────────────────────────────────────────────────────────────────────────────────────────────┘

┌──────────────────────────────────────────────────────────────────────────────────────────────────┐
│ DG{Float64}                                                                                      │
│ ═══════════                                                                                      │
│ basis: ……………………………………………………………… LobattoLegendreBasis{Float64}(polydeg=3)                         │
│ mortar: …………………………………………………………… LobattoLegendreMortarL2{Float64}(polydeg=3)                      │
│ surface integral: ………………………………… SurfaceIntegralWeakForm                                          │
│ │ surface flux: ……………………………………… FluxLaxFriedrichs(max_abs_speed)                                 │
│ volume integral: …………………………………… VolumeIntegralWeakForm                                           │
└──────────────────────────────────────────────────────────────────────────────────────────────────┘

┌──────────────────────────────────────────────────────────────────────────────────────────────────┐
│ AnalysisCallback                                                                                 │
│ ════════════════                                                                                 │
│ interval: ……………………………………………………… 100                                                              │
│ analyzer: ……………………………………………………… LobattoLegendreAnalyzer{Float64}(polydeg=6)                      │
│ │ error 1: …………………………………………………… l2_error                                                         │
│ │ error 2: …………………………………………………… linf_error                                                       │
│ │ integral 1: …………………………………………… entropy_timederivative                                           │
│ save analysis to file: …………………… no                                                               │
└──────────────────────────────────────────────────────────────────────────────────────────────────┘

┌──────────────────────────────────────────────────────────────────────────────────────────────────┐
│ SaveSolutionCallback                                                                             │
│ ════════════════════                                                                             │
│ interval: ……………………………………………………… 100                                                              │
│ solution variables: …………………………… cons2prim                                                        │
│ save initial solution: …………………… yes                                                              │
│ save final solution: ………………………… yes                                                              │
│ output directory: ………………………………… /home/runner/work/Ariadne.jl/Ariadne.jl/docs/build/generated/out │
└──────────────────────────────────────────────────────────────────────────────────────────────────┘

┌──────────────────────────────────────────────────────────────────────────────────────────────────┐
│ StepsizeCallback                                                                                 │
│ ════════════════                                                                                 │
│ CFL Advective: ………………………………………… Returns{Float64}(1.6)                                            │
│ CFL Diffusive: ………………………………………… Returns{Float64}(0.0)                                            │
│ Interval: ……………………………………………………… 1                                                                │
└──────────────────────────────────────────────────────────────────────────────────────────────────┘

┌──────────────────────────────────────────────────────────────────────────────────────────────────┐
│ Time integration                                                                                 │
│ ════════════════                                                                                 │
│ Start time: ………………………………………………… 0.0                                                              │
│ Final time: ………………………………………………… 1.0                                                              │
│ time integrator: …………………………………… TRBDF2                                                           │
│ adaptive: ……………………………………………………… false                                                            │
└──────────────────────────────────────────────────────────────────────────────────────────────────┘
┌──────────────────────────────────────────────────────────────────────────────────────────────────┐
│ Environment information                                                                          │
│ ═══════════════════════                                                                          │
│ #threads: ……………………………………………………… 1                                                                │
│ threading backend: ……………………………… static                                                           │
│ LoopVectorization: ……………………………… disabled                                                         │
└──────────────────────────────────────────────────────────────────────────────────────────────────┘

────────────────────────────────────────────────────────────────────────────────────────────────────
 Simulation running 'LinearScalarAdvectionEquation2D' with DGSEM(polydeg=3)
────────────────────────────────────────────────────────────────────────────────────────────────────
 #timesteps:                  0                run time:       8.62000000e-07 s
 Δt:             1.00000000e+00                └── GC time:    0.00000000e+00 s (0.000%)
 sim. time:      0.00000000e+00 (0.000%)       time/DOF/rhs!:  2.33884277e-08 s
                                               PID:                   Inf s
 #DOFs per field:          4096                alloc'd memory:        816.715 MiB
 #elements:                 256

 Variable:       scalar
 L2 error:       3.94120825e-06
 Linf error:     1.23670425e-05
 ∑∂S/∂U ⋅ Uₜ :  -2.19124035e-16
────────────────────────────────────────────────────────────────────────────────────────────────────


────────────────────────────────────────────────────────────────────────────────────────────────────
 Simulation running 'LinearScalarAdvectionEquation2D' with DGSEM(polydeg=3)
────────────────────────────────────────────────────────────────────────────────────────────────────
 #timesteps:                 18                run time:       3.90657590e+01 s
 Δt:             5.55555556e-02                └── GC time:    1.04686172e+00 s (2.680%)
 sim. time:      1.00000000e+00 (100.000%)     time/DOF/rhs!:  2.35496357e-08 s
                                               PID:            6.46447422e-08 s
 #DOFs per field:          4096                alloc'd memory:       1126.916 MiB
 #elements:                 256

 Variable:       scalar
 L2 error:       1.71161806e-04
 Linf error:     2.52798215e-04
 ∑∂S/∂U ⋅ Uₜ :  -1.97752414e-09
────────────────────────────────────────────────────────────────────────────────────────────────────

────────────────────────────────────────────────────────────────────────────────────
             Trixi.jl                      Time                    Allocations
                                  ───────────────────────   ────────────────────────
        Tot / % measured:              39.1s /  36.4%           4.74GiB /   0.1%

Section                   ncalls     time    %tot     avg     alloc    %tot      avg
────────────────────────────────────────────────────────────────────────────────────
rhs!                        148k    14.2s   99.9%  96.4μs   4.78KiB    0.1%    0.03B
  volume integral           148k    9.71s   68.3%  65.8μs     0.00B    0.0%    0.00B
  interface flux            148k    1.17s    8.2%  7.93μs     0.00B    0.0%    0.00B
  surface integral          148k    1.12s    7.8%  7.57μs     0.00B    0.0%    0.00B
  Jacobian                  148k    1.00s    7.0%  6.79μs     0.00B    0.0%    0.00B
  prolong2interfaces        148k    963ms    6.8%  6.52μs     0.00B    0.0%    0.00B
  ~rhs!~                    148k    125ms    0.9%   845ns   4.78KiB    0.1%    0.03B
  reset ∂u/∂t               148k    100ms    0.7%   679ns     0.00B    0.0%    0.00B
  prolong2mortars           148k   6.94ms    0.0%  47.0ns     0.00B    0.0%    0.00B
  mortar flux               148k   6.52ms    0.0%  44.2ns     0.00B    0.0%    0.00B
  prolong2boundaries        148k   5.61ms    0.0%  38.0ns     0.00B    0.0%    0.00B
  boundary flux             148k   4.38ms    0.0%  29.7ns     0.00B    0.0%    0.00B
  source terms              148k   4.36ms    0.0%  29.6ns     0.00B    0.0%    0.00B
analyze solution               2   9.66ms    0.1%  4.83ms   6.06MiB   96.7%  3.03MiB
I/O                            3   2.19ms    0.0%   731μs    203KiB    3.2%  67.8KiB
  save solution                2   2.19ms    0.0%  1.09ms    200KiB    3.1%   100KiB
  ~I/O~                        3   4.35μs    0.0%  1.45μs   1.70KiB    0.0%     581B
  get element variables        2   1.91μs    0.0%   956ns   1.41KiB    0.0%     720B
  get node variables           2   1.15μs    0.0%   576ns     0.00B    0.0%    0.00B
  save mesh                    2    100ns    0.0%  50.0ns     0.00B    0.0%    0.00B
calculate dt                  19   63.7μs    0.0%  3.35μs   2.67KiB    0.0%     144B
────────────────────────────────────────────────────────────────────────────────────

Plot the solution

plot(Trixi.PlotData2DTriangulated(sol_sdrik.u[end], sol_sdrik.prob.p))

Increase CFL numbers

trixi_include(@__MODULE__, joinpath(examples_dir(), "tree_2d_dgsem", "elixir_advection_basic.jl"), cfl = 10, sol = nothing);

sol = solve(
    ode, Theseus.ImplicitEuler();
    dt = 1.0, # solve needs some value here but it will be overwritten by the stepsize_callback
    ode_default_options()..., callback = callbacks,
    # verbose=1,
    krylov_algo = :gmres,
    # krylov_kwargs=(;verbose=1)
);

@show callbacks.discrete_callbacks[4]
┌──────────────────────────────────────────────────────────────────────────────────────────────────┐
│ StepsizeCallback                                                                                 │
│ ════════════════                                                                                 │
│ CFL Advective: ………………………………………… Returns{Int64}(10)                                               │
│ CFL Diffusive: ………………………………………… Returns{Float64}(0.0)                                            │
│ Interval: ……………………………………………………… 1                                                                │
└──────────────────────────────────────────────────────────────────────────────────────────────────┘

Plot the solution

plot(Trixi.PlotData2DTriangulated(sol.u[end], sol.prob.p))

This page was generated using Literate.jl.