Path: blob/main/examples/tree_1d_dgsem/elixir_diffusion_ldg.jl
2055 views
using OrdinaryDiffEqLowStorageRK1using Trixi23###############################################################################4# semidiscretization of the linear (advection) diffusion equation56advection_velocity = 0.0 # Note: This renders the equation mathematically purely parabolic7equations = LinearScalarAdvectionEquation1D(advection_velocity)8diffusivity() = 0.59equations_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)2223# Define initial condition if it is not defined already.24# For CI, the function is defined externally avoid "world age" issues that arise25# when running `Trixi.convergence_test`. The `isdefined` check is to allow the26# elixir to also be run outside of CI.27function initial_condition_pure_diffusion_1d_convergence_test(x, t,28equation)29nu = diffusivity()30c = 031A = 132omega = 133scalar = c + A * sin(omega * sum(x)) * exp(-nu * omega^2 * t)34return SVector(scalar)35end36initial_condition = initial_condition_pure_diffusion_1d_convergence_test3738# define periodic boundary conditions everywhere39boundary_conditions = boundary_condition_periodic40boundary_conditions_parabolic = boundary_condition_periodic4142# A semidiscretization collects data structures and functions for the spatial discretization43solver_parabolic = ViscousFormulationLocalDG()44semi = SemidiscretizationHyperbolicParabolic(mesh, (equations, equations_parabolic),45initial_condition,46solver; solver_parabolic,47boundary_conditions = (boundary_conditions,48boundary_conditions_parabolic))4950###############################################################################51# ODE solvers, callbacks etc.5253# Create ODE problem with time span from 0.0 to 0.154tspan = (0.0, 0.1)55ode = semidiscretize(semi, tspan)5657# At the beginning of the main loop, the SummaryCallback prints a summary of the simulation setup58# and resets the timers59summary_callback = SummaryCallback()6061# The AnalysisCallback allows to analyse the solution in regular intervals and prints the results62analysis_callback = AnalysisCallback(semi, interval = 100)6364# The AliveCallback prints short status information in regular intervals65alive_callback = AliveCallback(analysis_interval = 100)6667# The SaveRestartCallback allows to save a file from which a Trixi.jl simulation can be restarted68save_restart = SaveRestartCallback(interval = 100,69save_final_restart = true)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, save_restart)7374###############################################################################75# run the simulation7677# OrdinaryDiffEq's `solve` method evolves the solution in time and executes the passed callbacks78# For CI purposes, we use fixed time-stepping for this elixir.79sol = solve(ode, RDPK3SpFSAL35(); dt = 1.0e-3, adaptive = false,80ode_default_options()..., callback = callbacks)818283