1.2: First steps in Trixi.jl: Create your first setup

In this part of the introductory guide, we will create a first Trixi.jl setup as an extension of elixir_advection_basic.jl. Since Trixi.jl has a common basic structure for the setups, you can create your own by extending and modifying the following example.

Let's consider the linear advection equation for a state $u = u(x, y, t)$ on the two-dimensional spatial domain $[-1, 1] \times [-1, 1]$ with a source term

\[\frac{\partial}{\partial t}u + \frac{\partial}{\partial x} (0.2 u) - \frac{\partial}{\partial y} (0.7 u) = - 2 e^{-t} \sin\bigl(2 \pi (x - t) \bigr) \sin\bigl(2 \pi (y - t) \bigr),\]

with the initial condition

\[u(x, y, 0) = \sin\bigl(\pi x \bigr) \sin\bigl(\pi y \bigr),\]

and periodic boundary conditions.

The first step is to create and open a file with the .jl extension. You can do this with your favorite text editor (if you do not have one, we recommend VS Code). In this file, you will create your setup. The file can then be executed in Julia using, for example, trixi_include(). Alternatively, you can execute each line of the following code one by one in the Julia REPL. This will generate useful output for nearly every command and improve your comprehension of the process.

To be able to use functionalities of Trixi.jl, you always need to load Trixi.jl itself and the OrdinaryDiffEq.jl package.

using Trixi
using OrdinaryDiffEq

The next thing to do is to choose an equation that is suitable for your problem. To see all the currently implemented equations, take a look at src/equations. If you are interested in adding a new physics model that has not yet been implemented in Trixi.jl, take a look at the tutorials Adding a new scalar conservation law or Adding a non-conservative equation.

The linear scalar advection equation in two spatial dimensions

\[\frac{\partial}{\partial t}u + \frac{\partial}{\partial x} (a_1 u) + \frac{\partial}{\partial y} (a_2 u) = 0\]

is already implemented in Trixi.jl as LinearScalarAdvectionEquation2D, for which we need to define a two-dimensional parameter advection_velocity describing the parameters $a_1$ and $a_2$. Appropriate for our problem is (0.2, -0.7).

advection_velocity = (0.2, -0.7)
equations = LinearScalarAdvectionEquation2D(advection_velocity)
┌──────────────────────────────────────────────────────────────────────────────────────────────────┐
│ LinearScalarAdvectionEquation2D                                                                  │
│ ═══════════════════════════════                                                                  │
│ #variables: ………………………………………………… 1                                                                │
│ │ variable 1: …………………………………………… scalar                                                           │
└──────────────────────────────────────────────────────────────────────────────────────────────────┘

To solve our problem numerically using Trixi.jl, we have to discretize the spatial domain, for which we set up a mesh. One of the most used meshes in Trixi.jl is the TreeMesh. The spatial domain used is $[-1, 1] \times [-1, 1]$. We set an initial number of elements in the mesh using initial_refinement_level, which describes the initial number of hierarchical refinements. In this simple case, the total number of elements is 2^initial_refinement_level throughout the simulation. The variable n_cells_max is used to limit the number of elements in the mesh, which cannot be exceeded when using adaptive mesh refinement.

All minimum and all maximum coordinates must be combined into Tuples.

coordinates_min = (-1.0, -1.0)
coordinates_max = (1.0, 1.0)
mesh = TreeMesh(coordinates_min, coordinates_max,
                initial_refinement_level = 4,
                n_cells_max = 30_000)
┌──────────────────────────────────────────────────────────────────────────────────────────────────┐
│ TreeMesh{2, Trixi.SerialTree{2}}                                                                 │
│ ════════════════════════════════                                                                 │
│ center: …………………………………………………………… [0.0, 0.0]                                                       │
│ length: …………………………………………………………… 2.0                                                              │
│ periodicity: ……………………………………………… (true, true)                                                     │
│ current #cells: ……………………………………… 341                                                              │
│ #leaf-cells: ……………………………………………… 256                                                              │
│ maximum #cells: ……………………………………… 30000                                                            │
└──────────────────────────────────────────────────────────────────────────────────────────────────┘

To approximate the solution of the defined model, we create a DGSEM solver. The solution in each of the recently defined mesh elements will be approximated by a polynomial of degree polydeg. For more information about discontinuous Galerkin methods, check out the Introduction to DG methods tutorial. By default, in the weak formulation DGSEM initializes the surface flux as flux_central and uses the physical flux for the volume integral.

solver = DGSEM(polydeg = 3)
┌──────────────────────────────────────────────────────────────────────────────────────────────────┐
│ DG{Float64}                                                                                      │
│ ═══════════                                                                                      │
│ basis: ……………………………………………………………… LobattoLegendreBasis{Float64}(polydeg=3)                         │
│ mortar: …………………………………………………………… LobattoLegendreMortarL2{Float64}(polydeg=3)                      │
│ surface integral: ………………………………… SurfaceIntegralWeakForm                                          │
│ │ surface flux: ……………………………………… flux_central                                                     │
│ volume integral: …………………………………… VolumeIntegralWeakForm                                           │
└──────────────────────────────────────────────────────────────────────────────────────────────────┘

Now we need to define an initial condition for our problem. All the already implemented initial conditions for LinearScalarAdvectionEquation2D can be found in src/equations/linear_scalar_advection_2d.jl. If you want to use, for example, a Gaussian pulse, it can be used as follows:

initial_conditions = initial_condition_gauss

But to show you how an arbitrary initial condition can be implemented in a way suitable for Trixi.jl, we define our own initial conditions.

\[u(x, y, 0) = \sin\bigl(\pi x \bigr) \sin\bigl(\pi y \bigr).\]

The initial conditions function must take spatial coordinates, time and equation as arguments and returns an initial condition as a statically sized vector SVector. Following the same structure, you can define your own initial conditions. The time variable t can be unused in the initial condition, but might also be used to describe an analytical solution if known. If you use the initial condition as analytical solution, you can analyze your numerical solution by computing the error, see also the section about analyzing the solution.

function initial_condition_sinpi(x, t, equations::LinearScalarAdvectionEquation2D)
    u = sinpi(x[1]) * sinpi(x[2])
    return SVector(u)
end
initial_condition = initial_condition_sinpi
initial_condition_sinpi (generic function with 1 method)

The next step is to define a function of the source term corresponding to our problem.

\[f(u, x, y, t) = - 2 e^{-t} \sin\bigl(2 \pi (x - t) \bigr) \sin\bigl(2 \pi (y - t) \bigr)\]

This function must take the state variable, the spatial coordinates, the time and the equation itself as arguments and returns the source term as a static vector SVector.

function source_term_exp_sinpi(u, x, t, equations::LinearScalarAdvectionEquation2D)
    u = -2 * exp(-t) * sinpi(2 * (x[1] - t)) * sinpi(2 * (x[2] - t))
    return SVector(u)
end
source_term_exp_sinpi (generic function with 1 method)

Now we collect all the information that is necessary to define a spatial discretization,

semi = SemidiscretizationHyperbolic(mesh, equations, initial_condition, solver;
                                    source_terms = source_term_exp_sinpi)
┌──────────────────────────────────────────────────────────────────────────────────────────────────┐
│ SemidiscretizationHyperbolic                                                                     │
│ ════════════════════════════                                                                     │
│ #spatial dimensions: ………………………… 2                                                                │
│ mesh: ………………………………………………………………… TreeMesh{2, Trixi.SerialTree{2}} with length 341                 │
│ equations: …………………………………………………… LinearScalarAdvectionEquation2D                                  │
│ initial condition: ……………………………… initial_condition_sinpi                                          │
│ boundary conditions: ………………………… Trixi.BoundaryConditionPeriodic                                  │
│ source terms: …………………………………………… source_term_exp_sinpi                                            │
│ solver: …………………………………………………………… DG                                                               │
│ total #DOFs per field: …………………… 4096                                                             │
└──────────────────────────────────────────────────────────────────────────────────────────────────┘

which leaves us with an ODE problem in time with a span from 0.0 to 1.0. This approach is commonly referred to as the method of lines.

tspan = (0.0, 1.0)
ode = semidiscretize(semi, tspan)
ODEProblem with uType Vector{Float64} and tType Float64. In-place: true
timespan: (0.0, 1.0)
u0: 4096-element Vector{Float64}:
 0.0
 0.0
 0.0
 0.0
 0.0
 0.011734602208786986
 0.030369402873034978
 0.04145470668348209
 0.0
 0.030369402873034978
 ⋮
 0.0
 0.04145470668348209
 0.030369402873034978
 0.011734602208786986
 0.0
 0.0
 0.0
 0.0
 0.0

At this point, our problem is defined. We will use the solve function defined in OrdinaryDiffEq.jl to get the solution. OrdinaryDiffEq.jl gives us the ability to customize the solver using callbacks without actually modifying it. Trixi.jl already has some implemented Callbacks. The most widely used callbacks in Trixi.jl are step control callbacks that are activated at the end of each time step to perform some actions, e.g. to print statistics. We will show you how to use some of the common callbacks.

To print a summary of the simulation setup at the beginning and to reset timers to zero, we use the SummaryCallback.

summary_callback = SummaryCallback()
SummaryCallback

We also want to analyze the current state of the solution in regular intervals. The AnalysisCallback outputs some useful statistical information during the simulation every interval time steps.

analysis_callback = AnalysisCallback(semi, interval = 20)
┌──────────────────────────────────────────────────────────────────────────────────────────────────┐
│ AnalysisCallback                                                                                 │
│ ════════════════                                                                                 │
│ interval: ……………………………………………………… 20                                                               │
│ analyzer: ……………………………………………………… LobattoLegendreAnalyzer{Float64}(polydeg=6)                      │
│ │ error 1: …………………………………………………… l2_error                                                         │
│ │ error 2: …………………………………………………… linf_error                                                       │
│ │ integral 1: …………………………………………… entropy_timederivative                                           │
│ save analysis to file: …………………… no                                                               │
└──────────────────────────────────────────────────────────────────────────────────────────────────┘

To indicate that a simulation is still running, we utilize the inexpensive AliveCallback to periodically print information to the screen, such as the current time, every alive_interval time steps.

alive_callback = AliveCallback(alive_interval = 10)
┌──────────────────────────────────────────────────────────────────────────────────────────────────┐
│ AliveCallback                                                                                    │
│ ═════════════                                                                                    │
│ interval: ……………………………………………………… 10                                                               │
└──────────────────────────────────────────────────────────────────────────────────────────────────┘

It is also possible to control the time step size using the StepsizeCallback if the time integration method isn't adaptive itself. To get more details, look at CFL based step size control.

stepsize_callback = StepsizeCallback(cfl = 0.9)
┌──────────────────────────────────────────────────────────────────────────────────────────────────┐
│ StepsizeCallback                                                                                 │
│ ════════════════                                                                                 │
│ CFL number: ………………………………………………… 0.9                                                              │
└──────────────────────────────────────────────────────────────────────────────────────────────────┘

To save the current solution in regular intervals we use the SaveSolutionCallback. We would like to save the initial and final solutions as well. The data will be saved as HDF5 files located in the out folder. Afterwards it is possible to visualize a solution from saved files using Trixi2Vtk.jl and ParaView, which is described below in the section Visualize the solution.

save_solution = SaveSolutionCallback(interval = 20,
                                     save_initial_solution = true,
                                     save_final_solution = true)
┌──────────────────────────────────────────────────────────────────────────────────────────────────┐
│ SaveSolutionCallback                                                                             │
│ ════════════════════                                                                             │
│ interval: ……………………………………………………… 20                                                               │
│ solution variables: …………………………… cons2prim                                                        │
│ save initial solution: …………………… yes                                                              │
│ save final solution: ………………………… yes                                                              │
│ output directory: ………………………………… /home/runner/work/Trixi.jl/Trixi.jl/docs/out                     │
└──────────────────────────────────────────────────────────────────────────────────────────────────┘

Alternatively, we have the option to print solution files at fixed time intervals.

save_solution = SaveSolutionCallback(dt = 0.1,
                                     save_initial_solution = true,
                                     save_final_solution = true)

Another useful callback is the SaveRestartCallback. It saves information for restarting in regular intervals. We are interested in saving a restart file for the final solution as well. To perform a restart, you need to configure the restart setup in a special way, which is described in the section Restart simulation.

save_restart = SaveRestartCallback(interval = 100, save_final_restart = true)
┌──────────────────────────────────────────────────────────────────────────────────────────────────┐
│ SaveRestartCallback                                                                              │
│ ═══════════════════                                                                              │
│ interval: ……………………………………………………… 100                                                              │
│ save final solution: ………………………… yes                                                              │
│ output directory: ………………………………… /home/runner/work/Trixi.jl/Trixi.jl/docs/out                     │
└──────────────────────────────────────────────────────────────────────────────────────────────────┘

Create a CallbackSet to collect all callbacks so that they can be passed to the solve function.

callbacks = CallbackSet(summary_callback, analysis_callback, alive_callback,
                        stepsize_callback,
                        save_solution, save_restart);

The last step is to choose the time integration method. OrdinaryDiffEq.jl defines a wide range of ODE solvers, including the three-stage, third-order strong stability preserving Runge-Kutta method SSPRK33. We will pass the ODE problem, the ODE solver and the callbacks to the solve function. Also, to use StepsizeCallback, we must explicitly specify the initial trial time step dt, the selected value is not important, because it will be overwritten by the StepsizeCallback. And there is no need to save every step of the solution, as we are only interested the output provided by our callback SaveSolutionCallback.

sol = solve(ode, SSPRK33(); dt = 1.0, save_everystep = false, callback = callbacks);

████████╗██████╗ ██╗██╗  ██╗██╗
╚══██╔══╝██╔══██╗██║╚██╗██╔╝██║
   ██║   ██████╔╝██║ ╚███╔╝ ██║
   ██║   ██╔══██╗██║ ██╔██╗ ██║
   ██║   ██║  ██║██║██╔╝ ██╗██║
   ╚═╝   ╚═╝  ╚═╝╚═╝╚═╝  ╚═╝╚═╝

┌──────────────────────────────────────────────────────────────────────────────────────────────────┐
│ SemidiscretizationHyperbolic                                                                     │
│ ════════════════════════════                                                                     │
│ #spatial dimensions: ………………………… 2                                                                │
│ mesh: ………………………………………………………………… TreeMesh{2, Trixi.SerialTree{2}} with length 341                 │
│ equations: …………………………………………………… LinearScalarAdvectionEquation2D                                  │
│ initial condition: ……………………………… initial_condition_sinpi                                          │
│ boundary conditions: ………………………… Trixi.BoundaryConditionPeriodic                                  │
│ source terms: …………………………………………… source_term_exp_sinpi                                            │
│ solver: …………………………………………………………… DG                                                               │
│ total #DOFs per field: …………………… 4096                                                             │
└──────────────────────────────────────────────────────────────────────────────────────────────────┘

┌──────────────────────────────────────────────────────────────────────────────────────────────────┐
│ TreeMesh{2, Trixi.SerialTree{2}}                                                                 │
│ ════════════════════════════════                                                                 │
│ center: …………………………………………………………… [0.0, 0.0]                                                       │
│ length: …………………………………………………………… 2.0                                                              │
│ periodicity: ……………………………………………… (true, true)                                                     │
│ current #cells: ……………………………………… 341                                                              │
│ #leaf-cells: ……………………………………………… 256                                                              │
│ maximum #cells: ……………………………………… 30000                                                            │
└──────────────────────────────────────────────────────────────────────────────────────────────────┘

┌──────────────────────────────────────────────────────────────────────────────────────────────────┐
│ LinearScalarAdvectionEquation2D                                                                  │
│ ═══════════════════════════════                                                                  │
│ #variables: ………………………………………………… 1                                                                │
│ │ variable 1: …………………………………………… scalar                                                           │
└──────────────────────────────────────────────────────────────────────────────────────────────────┘

┌──────────────────────────────────────────────────────────────────────────────────────────────────┐
│ DG{Float64}                                                                                      │
│ ═══════════                                                                                      │
│ basis: ……………………………………………………………… LobattoLegendreBasis{Float64}(polydeg=3)                         │
│ mortar: …………………………………………………………… LobattoLegendreMortarL2{Float64}(polydeg=3)                      │
│ surface integral: ………………………………… SurfaceIntegralWeakForm                                          │
│ │ surface flux: ……………………………………… flux_central                                                     │
│ volume integral: …………………………………… VolumeIntegralWeakForm                                           │
└──────────────────────────────────────────────────────────────────────────────────────────────────┘

┌──────────────────────────────────────────────────────────────────────────────────────────────────┐
│ AnalysisCallback                                                                                 │
│ ════════════════                                                                                 │
│ interval: ……………………………………………………… 20                                                               │
│ analyzer: ……………………………………………………… LobattoLegendreAnalyzer{Float64}(polydeg=6)                      │
│ │ error 1: …………………………………………………… l2_error                                                         │
│ │ error 2: …………………………………………………… linf_error                                                       │
│ │ integral 1: …………………………………………… entropy_timederivative                                           │
│ save analysis to file: …………………… no                                                               │
└──────────────────────────────────────────────────────────────────────────────────────────────────┘

┌──────────────────────────────────────────────────────────────────────────────────────────────────┐
│ AliveCallback                                                                                    │
│ ═════════════                                                                                    │
│ interval: ……………………………………………………… 10                                                               │
└──────────────────────────────────────────────────────────────────────────────────────────────────┘

┌──────────────────────────────────────────────────────────────────────────────────────────────────┐
│ StepsizeCallback                                                                                 │
│ ════════════════                                                                                 │
│ CFL number: ………………………………………………… 0.9                                                              │
└──────────────────────────────────────────────────────────────────────────────────────────────────┘

┌──────────────────────────────────────────────────────────────────────────────────────────────────┐
│ SaveSolutionCallback                                                                             │
│ ════════════════════                                                                             │
│ interval: ……………………………………………………… 20                                                               │
│ solution variables: …………………………… cons2prim                                                        │
│ save initial solution: …………………… yes                                                              │
│ save final solution: ………………………… yes                                                              │
│ output directory: ………………………………… /home/runner/work/Trixi.jl/Trixi…build/tutorials/first_steps/out │
└──────────────────────────────────────────────────────────────────────────────────────────────────┘

┌──────────────────────────────────────────────────────────────────────────────────────────────────┐
│ SaveRestartCallback                                                                              │
│ ═══════════════════                                                                              │
│ interval: ……………………………………………………… 100                                                              │
│ save final solution: ………………………… yes                                                              │
│ output directory: ………………………………… /home/runner/work/Trixi.jl/Trixi…build/tutorials/first_steps/out │
└──────────────────────────────────────────────────────────────────────────────────────────────────┘

┌──────────────────────────────────────────────────────────────────────────────────────────────────┐
│ Time integration                                                                                 │
│ ════════════════                                                                                 │
│ Start time: ………………………………………………… 0.0                                                              │
│ Final time: ………………………………………………… 1.0                                                              │
│ time integrator: …………………………………… SSPRK33                                                          │
│ adaptive: ……………………………………………………… false                                                            │
└──────────────────────────────────────────────────────────────────────────────────────────────────┘
┌──────────────────────────────────────────────────────────────────────────────────────────────────┐
│ Environment information                                                                          │
│ ═══════════════════════                                                                          │
│ #threads: ……………………………………………………… 1                                                                │
└──────────────────────────────────────────────────────────────────────────────────────────────────┘

────────────────────────────────────────────────────────────────────────────────────────────────────
 Simulation running 'LinearScalarAdvectionEquation2D' with DGSEM(polydeg=3)
────────────────────────────────────────────────────────────────────────────────────────────────────
 #timesteps:                  0                run time:       1.02200000e-06 s
 Δt:             1.00000000e+00                └── GC time:    0.00000000e+00 s (0.000%)
 sim. time:      0.00000000e+00 (0.000%)       time/DOF/rhs!:         NaN s
                                               PID:                   Inf s
 #DOFs per field:          4096                alloc'd memory:       1120.748 MiB
 #elements:                 256

 Variable:       scalar
 L2 error:       5.57371016e-06
 Linf error:     2.37926999e-05
 ∑∂S/∂U ⋅ Uₜ :   4.42993996e-17
────────────────────────────────────────────────────────────────────────────────────────────────────

#timesteps:     10 │ Δt: 3.1250e-02 │ sim. time: 3.1250e-01 (31.250%)  │ run time: 4.4068e-02 s

────────────────────────────────────────────────────────────────────────────────────────────────────
 Simulation running 'LinearScalarAdvectionEquation2D' with DGSEM(polydeg=3)
────────────────────────────────────────────────────────────────────────────────────────────────────
 #timesteps:                 20                run time:       5.61696650e-02 s
 Δt:             3.12500000e-02                └── GC time:    0.00000000e+00 s (0.000%)
 sim. time:      6.25000000e-01 (62.500%)      time/DOF/rhs!:  3.76525158e-08 s
                                               PID:            2.24665344e-07 s
 #DOFs per field:          4096                alloc'd memory:       1125.192 MiB
 #elements:                 256

 Variable:       scalar
 L2 error:       6.70015343e-01
 Linf error:     1.45520916e+00
 ∑∂S/∂U ⋅ Uₜ :  -3.74390676e-02
────────────────────────────────────────────────────────────────────────────────────────────────────

#timesteps:     20 │ Δt: 3.1250e-02 │ sim. time: 6.2500e-01 (62.500%)  │ run time: 4.9692e-02 s
#timesteps:     30 │ Δt: 3.1250e-02 │ sim. time: 9.3750e-01 (93.750%)  │ run time: 5.4877e-02 s

────────────────────────────────────────────────────────────────────────────────────────────────────
 Simulation running 'LinearScalarAdvectionEquation2D' with DGSEM(polydeg=3)
────────────────────────────────────────────────────────────────────────────────────────────────────
 #timesteps:                 32                run time:       6.32636970e-02 s
 Δt:             3.12500000e-02                └── GC time:    0.00000000e+00 s (0.000%)
 sim. time:      1.00000000e+00 (100.000%)     time/DOF/rhs!:  3.72549752e-08 s
                                               PID:            4.19028592e-08 s
 #DOFs per field:          4096                alloc'd memory:       1129.434 MiB
 #elements:                 256

 Variable:       scalar
 L2 error:       8.65828105e-01
 Linf error:     1.75847551e+00
 ∑∂S/∂U ⋅ Uₜ :  -1.70590047e-02
────────────────────────────────────────────────────────────────────────────────────────────────────

────────────────────────────────────────────────────────────────────────────────────────────────────
Trixi.jl simulation finished.  Final time: 1.0  Time steps: 32 (accepted), 32 (total)
────────────────────────────────────────────────────────────────────────────────────────────────────

Finally, we print the timer summary.

summary_callback()
────────────────────────────────────────────────────────────────────────────────
           Trixi.jl                    Time                    Allocations
                              ───────────────────────   ────────────────────────
      Tot / % measured:           73.1ms /  27.6%           15.1MiB /  77.0%

Section               ncalls     time    %tot     avg     alloc    %tot      avg
────────────────────────────────────────────────────────────────────────────────
rhs!                      96   14.7ms   72.8%   153μs   9.33KiB    0.1%     100B
  source terms            96   8.90ms   44.1%  92.7μs     0.00B    0.0%    0.00B
  volume integral         96   3.80ms   18.8%  39.6μs     0.00B    0.0%    0.00B
  interface flux          96   1.09ms    5.4%  11.4μs     0.00B    0.0%    0.00B
  surface integral        96    377μs    1.9%  3.93μs     0.00B    0.0%    0.00B
  prolong2interfaces      96    319μs    1.6%  3.32μs     0.00B    0.0%    0.00B
  ~rhs!~                  96   95.4μs    0.5%   993ns   9.33KiB    0.1%     100B
  Jacobian                96   55.4μs    0.3%   577ns     0.00B    0.0%    0.00B
  reset ∂u/∂t             96   46.2μs    0.2%   481ns     0.00B    0.0%    0.00B
  prolong2mortars         96   5.00μs    0.0%  52.1ns     0.00B    0.0%    0.00B
  prolong2boundaries      96   4.32μs    0.0%  45.0ns     0.00B    0.0%    0.00B
  mortar flux             96   4.27μs    0.0%  44.5ns     0.00B    0.0%    0.00B
  boundary flux           96   2.95μs    0.0%  30.8ns     0.00B    0.0%    0.00B
I/O                        6   2.76ms   13.7%   460μs    371KiB    3.1%  61.8KiB
  save solution            3   1.59ms    7.9%   530μs    302KiB    2.5%   101KiB
  ~I/O~                    6   1.17ms    5.8%   194μs   68.5KiB    0.6%  11.4KiB
  get element vari...      3    661ns    0.0%   220ns     0.00B    0.0%    0.00B
  save mesh                3    271ns    0.0%  90.3ns     0.00B    0.0%    0.00B
  get node variables       3    110ns    0.0%  36.7ns     0.00B    0.0%    0.00B
analyze solution           3   2.73ms   13.5%   909μs   11.3MiB   96.8%  3.76MiB
calculate dt              33   3.85μs    0.0%   117ns     0.00B    0.0%    0.00B
────────────────────────────────────────────────────────────────────────────────

Now you can plot the solution as shown below, analyze it and improve the stability, accuracy or efficiency of your setup.

Visualize the solution

In the previous part of the tutorial, we calculated the final solution of the given problem, now we want to visualize it. A more detailed explanation of visualization methods can be found in the section Visualization.

Using Plots.jl

The first option is to use the Plots.jl package directly after calculations, when the solution is saved in the sol variable.

using Plots

As was shown in the Getting started section, you can plot all variables from the system of equations by executing the following.

plot(sol)

Alternatively, you can configure the plot more precisely. Trixi.jl provides a special data type, PlotData2D, to extract the visualization data from the solution.

pd = PlotData2D(sol);

You can plot specific variables from the system of equations by referring to their names. To obtain the names of all variables, execute the following.

@show pd.variable_names;
pd.variable_names = ["scalar"]

Plot the variable named "scalar" (which is the name of the variable for the linear advection equation in Trixi.jl).

plot(pd["scalar"])
Example block output

Mesh extraction is possible using the getmesh function. Plots.jl has the plot! function that allows you to modify an already built graph.

plot!(getmesh(pd))
Example block output

Using Trixi2Vtk.jl

Another way to visualize a solution is to extract it from a saved HDF5 file. After we used the solve function with SaveSolutionCallback there is a file with the final solution. It is located in the out folder and is named as follows: solution_index.h5. The index is the final time step of the solution that is padded to 6 digits with zeros from the beginning. With Trixi2Vtk you can convert the HDF5 output file generated by Trixi.jl into a VTK/VTU files. VTK/VTU are specialized formats designed to store structured data required for visualization purposes. This can be used in visualization tools such as ParaView or VisIt to plot the solution.

If you haven't added Trixi2Vtk.jl to your project yet, you can add it as follows.

import Pkg
Pkg.add(["Trixi2Vtk"])

Now we load the Trixi2Vtk.jl package and convert the file out/solution_000000032.h5 with the final solution using the trixi2vtk function saving the resulting file in the out folder.

using Trixi2Vtk
trixi2vtk(joinpath("out", "solution_000000032.h5"), output_directory = "out")
────────────────────────────────────────────────────────────────────────────────
                                       Time                    Allocations
                              ───────────────────────   ────────────────────────
      Tot / % measured:            1.15s /  62.6%           65.1MiB /  72.1%

Section               ncalls     time    %tot     avg     alloc    %tot      avg
────────────────────────────────────────────────────────────────────────────────
interpolate data           1    311ms   43.1%   311ms   7.76MiB   16.5%  7.76MiB
add data to VTK file       1    274ms   38.0%   274ms   9.09MiB   19.3%  9.09MiB
  scalar                   1    113ms   15.6%   113ms   5.50MiB   11.7%  5.50MiB
  add data to VTK ...      1    526μs    0.1%   526μs    118KiB    0.2%   118KiB
    cell_ids               1    267μs    0.0%   267μs   35.6KiB    0.1%  35.6KiB
    element_ids            1    197μs    0.0%   197μs   37.8KiB    0.1%  37.8KiB
    levels                 1   51.5μs    0.0%  51.5μs   37.8KiB    0.1%  37.8KiB
read mesh                  1   76.1ms   10.5%  76.1ms   9.00MiB   19.2%  9.00MiB
build VTK grid (no...      1   51.7ms    7.2%  51.7ms   15.6MiB   33.2%  15.6MiB
prepare VTK cells ...      1   7.03ms    1.0%  7.03ms   5.04MiB   10.7%  5.04MiB
build VTK grid (ce...      1    915μs    0.1%   915μs    353KiB    0.7%   353KiB
save VTK file              2    492μs    0.1%   246μs   4.91KiB    0.0%  2.45KiB
read data                  1    354μs    0.0%   354μs   69.6KiB    0.1%  69.6KiB
prepare VTK cells ...      1   99.4μs    0.0%  99.4μs    100KiB    0.2%   100KiB
────────────────────────────────────────────────────────────────────────────────

Now two files solution_000000032.vtu and solution_000000032_celldata.vtu have been generated in the out folder. The first one contains all the information for visualizing the solution, the second one contains all the cell-based or discretization-based information.

Now let's visualize the solution from the generated files in ParaView. Follow this short instruction to get the visualization.

  • Download, install and open ParaView.
  • Press Ctrl+O and select the generated files solution_000000032.vtu and solution_000000032_celldata.vtu from the out folder.
  • In the upper-left corner in the Pipeline Browser window, left-click on the eye-icon near solution_000000032.vtu.
  • In the lower-left corner in the Properties window, change the Coloring from Solid Color to scalar. This already generates the visualization of the final solution.
  • Now let's add the mesh to the visualization. In the upper-left corner in the Pipeline Browser window, left-click on the eye-icon near solution_000000032_celldata.vtu.
  • In the lower-left corner in the Properties window, change the Representation from Surface to Wireframe. Then a white grid should appear on the visualization.

Now, if you followed the instructions exactly, you should get a similar image as shown in the section Using Plots.jl:

paraview_trixi2vtk_example

After completing this tutorial you are able to set up your own simulations with Trixi.jl. If you have an interest in contributing to Trixi.jl as a developer, refer to the third part of the introduction titled Changing Trixi.jl itself.


This page was generated using Literate.jl.