Path: blob/main/examples/structured_2d_dgsem/elixir_advection_meshview.jl
2055 views
using OrdinaryDiffEqLowStorageRK1using Trixi23###############################################################################4# Coupled semidiscretization of two linear advection systems using converter functions5# and mesh views for the semidiscretizations. First we define a parent mesh6# for the entire physical domain, then we define the two mesh views on this parent.7#8# In this elixir, we have a square domain that is divided into left and right subdomains.9# On each half of the domain, a completely independent `SemidiscretizationHyperbolic`10# is created for the linear advection equations. The two systems are coupled in the11# x-direction.12# For a high-level overview, see also the figure below:13#14# (-1, 1) ( 1, 1)15# ┌────────────────────┬────────────────────┐16# │ ↑ periodic ↑ │ ↑ periodic ↑ │17# │ │ │18# │ ========= │ ========= │19# │ system #1 │ system #2 │20# │ ========= │ ========= │21# │ │ │22# │<-- coupled │<-- coupled │23# │ coupled -->│ coupled -->│24# │ │ │25# │ ↓ periodic ↓ │ ↓ periodic ↓ │26# └────────────────────┴────────────────────┘27# (-1, -1) ( 1, -1)2829advection_velocity = (0.2, -0.7)30equations = LinearScalarAdvectionEquation2D(advection_velocity)3132# Create DG solver with polynomial degree = 333solver = DGSEM(polydeg = 3, surface_flux = flux_lax_friedrichs)3435# Domain size of the parent mesh.36coordinates_min = (-1.0, -1.0)37coordinates_max = (1.0, 1.0)3839# Cell dimensions of the parent mesh.40cells_per_dimension = (16, 16)4142# Create parent mesh with 16 x 16 elements43parent_mesh = StructuredMesh(cells_per_dimension, coordinates_min, coordinates_max)4445# Create the two mesh views, each of which takes half of the parent mesh.46mesh1 = StructuredMeshView(parent_mesh; indices_min = (1, 1), indices_max = (8, 16))47mesh2 = StructuredMeshView(parent_mesh; indices_min = (9, 1), indices_max = (16, 16))4849# The coupling function is simply the identity, as we are dealing with two identical systems.50coupling_function = (x, u, equations_other, equations_own) -> u5152# Define the coupled boundary conditions53# The indices (:end, :i_forward) and (:begin, :i_forward) denote the interface indexing.54# For a system with coupling in x and y see examples/structured_2d_dgsem/elixir_advection_coupled.jl.55boundary_conditions1 = (56# Connect left boundary with right boundary of left mesh57x_neg = BoundaryConditionCoupled(2, (:end, :i_forward), Float64,58coupling_function),59x_pos = BoundaryConditionCoupled(2, (:begin, :i_forward), Float64,60coupling_function),61y_neg = boundary_condition_periodic,62y_pos = boundary_condition_periodic)63boundary_conditions2 = (64# Connect left boundary with right boundary of left mesh65x_neg = BoundaryConditionCoupled(1, (:end, :i_forward), Float64,66coupling_function),67x_pos = BoundaryConditionCoupled(1, (:begin, :i_forward), Float64,68coupling_function),69y_neg = boundary_condition_periodic,70y_pos = boundary_condition_periodic)7172# A semidiscretization collects data structures and functions for the spatial discretization73semi1 = SemidiscretizationHyperbolic(mesh1, equations, initial_condition_convergence_test,74solver,75boundary_conditions = boundary_conditions1)76semi2 = SemidiscretizationHyperbolic(mesh2, equations, initial_condition_convergence_test,77solver,78boundary_conditions = boundary_conditions2)79semi = SemidiscretizationCoupled(semi1, semi2)8081###############################################################################82# ODE solvers, callbacks etc.8384# Create ODE problem with time span from 0.0 to 1.085ode = semidiscretize(semi, (0.0, 1.0))8687# At the beginning of the main loop, the SummaryCallback prints a summary of the simulation setup88# and resets the timers89summary_callback = SummaryCallback()9091# The AnalysisCallback allows to analyse the solution in regular intervals and prints the results92analysis_callback1 = AnalysisCallback(semi1, interval = 100)93analysis_callback2 = AnalysisCallback(semi2, interval = 100)94analysis_callback = AnalysisCallbackCoupled(semi, analysis_callback1, analysis_callback2)9596analysis_interval = 10097alive_callback = AliveCallback(analysis_interval = analysis_interval)9899# The SaveSolutionCallback allows to save the solution to a file in regular intervals100save_solution = SaveSolutionCallback(interval = 100,101save_initial_solution = true,102save_final_solution = true,103solution_variables = cons2prim)104105# The StepsizeCallback handles the re-calculation of the maximum Δt after each time step106stepsize_callback = StepsizeCallback(cfl = 1.6)107108# Create a CallbackSet to collect all callbacks such that they can be passed to the ODE solver109callbacks = CallbackSet(summary_callback, analysis_callback, save_solution,110stepsize_callback)111112###############################################################################113# run the simulation114115# OrdinaryDiffEq's `solve` method evolves the solution in time and executes the passed callbacks116sol = solve(ode, CarpenterKennedy2N54(williamson_condition = false);117dt = 5.0e-2, # solve needs some value here but it will be overwritten by the stepsize_callback118ode_default_options()..., callback = callbacks);119120121