Path: blob/main/examples/t8code_3d_dgsem/elixir_advection_unstructured_curved.jl
2055 views
using OrdinaryDiffEqLowStorageRK1using Trixi23###############################################################################4# semidiscretization of the linear advection equation56advection_velocity = (0.2, -0.7, 0.5)7equations = LinearScalarAdvectionEquation3D(advection_velocity)89# Create DG solver with polynomial degree = 3 and (local) Lax-Friedrichs/Rusanov flux as surface flux10solver = DGSEM(polydeg = 3, surface_flux = flux_lax_friedrichs)1112initial_condition = initial_condition_convergence_test1314boundary_condition = BoundaryConditionDirichlet(initial_condition)15boundary_conditions = Dict(:all => boundary_condition)1617# Mapping as described in https://arxiv.org/abs/2012.12040 but with less warping.18# The mapping will be interpolated at tree level, and then refined without changing19# the geometry interpolant. The original mapping applied to this unstructured mesh20# causes some Jacobians to be negative, which makes the mesh invalid.21function mapping(xi, eta, zeta)22# Don't transform input variables between -1 and 1 onto [0,3] to obtain curved boundaries23# xi = 1.5 * xi_ + 1.524# eta = 1.5 * eta_ + 1.525# zeta = 1.5 * zeta_ + 1.52627y = eta +281 / 6 * (cos(1.5 * pi * (2 * xi - 3) / 3) *29cos(0.5 * pi * (2 * eta - 3) / 3) *30cos(0.5 * pi * (2 * zeta - 3) / 3))3132x = xi +331 / 6 * (cos(0.5 * pi * (2 * xi - 3) / 3) *34cos(2 * pi * (2 * y - 3) / 3) *35cos(0.5 * pi * (2 * zeta - 3) / 3))3637z = zeta +381 / 6 * (cos(0.5 * pi * (2 * x - 3) / 3) *39cos(pi * (2 * y - 3) / 3) *40cos(0.5 * pi * (2 * zeta - 3) / 3))4142return SVector(x, y, z)43end4445# Unstructured mesh with 68 cells of the cube domain [-1, 1]^346mesh_file = Trixi.download("https://gist.githubusercontent.com/efaulhaber/d45c8ac1e248618885fa7cc31a50ab40/raw/37fba24890ab37cfa49c39eae98b44faf4502882/cube_unstructured_1.inp",47joinpath(@__DIR__, "cube_unstructured_1.inp"))4849mesh = T8codeMesh(mesh_file, 3; polydeg = 3,50mapping = mapping,51initial_refinement_level = 2)5253# A semidiscretization collects data structures and functions for the spatial discretization54semi = SemidiscretizationHyperbolic(mesh, equations, initial_condition, solver,55boundary_conditions = boundary_conditions)5657###############################################################################58# ODE solvers, callbacks etc.5960# Create ODE problem with time span from 0.0 to 0.161ode = semidiscretize(semi, (0.0, 0.1))6263# At the beginning of the main loop, the SummaryCallback prints a summary of the simulation setup64# and resets the timers65summary_callback = SummaryCallback()6667# The AnalysisCallback allows to analyse the solution in regular intervals and prints the results68analysis_interval = 10069analysis_callback = AnalysisCallback(semi, interval = analysis_interval)7071alive_callback = AliveCallback(analysis_interval = analysis_interval)7273# The SaveRestartCallback allows to save a file from which a Trixi.jl simulation can be restarted74save_restart = SaveRestartCallback(interval = 100,75save_final_restart = true)7677# The SaveSolutionCallback allows to save the solution to a file in regular intervals78save_solution = SaveSolutionCallback(interval = 100,79solution_variables = cons2prim)8081# The StepsizeCallback handles the re-calculation of the maximum Δt after each time step82stepsize_callback = StepsizeCallback(cfl = 1.2)8384# Create a CallbackSet to collect all callbacks such that they can be passed to the ODE solver85callbacks = CallbackSet(summary_callback, analysis_callback, alive_callback, save_restart,86save_solution, stepsize_callback)8788###############################################################################89# run the simulation9091# OrdinaryDiffEq's `solve` method evolves the solution in time and executes the passed callbacks92sol = solve(ode, CarpenterKennedy2N54(williamson_condition = false);93dt = 1.0, # solve needs some value here but it will be overwritten by the stepsize_callback94ode_default_options()..., callback = callbacks);959697