Path: blob/main/examples/structured_2d_dgsem/elixir_mhd_coupled.jl
2055 views
using OrdinaryDiffEqLowStorageRK1using Trixi23###############################################################################4# Two semidiscretizations of the ideal GLM-MHD systems using converter functions such that5# they are coupled across the domain boundaries to generate a periodic system.6#7# In this elixir, we have a square domain that is divided into a left and right half.8# On each half of the domain, an independent SemidiscretizationHyperbolic is created for9# each set of ideal GLM-MHD equations. The two systems are coupled in the x-direction10# and are periodic in the y-direction.11# For a high-level overview, see also the figure below:12#13# (-2, 2) ( 2, 2)14# ┌────────────────────┬────────────────────┐15# │ ↑ periodic ↑ │ ↑ periodic ↑ │16# │ │ │17# │ ========= │ ========= │18# │ system #1 │ system #2 │19# │ ========= │ ========= │20# │ │ │21# │<-- coupled │<-- coupled │22# │ coupled -->│ coupled -->│23# │ │ │24# │ ↓ periodic ↓ │ ↓ periodic ↓ │25# └────────────────────┴────────────────────┘26# (-2, -2) ( 2, -2)2728gamma = 5 / 329equations = IdealGlmMhdEquations2D(gamma)3031cells_per_dimension = (32, 64)3233# Up to version 0.13.0, `max_abs_speed_naive` was used as the default wave speed estimate of34# `const flux_lax_friedrichs = FluxLaxFriedrichs(), i.e., `FluxLaxFriedrichs(max_abs_speed = max_abs_speed_naive)`.35# In the `StepsizeCallback`, though, the less diffusive `max_abs_speeds` is employed which is consistent with `max_abs_speed`.36# Thus, we exchanged in PR#2458 the default wave speed used in the LLF flux to `max_abs_speed`.37# To ensure that every example still runs we specify explicitly `FluxLaxFriedrichs(max_abs_speed_naive)`.38# We remark, however, that the now default `max_abs_speed` is in general recommended due to compliance with the39# `StepsizeCallback` (CFL-Condition) and less diffusion.40surface_flux = (FluxLaxFriedrichs(max_abs_speed_naive), flux_nonconservative_powell)41volume_flux = (flux_hindenlang_gassner, flux_nonconservative_powell)42solver = DGSEM(polydeg = 3, surface_flux = surface_flux,43volume_integral = VolumeIntegralFluxDifferencing(volume_flux))4445###########46# system #147###########4849initial_condition1 = initial_condition_convergence_test50coordinates_min1 = (-1 / sin(pi / 4), -1 / sin(pi / 4))51coordinates_max1 = (0.0, 1 / sin(pi / 4))52mesh1 = StructuredMesh(cells_per_dimension,53coordinates_min1,54coordinates_max1,55periodicity = (false, true))5657coupling_function1 = (x, u, equations_other, equations_own) -> u58boundary_conditions1 = (x_neg = BoundaryConditionCoupled(2, (:end, :i_forward), Float64,59coupling_function1),60x_pos = BoundaryConditionCoupled(2, (:begin, :i_forward), Float64,61coupling_function1),62y_neg = boundary_condition_periodic,63y_pos = boundary_condition_periodic)6465semi1 = SemidiscretizationHyperbolic(mesh1, equations, initial_condition1, solver,66boundary_conditions = boundary_conditions1)6768###########69# system #270###########7172initial_condition2 = initial_condition_convergence_test73coordinates_min2 = (0.0, -1 / sin(pi / 4))74coordinates_max2 = (1 / sin(pi / 4), 1 / sin(pi / 4))75mesh2 = StructuredMesh(cells_per_dimension,76coordinates_min2,77coordinates_max2,78periodicity = (false, true))7980coupling_function2 = (x, u, equations_other, equations_own) -> u81boundary_conditions2 = (x_neg = BoundaryConditionCoupled(1, (:end, :i_forward), Float64,82coupling_function2),83x_pos = BoundaryConditionCoupled(1, (:begin, :i_forward), Float64,84coupling_function2),85y_neg = boundary_condition_periodic,86y_pos = boundary_condition_periodic)8788semi2 = SemidiscretizationHyperbolic(mesh2, equations, initial_condition2, solver,89boundary_conditions = boundary_conditions2)9091# Create a semidiscretization that bundles all the semidiscretizations.92semi = SemidiscretizationCoupled(semi1, semi2)9394###############################################################################95# ODE solvers, callbacks etc.9697tspan = (0.0, 0.1)98ode = semidiscretize(semi, tspan)99100summary_callback = SummaryCallback()101102analysis_interval = 100103104analysis_callback1 = AnalysisCallback(semi1, interval = 100)105analysis_callback2 = AnalysisCallback(semi2, interval = 100)106analysis_callback = AnalysisCallbackCoupled(semi, analysis_callback1, analysis_callback2)107108alive_callback = AliveCallback(analysis_interval = analysis_interval)109110save_solution = SaveSolutionCallback(interval = 50,111save_initial_solution = true,112save_final_solution = true,113solution_variables = cons2prim)114115cfl = 1.0116117stepsize_callback = StepsizeCallback(cfl = cfl)118119glm_speed_callback = GlmSpeedCallback(glm_scale = 0.5, cfl = cfl,120semi_indices = [1, 2])121122callbacks = CallbackSet(summary_callback,123analysis_callback, alive_callback,124save_solution,125stepsize_callback,126glm_speed_callback)127128###############################################################################129# run the simulation130131sol = solve(ode, CarpenterKennedy2N54(williamson_condition = false);132dt = 0.01, # solve needs some value here but it will be overwritten by the stepsize_callback133ode_default_options()..., callback = callbacks);134135136