Path: blob/main/examples/tree_2d_dgsem/elixir_advection_amr_visualization.jl
2055 views
using OrdinaryDiffEqLowStorageRK1using Trixi2using Plots34###############################################################################5# semidiscretization of the linear advection equation67advection_velocity = (0.2, -0.7)8equations = LinearScalarAdvectionEquation2D(advection_velocity)910function initial_condition_gauss_largedomain(x, t,11equation::LinearScalarAdvectionEquation2D)12# Store translated coordinate for easy use of exact solution13domain_length = SVector(10, 10)14x_trans = Trixi.x_trans_periodic_2d(x - equation.advection_velocity * t, domain_length)1516return SVector(exp(-(x_trans[1]^2 + x_trans[2]^2)))17end18initial_condition = initial_condition_gauss_largedomain1920solver = DGSEM(polydeg = 3, surface_flux = flux_lax_friedrichs)2122coordinates_min = (-5.0, -5.0)23coordinates_max = (5.0, 5.0)24mesh = TreeMesh(coordinates_min, coordinates_max,25initial_refinement_level = 3,26n_cells_max = 30_000)2728semi = SemidiscretizationHyperbolic(mesh, equations, initial_condition, solver)2930###############################################################################31# ODE solvers, callbacks etc.3233tspan = (0.0, 20.0)34ode = semidiscretize(semi, tspan)3536summary_callback = SummaryCallback()3738analysis_interval = 10039analysis_callback = AnalysisCallback(semi, interval = analysis_interval,40extra_analysis_integrals = (entropy,))4142alive_callback = AliveCallback(analysis_interval = analysis_interval)4344save_solution = SaveSolutionCallback(interval = 100,45save_initial_solution = true,46save_final_solution = true,47solution_variables = cons2prim)4849# Enable in-situ visualization with a new plot generated every 20 time steps50# and additional plotting options passed as keyword arguments51visualization = VisualizationCallback(semi; interval = 20, clims = (0, 1))5253amr_controller = ControllerThreeLevel(semi, IndicatorMax(semi, variable = first),54base_level = 3,55med_level = 4, med_threshold = 0.1,56max_level = 5, max_threshold = 0.6)57amr_callback = AMRCallback(semi, amr_controller,58interval = 5,59adapt_initial_condition = true,60adapt_initial_condition_only_refine = true)6162stepsize_callback = StepsizeCallback(cfl = 1.6)6364callbacks = CallbackSet(summary_callback,65analysis_callback, alive_callback,66save_solution, visualization,67amr_callback, stepsize_callback);6869###############################################################################70# run the simulation7172sol = solve(ode, CarpenterKennedy2N54(williamson_condition = false);73dt = 1.0, # solve needs some value here but it will be overwritten by the stepsize_callback74ode_default_options()..., callback = callbacks);757677