Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

WIP: In-situ visualization #1: Make Makie available in VisualizationCallback #2199

Draft
wants to merge 7 commits into
base: main
Choose a base branch
from
Draft
3 changes: 2 additions & 1 deletion examples/tree_2d_dgsem/elixir_advection_amr_visualization.jl
Original file line number Diff line number Diff line change
Expand Up @@ -2,6 +2,7 @@
using OrdinaryDiffEq
using Trixi
using Plots
using GLMakie

###############################################################################
# semidiscretization of the linear advection equation
Expand Down Expand Up @@ -50,7 +51,7 @@ save_solution = SaveSolutionCallback(interval = 100,

# Enable in-situ visualization with a new plot generated every 20 time steps
# and additional plotting options passed as keyword arguments
visualization = VisualizationCallback(interval = 20, clims = (0, 1))
visualization = VisualizationCallback(interval = 20, clims = (0, 1), plot_creator = Trixi.show_plot_makie)

amr_controller = ControllerThreeLevel(semi, IndicatorMax(semi, variable = first),
base_level = 3,
Expand Down
109 changes: 109 additions & 0 deletions examples/tree_2d_dgsem/elixir_euler_blast_wave_amr_visualization.jl
Original file line number Diff line number Diff line change
@@ -0,0 +1,109 @@

using OrdinaryDiffEq
using Trixi
using GLMakie

###############################################################################
# semidiscretization of the compressible Euler equations

equations = CompressibleEulerEquations2D(1.4)

"""
initial_condition_blast_wave(x, t, equations::CompressibleEulerEquations2D)

A medium blast wave taken from
- Sebastian Hennemann, Gregor J. Gassner (2020)
A provably entropy stable subcell shock capturing approach for high order split form DG
[arXiv: 2008.12044](https://arxiv.org/abs/2008.12044)
"""
function initial_condition_blast_wave(x, t, equations::CompressibleEulerEquations2D)
# Modified From Hennemann & Gassner JCP paper 2020 (Sec. 6.3) -> "medium blast wave"
# Set up polar coordinates
inicenter = SVector(0.0, 0.0)
x_norm = x[1] - inicenter[1]
y_norm = x[2] - inicenter[2]
r = sqrt(x_norm^2 + y_norm^2)
phi = atan(y_norm, x_norm)
sin_phi, cos_phi = sincos(phi)

# Calculate primitive variables
rho = r > 0.5 ? 1.0 : 1.1691
v1 = r > 0.5 ? 0.0 : 0.1882 * cos_phi
v2 = r > 0.5 ? 0.0 : 0.1882 * sin_phi
p = r > 0.5 ? 1.0E-3 : 1.245

return prim2cons(SVector(rho, v1, v2, p), equations)
end
initial_condition = initial_condition_blast_wave

surface_flux = flux_lax_friedrichs
volume_flux = flux_ranocha
basis = LobattoLegendreBasis(3)
indicator_sc = IndicatorHennemannGassner(equations, basis,
alpha_max = 0.5,
alpha_min = 0.001,
alpha_smooth = true,
variable = density_pressure)
volume_integral = VolumeIntegralShockCapturingHG(indicator_sc;
volume_flux_dg = volume_flux,
volume_flux_fv = surface_flux)
solver = DGSEM(basis, surface_flux, volume_integral)

coordinates_min = (-2.0, -2.0)
coordinates_max = (2.0, 2.0)
mesh = TreeMesh(coordinates_min, coordinates_max,
initial_refinement_level = 6,
n_cells_max = 10_000)

semi = SemidiscretizationHyperbolic(mesh, equations, initial_condition, solver)

###############################################################################
# ODE solvers, callbacks etc.

tspan = (0.0, 12.5)
ode = semidiscretize(semi, tspan)

summary_callback = SummaryCallback()

analysis_interval = 100
analysis_callback = AnalysisCallback(semi, interval = analysis_interval)

alive_callback = AliveCallback(analysis_interval = analysis_interval)

save_solution = SaveSolutionCallback(interval = 100,
save_initial_solution = true,
save_final_solution = true,
solution_variables = cons2prim)

amr_indicator = IndicatorHennemannGassner(semi,
alpha_max = 0.5,
alpha_min = 0.001,
alpha_smooth = true,
variable = density_pressure)


visualization = VisualizationCallback(interval = 20, clims = (0, 1), plot_creator = Trixi.show_plot_makie)

amr_controller = ControllerThreeLevel(semi, amr_indicator,
base_level = 4,
max_level = 6, max_threshold = 0.01)
amr_callback = AMRCallback(semi, amr_controller,
interval = 5,
adapt_initial_condition = true,
adapt_initial_condition_only_refine = true)

stepsize_callback = StepsizeCallback(cfl = 0.9)

callbacks = CallbackSet(summary_callback,
analysis_callback, alive_callback,
save_solution,
visualization,
amr_callback, stepsize_callback)

###############################################################################
# run the simulation

sol = solve(ode, CarpenterKennedy2N54(williamson_condition = false),
dt = 1.0, # solve needs some value here but it will be overwritten by the stepsize_callback
save_everystep = false, callback = callbacks);
summary_callback() # print the timer summary
71 changes: 71 additions & 0 deletions examples/tree_3d_dgsem/elixir_advection_amr_visualization.jl
Original file line number Diff line number Diff line change
@@ -0,0 +1,71 @@

using OrdinaryDiffEq
using Trixi
using Plots
using GLMakie

###############################################################################
# semidiscretization of the linear advection equation

advection_velocity = (0.2, -0.7, 0.5)
equations = LinearScalarAdvectionEquation3D(advection_velocity)

initial_condition = initial_condition_gauss
solver = DGSEM(polydeg = 3, surface_flux = flux_lax_friedrichs)

coordinates_min = (-5.0, -5.0, -5.0)
coordinates_max = (5.0, 5.0, 5.0)
mesh = TreeMesh(coordinates_min, coordinates_max,
initial_refinement_level = 4,
n_cells_max = 30_000)

semi = SemidiscretizationHyperbolic(mesh, equations, initial_condition, solver)

###############################################################################
# ODE solvers, callbacks etc.

tspan = (0.0, 0.3)
ode = semidiscretize(semi, tspan)

summary_callback = SummaryCallback()

analysis_interval = 100
analysis_callback = AnalysisCallback(semi, interval = analysis_interval,
extra_analysis_integrals = (entropy,))

alive_callback = AliveCallback(analysis_interval = analysis_interval)

save_solution = SaveSolutionCallback(interval = 100,
save_initial_solution = true,
save_final_solution = true,
solution_variables = cons2prim)


visualization = VisualizationCallback(interval = 20, clims = (0, 1), plot_data_creator = Trixi.PlotData3D, plot_creator = Trixi.show_plot_makie)

amr_controller = ControllerThreeLevel(semi, IndicatorMax(semi, variable = first),
base_level = 4,
med_level = 5, med_threshold = 0.1,
max_level = 6, max_threshold = 0.6)
amr_callback = AMRCallback(semi, amr_controller,
interval = 5,
adapt_initial_condition = true,
adapt_initial_condition_only_refine = true)

stepsize_callback = StepsizeCallback(cfl = 1.2)

callbacks = CallbackSet(summary_callback,
analysis_callback,
alive_callback,
save_solution,
visualization,
amr_callback,
stepsize_callback)

###############################################################################
# run the simulation

sol = solve(ode, CarpenterKennedy2N54(williamson_condition = false),
dt = 1.0, # solve needs some value here but it will be overwritten by the stepsize_callback
save_everystep = false, callback = callbacks);
summary_callback() # print the timer summary
87 changes: 87 additions & 0 deletions examples/tree_3d_dgsem/elixir_euler_amr_visualization.jl
Original file line number Diff line number Diff line change
@@ -0,0 +1,87 @@

using OrdinaryDiffEq
using Trixi
using Plots
using GLMakie

###############################################################################
# semidiscretization of the compressible Euler equations

equations = CompressibleEulerEquations3D(1.4)

"""
initial_condition_density_pulse(x, t, equations::CompressibleEulerEquations3D)

A Gaussian pulse in the density with constant velocity and pressure; reduces the
compressible Euler equations to the linear advection equations.
"""
function initial_condition_density_pulse(x, t, equations::CompressibleEulerEquations3D)
rho = 1 + exp(-(x[1]^2 + x[2]^2 + x[3]^2)) / 2
v1 = 1
v2 = 1
v3 = 1
rho_v1 = rho * v1
rho_v2 = rho * v2
rho_v3 = rho * v3
p = 1
rho_e = p / (equations.gamma - 1) + 1 / 2 * rho * (v1^2 + v2^2 + v3^2)
return SVector(rho, rho_v1, rho_v2, rho_v3, rho_e)
end
initial_condition = initial_condition_density_pulse
solver = DGSEM(polydeg = 3, surface_flux = flux_lax_friedrichs)

coordinates_min = (-5.0, -5.0, -5.0)
coordinates_max = (5.0, 5.0, 5.0)
mesh = TreeMesh(coordinates_min, coordinates_max,
initial_refinement_level = 4,
n_cells_max = 10_000)

semi = SemidiscretizationHyperbolic(mesh, equations, initial_condition, solver)

###############################################################################
# ODE solvers, callbacks etc.

tspan = (0.0, 10.0)
ode = semidiscretize(semi, tspan)

summary_callback = SummaryCallback()

analysis_interval = 100
analysis_callback = AnalysisCallback(semi, interval = analysis_interval,
extra_analysis_integrals = (entropy,))

alive_callback = AliveCallback(analysis_interval = analysis_interval)

save_restart = SaveRestartCallback(interval = 100,
save_final_restart = true)

save_solution = SaveSolutionCallback(interval = 100,
save_initial_solution = true,
save_final_solution = true,
solution_variables = cons2prim)

visualization = VisualizationCallback(interval = 20, clims = (0, 1), plot_data_creator = Trixi.PlotData3D, plot_creator = Trixi.show_plot_makie)

amr_controller = ControllerThreeLevel(semi, IndicatorMax(semi, variable = first),
base_level = 1,
med_level = 2, med_threshold = 1.05,
max_level = 3, max_threshold = 1.3)
amr_callback = AMRCallback(semi, amr_controller,
interval = 5,
adapt_initial_condition = true,
adapt_initial_condition_only_refine = true)

stepsize_callback = StepsizeCallback(cfl = 0.9)

callbacks = CallbackSet(summary_callback,
analysis_callback, alive_callback,
save_restart, save_solution, visualization,
amr_callback, stepsize_callback);

###############################################################################
# run the simulation

sol = solve(ode, CarpenterKennedy2N54(williamson_condition = false),
dt = 1.0, # solve needs some value here but it will be overwritten by the stepsize_callback
save_everystep = false, callback = callbacks);
summary_callback() # print the timer summary
8 changes: 8 additions & 0 deletions src/Trixi.jl
Original file line number Diff line number Diff line change
Expand Up @@ -314,6 +314,14 @@ function __init__()
end
end

@require GLMakie="e9467ef8-e4e7-5192-8a1a-b1aee30e663a" begin
using .GLMakie: GLMakie
end

@require CairoMakie="13f3f980-e62b-5c42-98c6-ff1f3baf88f0" begin
using .CairoMakie: CairoMakie
end

@static if !isdefined(Base, :get_extension)
@require Convex="f65535da-76fb-5f13-bab9-19810c17039a" begin
@require ECOS="e2685f51-7e38-5353-a97d-a921fd2c8199" begin
Expand Down
Loading