Path: blob/main/examples/structured_2d_dgsem/elixir_mhdmultiion_ec.jl
5586 views
using OrdinaryDiffEqLowStorageRK1using Trixi23###############################################################################4# semidiscretization of the ideal multi-ion MHD equations5equations = IdealGlmMhdMultiIonEquations2D(gammas = (1.4, 1.667),6charge_to_mass = (1.0, 2.0))78initial_condition = initial_condition_weak_blast_wave910# Entropy conservative numerical fluxes11volume_flux = (flux_ruedaramirez_etal, flux_nonconservative_ruedaramirez_etal)12surface_flux = (flux_ruedaramirez_etal, flux_nonconservative_ruedaramirez_etal)1314solver = DGSEM(polydeg = 3, surface_flux = surface_flux,15volume_integral = VolumeIntegralFluxDifferencing(volume_flux))1617coordinates_min = (-2.0, -2.0)18coordinates_max = (2.0, 2.0)19cells_per_dimension = (100, 100)20mesh = StructuredMesh(cells_per_dimension, coordinates_min,21coordinates_max, periodicity = true)2223# The multi-ion GLM-MHD equations require the inclusion of source_terms_lorentz24# whenever multiple ion species are present25semi = SemidiscretizationHyperbolic(mesh, equations, initial_condition, solver,26source_terms = source_terms_lorentz,27boundary_conditions = boundary_condition_periodic)2829###############################################################################30# ODE solvers, callbacks etc.3132tspan = (0.0, 0.4)33ode = semidiscretize(semi, tspan)3435summary_callback = SummaryCallback()3637analysis_interval = 1038analysis_callback = AnalysisCallback(semi, interval = analysis_interval)39alive_callback = AliveCallback(analysis_interval = analysis_interval)4041save_solution = SaveSolutionCallback(dt = 0.1, # interval=100,42save_initial_solution = true,43save_final_solution = true,44solution_variables = cons2prim)4546cfl = 0.54748stepsize_callback = StepsizeCallback(cfl = cfl)4950glm_speed_callback = GlmSpeedCallback(glm_scale = 0.5, cfl = cfl)5152callbacks = CallbackSet(summary_callback,53analysis_callback, alive_callback,54save_solution,55stepsize_callback,56glm_speed_callback)5758###############################################################################59# run the simulation6061sol = solve(ode, CarpenterKennedy2N54(williamson_condition = false);62dt = 1.0, # solve needs some value here but it will be overwritten by the stepsize_callback63ode_default_options()..., callback = callbacks);646566