Path: blob/main/examples/tree_2d_dgsem/elixir_advection_diffusion_nonperiodic.jl
2055 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))1617# Create a uniformly refined mesh with periodic boundaries18mesh = TreeMesh(coordinates_min, coordinates_max,19initial_refinement_level = 4,20periodicity = false,21n_cells_max = 30_000) # set maximum capacity of tree data structure2223# Example setup taken from24# - Truman Ellis, Jesse Chan, and Leszek Demkowicz (2016).25# Robust DPG methods for transient convection-diffusion.26# In: Building bridges: connections and challenges in modern approaches27# to numerical partial differential equations.28# [DOI](https://doi.org/10.1007/978-3-319-41640-3_6).29function initial_condition_eriksson_johnson(x, t, equations)30RealT = eltype(x)31l = 432epsilon = diffusivity() # TODO: this requires epsilon < .6 due to sqrt33lambda_1 = (-1 + sqrt(1 - 4 * epsilon * l)) / (-2 * epsilon)34lambda_2 = (-1 - sqrt(1 - 4 * epsilon * l)) / (-2 * epsilon)35r1 = (1 + sqrt(1 + 4 * convert(RealT, pi)^2 * epsilon^2)) / (2 * epsilon)36s1 = (1 - sqrt(1 + 4 * convert(RealT, pi)^2 * epsilon^2)) / (2 * epsilon)37u = exp(-l * t) * (exp(lambda_1 * x[1]) - exp(lambda_2 * x[1])) +38cospi(x[2]) * (exp(s1 * x[1]) - exp(r1 * x[1])) / (exp(-s1) - exp(-r1))39return SVector{1}(u)40end41initial_condition = initial_condition_eriksson_johnson4243boundary_conditions = (; x_neg = BoundaryConditionDirichlet(initial_condition),44y_neg = BoundaryConditionDirichlet(initial_condition),45y_pos = BoundaryConditionDirichlet(initial_condition),46x_pos = boundary_condition_do_nothing)4748boundary_conditions_parabolic = BoundaryConditionDirichlet(initial_condition)4950# A semidiscretization collects data structures and functions for the spatial discretization51semi = SemidiscretizationHyperbolicParabolic(mesh,52(equations, equations_parabolic),53initial_condition, solver;54solver_parabolic = ViscousFormulationBassiRebay1(),55boundary_conditions = (boundary_conditions,56boundary_conditions_parabolic))5758###############################################################################59# ODE solvers, callbacks etc.6061# Create ODE problem with time span `tspan`62tspan = (0.0, 1.5)63ode = semidiscretize(semi, tspan)6465# At the beginning of the main loop, the SummaryCallback prints a summary of the simulation setup66# and resets the timers67summary_callback = SummaryCallback()6869# The AnalysisCallback allows to analyse the solution in regular intervals and prints the results70analysis_interval = 10071analysis_callback = AnalysisCallback(semi, interval = analysis_interval)7273# The AliveCallback prints short status information in regular intervals74alive_callback = AliveCallback(analysis_interval = analysis_interval)7576# Create a CallbackSet to collect all callbacks such that they can be passed to the ODE solver77callbacks = CallbackSet(summary_callback, analysis_callback, alive_callback)7879###############################################################################80# run the simulation8182# OrdinaryDiffEq's `solve` method evolves the solution in time and executes the passed callbacks83time_int_tol = 1.0e-1184sol = solve(ode, RDPK3SpFSAL49(); abstol = time_int_tol, reltol = time_int_tol,85ode_default_options()..., callback = callbacks)868788