Path: blob/main/examples/structured_2d_dgsem/elixir_advection_parallelogram.jl
2055 views
# This elixir transforms the setup of elixir_advection_basic to a parallelogram.1# The nodal values of the initial condition and the exact solution are the same as2# in elixir_advection_basic.3# However, on this non-rectangular mesh, the metric terms are non-trivial.4# The same errors as with elixir_advection_basic are expected.56using OrdinaryDiffEqLowStorageRK7using Trixi89# initial_condition_convergence_test transformed to the parallelogram10function initial_condition_parallelogram(x, t, equation::LinearScalarAdvectionEquation2D)11# Transform back to unit square12x_transformed = SVector(x[1] - x[2], x[2])13a = equation.advection_velocity14a_transformed = SVector(a[1] - a[2], a[2])1516# Store translated coordinate for easy use of exact solution17x_translated = x_transformed - a_transformed * t1819c = 1.020A = 0.521L = 222f = 1 / L23omega = 2 * pi * f24scalar = c + A * sin(omega * sum(x_translated))25return SVector(scalar)26end2728###############################################################################29# semidiscretization of the linear advection equation3031# Transformed advection_velocity = (0.2, -0.7) by transformation mapping32advection_velocity = (-0.5, -0.7)33equations = LinearScalarAdvectionEquation2D(advection_velocity)3435# Create DG solver with polynomial degree = 3 and (local) Lax-Friedrichs/Rusanov flux as surface flux36solver = DGSEM(polydeg = 3, surface_flux = flux_lax_friedrichs)3738# Define faces for a parallelogram that looks like this39#40# (0,1) __________ (2, 1)41# ⟋ ⟋42# ⟋ ⟋43# ⟋ ⟋44# (-2,-1) ‾‾‾‾‾‾‾‾‾‾ (0,-1)45mapping(xi, eta) = SVector(xi + eta, eta)4647cells_per_dimension = (16, 16)4849# Create curved mesh with 16 x 16 elements, periodic in both dimensions50mesh = StructuredMesh(cells_per_dimension, mapping; periodicity = (true, true))5152# A semidiscretization collects data structures and functions for the spatial discretization53semi = SemidiscretizationHyperbolic(mesh, equations, initial_condition_parallelogram,54solver)5556###############################################################################57# ODE solvers, callbacks etc.5859# Create ODE problem with time span from 0.0 to 1.060ode = semidiscretize(semi, (0.0, 1.0))6162# At the beginning of the main loop, the SummaryCallback prints a summary of the simulation setup63# and resets the timers64summary_callback = SummaryCallback()6566# The AnalysisCallback allows to analyse the solution in regular intervals and prints the results67analysis_callback = AnalysisCallback(semi, interval = 100)6869# The SaveSolutionCallback allows to save the solution to a file in regular intervals70save_solution = SaveSolutionCallback(interval = 100,71solution_variables = cons2prim)7273# The StepsizeCallback handles the re-calculation of the maximum Δt after each time step74stepsize_callback = StepsizeCallback(cfl = 1.6)7576# Create a CallbackSet to collect all callbacks such that they can be passed to the ODE solver77callbacks = CallbackSet(summary_callback, analysis_callback, save_solution,78stepsize_callback)7980###############################################################################81# run the simulation8283# OrdinaryDiffEq's `solve` method evolves the solution in time and executes the passed callbacks84sol = solve(ode, CarpenterKennedy2N54(williamson_condition = false);85dt = 1.0, # solve needs some value here but it will be overwritten by the stepsize_callback86ode_default_options()..., callback = callbacks);878889