Path: blob/main/examples/tree_2d_dgsem/elixir_advection_diffusion.jl
2055 views
using OrdinaryDiffEqLowStorageRK1using Trixi23###############################################################################4# semidiscretization of the linear advection-diffusion equation56advection_velocity = (1.5, 1.0)7equations = LinearScalarAdvectionEquation2D(advection_velocity)8diffusivity() = 5.0e-29equations_parabolic = LaplaceDiffusion2D(diffusivity(), equations)1011# Create DG solver with polynomial degree = 3 and (local) Lax-Friedrichs/Rusanov flux as surface flux12solver = DGSEM(polydeg = 3, surface_flux = flux_lax_friedrichs)1314coordinates_min = (-1.0, -1.0) # minimum coordinates (min(x), min(y))15coordinates_max = (1.0, 1.0) # maximum coordinates (max(x), max(y))1617# Create a uniformly refined mesh with periodic boundaries18mesh = TreeMesh(coordinates_min, coordinates_max,19initial_refinement_level = 4,20periodicity = true,21n_cells_max = 30_000) # set maximum capacity of tree data structure2223# Define initial condition24function initial_condition_diffusive_convergence_test(x, t,25equation::LinearScalarAdvectionEquation2D)26# Store translated coordinate for easy use of exact solution27RealT = eltype(x)28x_trans = x - equation.advection_velocity * t2930nu = diffusivity()31c = 132A = 0.5f033L = 234f = 1.0f0 / L35omega = 2 * convert(RealT, pi) * f36scalar = c + A * sin(omega * sum(x_trans)) * exp(-2 * nu * omega^2 * t)37return SVector(scalar)38end39initial_condition = initial_condition_diffusive_convergence_test4041# define periodic boundary conditions everywhere42boundary_conditions = boundary_condition_periodic43boundary_conditions_parabolic = boundary_condition_periodic4445# A semidiscretization collects data structures and functions for the spatial discretization46semi = SemidiscretizationHyperbolicParabolic(mesh,47(equations, equations_parabolic),48initial_condition, solver;49solver_parabolic = ViscousFormulationBassiRebay1(),50boundary_conditions = (boundary_conditions,51boundary_conditions_parabolic))5253###############################################################################54# ODE solvers, callbacks etc.5556# Create ODE problem with time span from 0.0 to 1.557tspan = (0.0, 1.5)58ode = semidiscretize(semi, tspan)5960# At the beginning of the main loop, the SummaryCallback prints a summary of the simulation setup61# and resets the timers62summary_callback = SummaryCallback()6364# The AnalysisCallback allows to analyse the solution in regular intervals and prints the results65analysis_interval = 10066analysis_callback = AnalysisCallback(semi, interval = analysis_interval)6768# The AliveCallback prints short status information in regular intervals69alive_callback = AliveCallback(analysis_interval = analysis_interval)7071# Create a CallbackSet to collect all callbacks such that they can be passed to the ODE solver72callbacks = CallbackSet(summary_callback, analysis_callback, alive_callback)7374###############################################################################75# run the simulation7677# OrdinaryDiffEq's `solve` method evolves the solution in time and executes the passed callbacks78alg = RDPK3SpFSAL49()79time_int_tol = 1.0e-1180sol = solve(ode, alg; abstol = time_int_tol, reltol = time_int_tol,81ode_default_options()..., callback = callbacks)828384