Using the an implicit solver based on Ariadne with Trixi.jl
using Trixi
using Theseus
using CairoMakie
using LinearAlgebra
import Ariadne: JacobianOperatorNotes: 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_LOOPVECTORIZATIONLoad 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.000213 seconds (86 allocations: 8.844 KiB)
0.000132 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.19077413e+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.76193848e-08 s
PID: Inf s
#DOFs per field: 4096 alloc'd memory: 545.031 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.34022518e+01 s
Δt: 5.55555556e-02 └── GC time: 2.49159323e-01 s (0.746%)
sim. time: 1.00000000e+00 (100.000%) time/DOF/rhs!: 2.10988975e-08 s
PID: 1.46577978e-05 s
#DOFs per field: 4096 alloc'd memory: 486.662 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: 33.7s / 7.2% 1.21GiB / 3.8%
Section ncalls time %tot avg alloc %tot avg
────────────────────────────────────────────────────────────────────────────────────
I/O 3 2.16s 89.3% 719ms 28.8MiB 60.7% 9.60MiB
~I/O~ 3 1.67s 69.2% 557ms 17.2MiB 36.3% 5.73MiB
save solution 2 479ms 19.8% 239ms 11.0MiB 23.2% 5.50MiB
get element variables 2 6.67ms 0.3% 3.34ms 620KiB 1.3% 310KiB
get node variables 2 1.05μs 0.0% 526ns 0.00B 0.0% 0.00B
save mesh 2 430ns 0.0% 215ns 0.00B 0.0% 0.00B
analyze solution 2 212ms 8.8% 106ms 18.6MiB 39.3% 9.31MiB
rhs! 533 46.0ms 1.9% 86.2μs 4.78KiB 0.0% 9.19B
volume integral 533 30.6ms 1.3% 57.5μs 0.00B 0.0% 0.00B
interface flux 533 4.40ms 0.2% 8.26μs 0.00B 0.0% 0.00B
surface integral 533 3.49ms 0.1% 6.55μs 0.00B 0.0% 0.00B
prolong2interfaces 533 3.33ms 0.1% 6.25μs 0.00B 0.0% 0.00B
Jacobian 533 2.80ms 0.1% 5.26μs 0.00B 0.0% 0.00B
reset ∂u/∂t 533 591μs 0.0% 1.11μs 0.00B 0.0% 0.00B
~rhs!~ 533 561μs 0.0% 1.05μs 4.78KiB 0.0% 9.19B
prolong2mortars 533 42.4μs 0.0% 79.6ns 0.00B 0.0% 0.00B
prolong2boundaries 533 31.8μs 0.0% 59.7ns 0.00B 0.0% 0.00B
mortar flux 533 29.3μs 0.0% 54.9ns 0.00B 0.0% 0.00B
boundary flux 533 15.7μs 0.0% 29.5ns 0.00B 0.0% 0.00B
source terms 533 15.5μs 0.0% 29.1ns 0.00B 0.0% 0.00B
calculate dt 19 17.3μs 0.0% 910ns 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: 6.11000000e-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.34567871e-08 s
PID: Inf s
#DOFs per field: 4096 alloc'd memory: 720.458 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.89417684e+01 s
Δt: 5.55555556e-02 └── GC time: 8.31005477e-01 s (2.134%)
sim. time: 1.00000000e+00 (100.000%) time/DOF/rhs!: 2.35538028e-08 s
PID: 6.44394976e-08 s
#DOFs per field: 4096 alloc'd memory: 909.094 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: 38.9s / 36.5% 4.62GiB / 0.1%
Section ncalls time %tot avg alloc %tot avg
────────────────────────────────────────────────────────────────────────────────────
rhs! 148k 14.2s 100.0% 96.4μs 4.78KiB 0.1% 0.03B
volume integral 148k 9.72s 68.3% 65.9μs 0.00B 0.0% 0.00B
interface flux 148k 1.15s 8.1% 7.80μs 0.00B 0.0% 0.00B
surface integral 148k 1.11s 7.8% 7.55μs 0.00B 0.0% 0.00B
Jacobian 148k 995ms 7.0% 6.75μs 0.00B 0.0% 0.00B
prolong2interfaces 148k 982ms 6.9% 6.65μs 0.00B 0.0% 0.00B
~rhs!~ 148k 126ms 0.9% 853ns 4.78KiB 0.1% 0.03B
reset ∂u/∂t 148k 103ms 0.7% 697ns 0.00B 0.0% 0.00B
prolong2mortars 148k 6.54ms 0.0% 44.3ns 0.00B 0.0% 0.00B
prolong2boundaries 148k 6.33ms 0.0% 42.9ns 0.00B 0.0% 0.00B
mortar flux 148k 5.50ms 0.0% 37.3ns 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.32ms 0.0% 29.3ns 0.00B 0.0% 0.00B
I/O 3 2.15ms 0.0% 717μs 205KiB 3.9% 68.2KiB
save solution 2 2.14ms 0.0% 1.07ms 202KiB 3.9% 101KiB
~I/O~ 3 4.30μs 0.0% 1.43μs 1.70KiB 0.0% 581B
get element variables 2 1.37μs 0.0% 686ns 1.41KiB 0.0% 720B
get node variables 2 1.10μs 0.0% 551ns 0.00B 0.0% 0.00B
save mesh 2 131ns 0.0% 65.5ns 0.00B 0.0% 0.00B
analyze solution 2 1.76ms 0.0% 878μs 4.87MiB 95.9% 2.43MiB
calculate dt 19 71.8μs 0.0% 3.78μs 2.67KiB 0.1% 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.