Skip to content

Commit

Permalink
LinearizedEuler3D (#1903)
Browse files Browse the repository at this point in the history
* LinearizedEuler3D

* fmt

* Apply suggestions from code review

Co-authored-by: Lars Christmann <account-github@l12n.eu>

* nicer IC

* test vals

* fmt

* unit tests for full coverage

---------

Co-authored-by: Lars Christmann <account-github@l12n.eu>
  • Loading branch information
DanielDoehring and lchristm committed Apr 16, 2024
1 parent b1a84a6 commit 3cd1266
Show file tree
Hide file tree
Showing 11 changed files with 501 additions and 3 deletions.
64 changes: 64 additions & 0 deletions examples/p4est_3d_dgsem/elixir_linearizedeuler_convergence.jl
Original file line number Diff line number Diff line change
@@ -0,0 +1,64 @@
using OrdinaryDiffEq
using Trixi

###############################################################################
# semidiscretization of the linearized Euler equations

equations = LinearizedEulerEquations3D(v_mean_global = (0.0, 0.0, 0.0), c_mean_global = 1.0,
rho_mean_global = 1.0)

initial_condition = initial_condition_convergence_test

solver = DGSEM(polydeg = 3, surface_flux = flux_hll)

coordinates_min = (-1.0, -1.0, -1.0) # minimum coordinates (min(x), min(y), min(z))
coordinates_max = (1.0, 1.0, 1.0) # maximum coordinates (max(x), max(y), max(z))

# `initial_refinement_level` is provided here to allow for a
# convenient convergence test, see
# https://trixi-framework.github.io/Trixi.jl/stable/#Performing-a-convergence-analysis
trees_per_dimension = (4, 4, 4)
mesh = P4estMesh(trees_per_dimension, polydeg = 3,
coordinates_min = coordinates_min,
coordinates_max = coordinates_max,
initial_refinement_level = 0)

# A semidiscretization collects data structures and functions for the spatial discretization
semi = SemidiscretizationHyperbolic(mesh, equations, initial_condition, solver)

###############################################################################
# ODE solvers, callbacks etc.

# Create ODE problem with time span from 0.0 to 0.2
tspan = (0.0, 0.2)
ode = semidiscretize(semi, tspan)

# At the beginning of the main loop, the SummaryCallback prints a summary of the simulation setup
# and resets the timers
summary_callback = SummaryCallback()

analysis_interval = 100

# The AnalysisCallback allows to analyse the solution in regular intervals and prints the results
analysis_callback = AnalysisCallback(semi, interval = analysis_interval)

# The AliveCallback prints short status information in regular intervals
alive_callback = AliveCallback(analysis_interval = analysis_interval)

# The StepsizeCallback handles the re-calculation of the maximum Δt after each time step
stepsize_callback = StepsizeCallback(cfl = 0.8)

# Create a CallbackSet to collect all callbacks such that they can be passed to the ODE solver
callbacks = CallbackSet(summary_callback, analysis_callback, alive_callback,
stepsize_callback)

###############################################################################
# run the simulation

# OrdinaryDiffEq's `solve` method evolves the solution in time and executes the passed callbacks
sol = solve(ode, CarpenterKennedy2N54(williamson_condition = false),
dt = 1.0, # solve needs some value here but it will be overwritten by the stepsize_callback
save_everystep = false, callback = callbacks);

# print the timer summary
summary_callback() # print the timer summary
65 changes: 65 additions & 0 deletions examples/tree_3d_dgsem/elixir_linearizedeuler_gauss_wall.jl
Original file line number Diff line number Diff line change
@@ -0,0 +1,65 @@
using OrdinaryDiffEq
using Trixi

###############################################################################
# semidiscretization of the linearized Euler equations

equations = LinearizedEulerEquations3D(v_mean_global = (0.5, 0.5, 0.5), c_mean_global = 1.0,
rho_mean_global = 1.0)

solver = DGSEM(polydeg = 5, surface_flux = flux_lax_friedrichs)

coordinates_min = (0.0, 0.0, 0.0)
coordinates_max = (90.0, 90.0, 90.0)

mesh = TreeMesh(coordinates_min, coordinates_max,
initial_refinement_level = 4,
n_cells_max = 100_000,
periodicity = false)

# Initialize density and pressure perturbation with a Gaussian bump
# that splits into radial waves which are advected with v - c and v + c.
function initial_condition_gauss_wall(x, t, equations::LinearizedEulerEquations3D)
v1_prime = 0.0
v2_prime = 0.0
v3_prime = 0.0
rho_prime = p_prime = 2 * exp(-((x[1] - 45)^2 + (x[2] - 45)^2) / 25)
return SVector(rho_prime, v1_prime, v2_prime, v3_prime, p_prime)
end
initial_condition = initial_condition_gauss_wall

# A semidiscretization collects data structures and functions for the spatial discretization
semi = SemidiscretizationHyperbolic(mesh, equations, initial_condition, solver,
boundary_conditions = boundary_condition_wall)

###############################################################################
# ODE solvers, callbacks etc.

# At t = 30, the wave moving with v + c crashes into the wall
tspan = (0.0, 30.0)
ode = semidiscretize(semi, tspan)

# At the beginning of the main loop, the SummaryCallback prints a summary of the simulation setup
# and resets the timers
summary_callback = SummaryCallback()

# The AnalysisCallback allows to analyse the solution in regular intervals and prints the results
analysis_callback = AnalysisCallback(semi, interval = 100)

# The StepsizeCallback handles the re-calculation of the maximum Δt after each time step
stepsize_callback = StepsizeCallback(cfl = 0.9)

# Create a CallbackSet to collect all callbacks such that they can be passed to the ODE solver
callbacks = CallbackSet(summary_callback, analysis_callback,
stepsize_callback)

###############################################################################
# run the simulation

# OrdinaryDiffEq's `solve` method evolves the solution in time and executes the passed callbacks
sol = solve(ode, CarpenterKennedy2N54(williamson_condition = false),
dt = 1.0, # solve needs some value here but it will be overwritten by the stepsize_callback
save_everystep = false, callback = callbacks)

# Print the timer summary
summary_callback()
2 changes: 1 addition & 1 deletion src/Trixi.jl
Original file line number Diff line number Diff line change
Expand Up @@ -160,7 +160,7 @@ export AcousticPerturbationEquations2D,
LatticeBoltzmannEquations2D, LatticeBoltzmannEquations3D,
ShallowWaterEquations1D, ShallowWaterEquations2D,
ShallowWaterEquationsQuasi1D,
LinearizedEulerEquations1D, LinearizedEulerEquations2D,
LinearizedEulerEquations1D, LinearizedEulerEquations2D, LinearizedEulerEquations3D,
PolytropicEulerEquations2D,
TrafficFlowLWREquations1D

Expand Down
1 change: 1 addition & 0 deletions src/equations/equations.jl
Original file line number Diff line number Diff line change
Expand Up @@ -503,6 +503,7 @@ abstract type AbstractLinearizedEulerEquations{NDIMS, NVARS} <:
AbstractEquations{NDIMS, NVARS} end
include("linearized_euler_1d.jl")
include("linearized_euler_2d.jl")
include("linearized_euler_3d.jl")

abstract type AbstractEquationsParabolic{NDIMS, NVARS, GradientVariables} <:
AbstractEquations{NDIMS, NVARS} end
Expand Down
2 changes: 1 addition & 1 deletion src/equations/linearized_euler_1d.jl
Original file line number Diff line number Diff line change
Expand Up @@ -8,7 +8,7 @@
@doc raw"""
LinearizedEulerEquations1D(v_mean_global, c_mean_global, rho_mean_global)
Linearized euler equations in one space dimension. The equations are given by
Linearized Euler equations in one space dimension. The equations are given by
```math
\partial_t
\begin{pmatrix}
Expand Down
2 changes: 1 addition & 1 deletion src/equations/linearized_euler_2d.jl
Original file line number Diff line number Diff line change
Expand Up @@ -8,7 +8,7 @@
@doc raw"""
LinearizedEulerEquations2D(v_mean_global, c_mean_global, rho_mean_global)
Linearized euler equations in two space dimensions. The equations are given by
Linearized Euler equations in two space dimensions. The equations are given by
```math
\partial_t
\begin{pmatrix}
Expand Down
Loading

0 comments on commit 3cd1266

Please sign in to comment.