Path: blob/main/examples/tree_2d_dgsem/elixir_advection_diffusion_nonperiodic.jl
5586 views
using OrdinaryDiffEqLowStorageRK1using Trixi23###############################################################################4# semidiscretization of the linear advection-diffusion equation56diffusivity() = 5.0e-27advection_velocity = (1.0, 0.0)8equations = LinearScalarAdvectionEquation2D(advection_velocity)9equations_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, -0.5) # minimum coordinates (min(x), min(y))15coordinates_max = (0.0, 0.5) # maximum coordinates (max(x), max(y))1617mesh = TreeMesh(coordinates_min, coordinates_max,18initial_refinement_level = 4,19periodicity = false,20n_cells_max = 30_000) # set maximum capacity of tree data structure2122# Example setup taken from23# - Truman Ellis, Jesse Chan, and Leszek Demkowicz (2016).24# Robust DPG methods for transient convection-diffusion.25# In: Building bridges: connections and challenges in modern approaches26# to numerical partial differential equations.27# [DOI](https://doi.org/10.1007/978-3-319-41640-3_6).28function initial_condition_eriksson_johnson(x, t, equations)29RealT = eltype(x)30l = 431epsilon = diffusivity() # NOTE: this requires epsilon <= 1/16 due to sqrt32lambda_1 = (-1 + sqrt(1 - 4 * epsilon * l)) / (-2 * epsilon)33lambda_2 = (-1 - sqrt(1 - 4 * epsilon * l)) / (-2 * epsilon)34r1 = (1 + sqrt(1 + 4 * convert(RealT, pi)^2 * epsilon^2)) / (2 * epsilon)35s1 = (1 - sqrt(1 + 4 * convert(RealT, pi)^2 * epsilon^2)) / (2 * epsilon)36u = exp(-l * t) * (exp(lambda_1 * x[1]) - exp(lambda_2 * x[1])) +37cospi(x[2]) * (exp(s1 * x[1]) - exp(r1 * x[1])) / (exp(-s1) - exp(-r1))38return SVector{1}(u)39end40initial_condition = initial_condition_eriksson_johnson4142boundary_conditions = (; x_neg = BoundaryConditionDirichlet(initial_condition),43y_neg = BoundaryConditionDirichlet(initial_condition),44y_pos = BoundaryConditionDirichlet(initial_condition),45x_pos = boundary_condition_do_nothing)4647boundary_conditions_parabolic = BoundaryConditionDirichlet(initial_condition)4849# A semidiscretization collects data structures and functions for the spatial discretization50semi = SemidiscretizationHyperbolicParabolic(mesh,51(equations, equations_parabolic),52initial_condition, solver;53solver_parabolic = ParabolicFormulationBassiRebay1(),54boundary_conditions = (boundary_conditions,55boundary_conditions_parabolic))5657###############################################################################58# ODE solvers, callbacks etc.5960# Create ODE problem with time span `tspan`61tspan = (0.0, 1.5)62ode = semidiscretize(semi, tspan)6364# At the beginning of the main loop, the SummaryCallback prints a summary of the simulation setup65# and resets the timers66summary_callback = SummaryCallback()6768# The AnalysisCallback allows to analyse the solution in regular intervals and prints the results69analysis_interval = 10070analysis_callback = AnalysisCallback(semi, interval = analysis_interval)7172# The AliveCallback prints short status information in regular intervals73alive_callback = AliveCallback(analysis_interval = analysis_interval)7475# Create a CallbackSet to collect all callbacks such that they can be passed to the ODE solver76callbacks = CallbackSet(summary_callback, analysis_callback, alive_callback)7778###############################################################################79# run the simulation8081# OrdinaryDiffEq's `solve` method evolves the solution in time and executes the passed callbacks82time_int_tol = 1.0e-1183sol = solve(ode, RDPK3SpFSAL49(); abstol = time_int_tol, reltol = time_int_tol,84ode_default_options()..., callback = callbacks)858687