Path: blob/main/examples/tree_3d_dgsem/elixir_mhdmultiion_ec.jl
5586 views
1using OrdinaryDiffEqLowStorageRK2using Trixi34###############################################################################5# semidiscretization of the ideal MHD equations6equations = IdealGlmMhdMultiIonEquations3D(gammas = (1.4, 1.667),7charge_to_mass = (1.0, 2.0))89initial_condition = initial_condition_weak_blast_wave1011volume_flux = (flux_ruedaramirez_etal, flux_nonconservative_ruedaramirez_etal)12surface_flux = (flux_ruedaramirez_etal, flux_nonconservative_ruedaramirez_etal)13solver = DGSEM(polydeg = 3, surface_flux = surface_flux,14volume_integral = VolumeIntegralFluxDifferencing(volume_flux))1516coordinates_min = (-2.0, -2.0, -2.0)17coordinates_max = (2.0, 2.0, 2.0)18mesh = TreeMesh(coordinates_min, coordinates_max,19initial_refinement_level = 4,20n_cells_max = 10_000, periodicity = true)2122# The multi-ion GLM-MHD equations require the inclusion of source_terms_lorentz23# whenever multiple ion species are present24semi = SemidiscretizationHyperbolic(mesh, equations, initial_condition, solver;25source_terms = source_terms_lorentz,26boundary_conditions = boundary_condition_periodic)2728###############################################################################29# ODE solvers, callbacks etc.3031tspan = (0.0, 0.4)32ode = semidiscretize(semi, tspan)3334summary_callback = SummaryCallback()3536analysis_interval = 1037analysis_callback = AnalysisCallback(semi, interval = analysis_interval)38alive_callback = AliveCallback(analysis_interval = analysis_interval)3940save_solution = SaveSolutionCallback(dt = 0.1,41save_initial_solution = true,42save_final_solution = true,43solution_variables = cons2prim)4445cfl = 0.54647stepsize_callback = StepsizeCallback(cfl = cfl)4849glm_speed_callback = GlmSpeedCallback(glm_scale = 0.5, cfl = cfl)5051callbacks = CallbackSet(summary_callback,52analysis_callback, alive_callback,53save_solution,54stepsize_callback,55glm_speed_callback)5657###############################################################################58# run the simulation5960sol = solve(ode, CarpenterKennedy2N54(williamson_condition = false);61dt = 1.0, # solve needs some value here but it will be overwritten by the stepsize_callback62ode_default_options()..., callback = callbacks);636465