Path: blob/main/examples/tree_1d_dgsem/elixir_advection_diffusion.jl
2055 views
using OrdinaryDiffEqSDIRK, ADTypes1using Trixi23###############################################################################4# semidiscretization of the linear advection diffusion equation56advection_velocity = 0.17equations = LinearScalarAdvectionEquation1D(advection_velocity)8diffusivity() = 0.19equations_parabolic = LaplaceDiffusion1D(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 = -convert(Float64, pi) # minimum coordinate15coordinates_max = convert(Float64, pi) # maximum coordinate1617# Create a uniformly refined mesh with periodic boundaries18mesh = TreeMesh(coordinates_min, coordinates_max,19initial_refinement_level = 4,20n_cells_max = 30_000, # set maximum capacity of tree data structure21periodicity = true)2223function x_trans_periodic(x, domain_length = SVector(oftype(x[1], 2 * pi)),24center = SVector(oftype(x[1], 0)))25x_normalized = x .- center26x_shifted = x_normalized .% domain_length27x_offset = ((x_shifted .< -0.5f0 * domain_length) -28(x_shifted .> 0.5f0 * domain_length)) .*29domain_length30return center + x_shifted + x_offset31end3233# Define initial condition34function initial_condition_diffusive_convergence_test(x, t,35equation::LinearScalarAdvectionEquation1D)36# Store translated coordinate for easy use of exact solution37x_trans = x_trans_periodic(x - equation.advection_velocity * t)3839nu = diffusivity()40c = 041A = 142omega = 143scalar = c + A * sin(omega * sum(x_trans)) * exp(-nu * omega^2 * t)44return SVector(scalar)45end46initial_condition = initial_condition_diffusive_convergence_test4748# define periodic boundary conditions everywhere49boundary_conditions = boundary_condition_periodic50boundary_conditions_parabolic = boundary_condition_periodic5152# A semidiscretization collects data structures and functions for the spatial discretization53semi = SemidiscretizationHyperbolicParabolic(mesh, (equations, equations_parabolic),54initial_condition,55solver;56boundary_conditions = (boundary_conditions,57boundary_conditions_parabolic))5859###############################################################################60# ODE solvers, callbacks etc.6162# Create ODE problem with time span from 0.0 to 1.063tspan = (0.0, 1.0)64ode = semidiscretize(semi, tspan)6566# At the beginning of the main loop, the SummaryCallback prints a summary of the simulation setup67# and resets the timers68summary_callback = SummaryCallback()6970# The AnalysisCallback allows to analyse the solution in regular intervals and prints the results71analysis_callback = AnalysisCallback(semi, interval = 100)7273# The AliveCallback prints short status information in regular intervals74alive_callback = AliveCallback(analysis_interval = 100)7576# The SaveRestartCallback allows to save a file from which a Trixi.jl simulation can be restarted77save_restart = SaveRestartCallback(interval = 100,78save_final_restart = true)7980# Create a CallbackSet to collect all callbacks such that they can be passed to the ODE solver81callbacks = CallbackSet(summary_callback, analysis_callback, alive_callback, save_restart)8283###############################################################################84# run the simulation8586# OrdinaryDiffEq's `solve` method evolves the solution in time and executes the passed callbacks87time_int_tol = 1.0e-1088time_abs_tol = 1.0e-1089sol = solve(ode, KenCarp4(autodiff = AutoFiniteDiff()); # This is an IMEX SDIRK method90abstol = time_abs_tol, reltol = time_int_tol,91ode_default_options()..., callback = callbacks)929394