Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
trixi-framework
GitHub Repository: trixi-framework/Trixi.jl
Path: blob/main/examples/p4est_2d_dgsem/elixir_advection_coupled.jl
5586 views
1
using OrdinaryDiffEqSSPRK, OrdinaryDiffEqLowStorageRK
2
using Trixi
3
4
###############################################################################
5
# Simplest coupled setup consisting of two non-trivial mesh views.
6
7
advection_velocity = (0.2, -0.7)
8
equations = LinearScalarAdvectionEquation2D(advection_velocity)
9
10
# Create DG solver with polynomial degree = 3 and (local) Lax-Friedrichs/Rusanov flux as surface flux
11
solver = DGSEM(polydeg = 3, surface_flux = flux_lax_friedrichs)
12
13
# Define the physical domain for the parent mesh.
14
coordinates_min = (-1.0, -1.0) # minimum coordinates (min(x), min(y))
15
coordinates_max = (1.0, 1.0) # maximum coordinates (max(x), max(y))
16
17
trees_per_dimension = (8, 8)
18
19
# Create parent P4estMesh with 8 x 8 trees and 8 x 8 elements
20
# Since we couple through the boundaries, the periodicity does not matter here,
21
# but it is to trigger parts of the code for the test.
22
parent_mesh = P4estMesh(trees_per_dimension, polydeg = 3,
23
coordinates_min = coordinates_min,
24
coordinates_max = coordinates_max,
25
initial_refinement_level = 0,
26
periodicity = false)
27
28
# Define the mesh views consisting of a small square in the center
29
# and a square ring around it.
30
cell_ids1 = vcat((1:18), (23:26), (31:34), (39:42), (47:64))
31
mesh1 = P4estMeshView(parent_mesh, cell_ids1)
32
cell_ids2 = vcat((19:22), (27:30), (35:38), (43:46))
33
mesh2 = P4estMeshView(parent_mesh, cell_ids2)
34
35
# Define a trivial coupling function.
36
coupling_function = (x, u, equations_other, equations_own) -> u
37
38
# The mesh is coupled across the physical boundaries, which makes this setup
39
# effectively double periodic.
40
boundary_conditions = (; x_neg = BoundaryConditionCoupledP4est(coupling_function),
41
y_neg = BoundaryConditionCoupledP4est(coupling_function),
42
y_pos = BoundaryConditionCoupledP4est(coupling_function),
43
x_pos = BoundaryConditionCoupledP4est(coupling_function))
44
45
semi1 = SemidiscretizationHyperbolic(mesh1, equations, initial_condition_convergence_test,
46
solver,
47
boundary_conditions = boundary_conditions)
48
semi2 = SemidiscretizationHyperbolic(mesh2, equations, initial_condition_convergence_test,
49
solver,
50
boundary_conditions = boundary_conditions)
51
52
# Create a semidiscretization that bundles semi1 and semi2
53
semi = SemidiscretizationCoupledP4est(semi1, semi2)
54
55
###############################################################################
56
# ODE solvers, callbacks etc.
57
58
# Create ODE problem with time span from 0.0 to 2.0
59
ode = semidiscretize(semi, (0.0, 2.0))
60
61
# At the beginning of the main loop, the SummaryCallback prints a summary of the simulation setup
62
# and resets the timers
63
summary_callback = SummaryCallback()
64
65
# The AnalysisCallback allows to analyse the solution in regular intervals and prints the results
66
# We require this definition for the test, even though we don't use it in the CallbackSet.
67
analysis_callback1 = AnalysisCallback(semi1, interval = 100)
68
analysis_callback2 = AnalysisCallback(semi2, interval = 100)
69
analysis_callback = AnalysisCallbackCoupledP4est(semi, analysis_callback1,
70
analysis_callback2)
71
72
# The SaveSolutionCallback allows to save the solution to a file in regular intervals
73
save_solution = SaveSolutionCallback(interval = 100,
74
solution_variables = cons2prim)
75
76
# The StepsizeCallback handles the re-calculation of the maximum Δt after each time step
77
stepsize_callback = StepsizeCallback(cfl = 1.6)
78
79
# Create a CallbackSet to collect all callbacks such that they can be passed to the ODE solver
80
callbacks = CallbackSet(summary_callback, save_solution, stepsize_callback)
81
82
###############################################################################
83
# run the simulation
84
85
# OrdinaryDiffEq's `solve` method evolves the solution in time and executes the passed callbacks
86
sol = solve(ode, CarpenterKennedy2N54(williamson_condition = false);
87
dt = 1.0, # solve needs some value here but it will be overwritten by the stepsize_callback
88
ode_default_options()..., callback = callbacks);
89
90