Path: blob/main/examples/p4est_3d_dgsem/elixir_mhdmultiion_ec.jl
5586 views
1using OrdinaryDiffEqLowStorageRK2using Trixi34###############################################################################5# Compute non-trivial electron pressure as 0.2 * total_pressure6function electron_pressure(u, equations::IdealGlmMhdMultiIonEquations3D)7prim = cons2prim(u, equations)8pe = zero(eltype(u))9for k in eachcomponent(equations)10rho, v1, v2, v3, p = get_component(k, prim, equations)11pe += 0.2 * p12end13return zero(u[1])14end1516# semidiscretization of the ideal GLM-MHD multi-ion equations17equations = IdealGlmMhdMultiIonEquations3D(gammas = (1.4, 1.667),18charge_to_mass = (1.0, 2.0))1920initial_condition = initial_condition_weak_blast_wave2122volume_flux = (flux_ruedaramirez_etal, flux_nonconservative_ruedaramirez_etal)23surface_flux = (flux_ruedaramirez_etal, flux_nonconservative_ruedaramirez_etal)24# For provably entropy-stable surface fluxes, use25# surface_flux = (FluxPlusDissipation(flux_ruedaramirez_etal, DissipationLaxFriedrichsEntropyVariables()),26# flux_nonconservative_ruedaramirez_etal)27solver = DGSEM(polydeg = 3, surface_flux = surface_flux,28volume_integral = VolumeIntegralFluxDifferencing(volume_flux))2930# Mapping as described in https://arxiv.org/abs/2012.1204031function mapping(xi_, eta_, zeta_)32# Transform input variables between -1 and 1 onto [0,3]33xi = 1.5 * xi_ + 1.534eta = 1.5 * eta_ + 1.535zeta = 1.5 * zeta_ + 1.53637y = eta +383 / 8 * (cos(1.5 * pi * (2 * xi - 3) / 3) *39cos(0.5 * pi * (2 * eta - 3) / 3) *40cos(0.5 * pi * (2 * zeta - 3) / 3))4142x = xi +433 / 8 * (cos(0.5 * pi * (2 * xi - 3) / 3) *44cos(2 * pi * (2 * y - 3) / 3) *45cos(0.5 * pi * (2 * zeta - 3) / 3))4647z = zeta +483 / 8 * (cos(0.5 * pi * (2 * x - 3) / 3) *49cos(pi * (2 * y - 3) / 3) *50cos(0.5 * pi * (2 * zeta - 3) / 3))5152return SVector(x, y, z)53end5455trees_per_dimension = (8, 8, 8)56mesh = P4estMesh(trees_per_dimension,57polydeg = 3, mapping = mapping,58periodicity = true)5960# The multi-ion GLM-MHD equations require the inclusion of source_terms_lorentz61# whenever multiple ion species are present62semi = SemidiscretizationHyperbolic(mesh, equations, initial_condition, solver,63source_terms = source_terms_lorentz,64boundary_conditions = boundary_condition_periodic)6566###############################################################################67# ODE solvers, callbacks etc.6869tspan = (0.0, 0.4)70ode = semidiscretize(semi, tspan)7172summary_callback = SummaryCallback()7374analysis_interval = 1075analysis_callback = AnalysisCallback(semi, interval = analysis_interval)76alive_callback = AliveCallback(analysis_interval = analysis_interval)7778save_solution = SaveSolutionCallback(dt = 0.1,79save_initial_solution = true,80save_final_solution = true,81solution_variables = cons2prim)8283cfl = 0.58485stepsize_callback = StepsizeCallback(cfl = cfl)8687glm_speed_callback = GlmSpeedCallback(glm_scale = 0.5, cfl = cfl)8889callbacks = CallbackSet(summary_callback,90analysis_callback, alive_callback,91save_solution,92stepsize_callback,93glm_speed_callback)9495###############################################################################96# run the simulation9798sol = solve(ode, CarpenterKennedy2N54(williamson_condition = false);99dt = 1.0, # solve needs some value here but it will be overwritten by the stepsize_callback100ode_default_options()..., callback = callbacks);101102103