Path: blob/main/examples/dgmulti_2d/elixir_euler_triangulate_pkg_mesh.jl
2055 views
using OrdinaryDiffEqLowStorageRK1using Trixi23dg = DGMulti(polydeg = 3, element_type = Tri(),4surface_integral = SurfaceIntegralWeakForm(flux_hll),5volume_integral = VolumeIntegralWeakForm())67equations = CompressibleEulerEquations2D(1.4)8initial_condition = initial_condition_convergence_test9source_terms = source_terms_convergence_test1011# use pre-defined Triangulate geometry in StartUpDG12meshIO = StartUpDG.triangulate_domain(StartUpDG.RectangularDomainWithHole())1314# the pre-defined Triangulate geometry in StartUpDG has integer boundary tags. this routine15# assigns boundary faces based on these integer boundary tags.16mesh = DGMultiMesh(dg, meshIO, Dict(:outer_boundary => 1, :inner_boundary => 2))1718boundary_condition_convergence_test = BoundaryConditionDirichlet(initial_condition)19boundary_conditions = (; :outer_boundary => boundary_condition_convergence_test,20:inner_boundary => boundary_condition_convergence_test)2122semi = SemidiscretizationHyperbolic(mesh, equations, initial_condition, dg,23source_terms = source_terms,24boundary_conditions = boundary_conditions)2526tspan = (0.0, 0.2)27ode = semidiscretize(semi, tspan)2829summary_callback = SummaryCallback()30alive_callback = AliveCallback(alive_interval = 10)31analysis_interval = 10032analysis_callback = AnalysisCallback(semi, interval = analysis_interval, uEltype = real(dg))33save_solution = SaveSolutionCallback(interval = analysis_interval,34solution_variables = cons2prim)35callbacks = CallbackSet(summary_callback, alive_callback, analysis_callback, save_solution)3637###############################################################################38# run the simulation3940sol = solve(ode, CarpenterKennedy2N54(williamson_condition = false);41dt = 0.5 * estimate_dt(mesh, dg),42ode_default_options()...,43callback = callbacks);444546