Path: blob/main/examples/tree_2d_dgsem/elixir_advection_amr_visualization.jl
5586 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, periodicity = true)2728semi = SemidiscretizationHyperbolic(mesh, equations, initial_condition, solver;29boundary_conditions = boundary_condition_periodic)3031###############################################################################32# ODE solvers, callbacks etc.3334tspan = (0.0, 20.0)35ode = semidiscretize(semi, tspan)3637summary_callback = SummaryCallback()3839analysis_interval = 10040analysis_callback = AnalysisCallback(semi, interval = analysis_interval,41extra_analysis_integrals = (entropy,))4243alive_callback = AliveCallback(analysis_interval = analysis_interval)4445save_solution = SaveSolutionCallback(interval = 100,46save_initial_solution = true,47save_final_solution = true,48solution_variables = cons2prim)4950# Enable in-situ visualization with a new plot generated every 20 time steps51# and additional plotting options passed as keyword arguments52visualization = VisualizationCallback(semi; interval = 20, clims = (0, 1))5354amr_controller = ControllerThreeLevel(semi, IndicatorMax(semi, variable = first),55base_level = 3,56med_level = 4, med_threshold = 0.1,57max_level = 5, max_threshold = 0.6)58amr_callback = AMRCallback(semi, amr_controller,59interval = 5,60adapt_initial_condition = true,61adapt_initial_condition_only_refine = true)6263stepsize_callback = StepsizeCallback(cfl = 1.6)6465callbacks = CallbackSet(summary_callback,66analysis_callback, alive_callback,67save_solution, visualization,68amr_callback, stepsize_callback);6970###############################################################################71# run the simulation7273sol = solve(ode, CarpenterKennedy2N54(williamson_condition = false);74dt = 1.0, # solve needs some value here but it will be overwritten by the stepsize_callback75ode_default_options()..., callback = callbacks);767778