{ "cells": [ { "cell_type": "markdown", "source": [ "# 1.2: First steps in Trixi.jl: Create your first setup" ], "metadata": {} }, { "cell_type": "markdown", "source": [ "**Note:** To improve responsiveness via caching, the notebooks are updated only once a week. They are only\n", "available for the latest stable release of Trixi.jl at the time of caching." ], "metadata": {} }, { "cell_type": "markdown", "source": [ "In this part of the introductory guide, we will create a first Trixi.jl setup as an extension of\n", "[`elixir_advection_basic.jl`](https://github.com/trixi-framework/Trixi.jl/blob/main/examples/tree_2d_dgsem/elixir_advection_basic.jl).\n", "Since Trixi.jl has a common basic structure for the setups, you can create your own by extending\n", "and modifying the following example." ], "metadata": {} }, { "cell_type": "markdown", "source": [ "Let's consider the linear advection equation for a state $u = u(x, y, t)$ on the two-dimensional spatial domain\n", "$[-1, 1] \\times [-1, 1]$ with a source term\n", "$$\n", "\\frac{\\partial}{\\partial t}u + \\frac{\\partial}{\\partial x} (0.2 u) - \\frac{\\partial}{\\partial y} (0.7 u) = - 2 e^{-t}\n", "\\sin\\bigl(2 \\pi (x - t) \\bigr) \\sin\\bigl(2 \\pi (y - t) \\bigr),\n", "$$\n", "with the initial condition\n", "$$\n", "u(x, y, 0) = \\sin\\bigl(\\pi x \\bigr) \\sin\\bigl(\\pi y \\bigr),\n", "$$\n", "and periodic boundary conditions." ], "metadata": {} }, { "cell_type": "markdown", "source": [ "The first step is to create and open a file with the .jl extension. You can do this with your\n", "favorite text editor (if you do not have one, we recommend [VS Code](https://code.visualstudio.com/)).\n", "In this file, you will create your setup. The file can then be executed in Julia using, for example, `trixi_include()`.\n", "Alternatively, you can execute each line of the following code one by one in the\n", "Julia REPL. This will generate useful output for nearly every\n", "command and improve your comprehension of the process." ], "metadata": {} }, { "cell_type": "markdown", "source": [ "To be able to use functionalities of Trixi.jl, you always need to load Trixi.jl itself\n", "and the [OrdinaryDiffEq.jl](https://github.com/SciML/OrdinaryDiffEq.jl) package." ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "using Trixi\n", "using OrdinaryDiffEq" ], "metadata": {}, "execution_count": null }, { "cell_type": "markdown", "source": [ "The next thing to do is to choose an equation that is suitable for your problem. To see all the\n", "currently implemented equations, take a look at\n", "[`src/equations`](https://github.com/trixi-framework/Trixi.jl/tree/main/src/equations).\n", "If you are interested in adding a new physics model that has not yet been implemented in\n", "Trixi.jl, take a look at the tutorials\n", "Adding a new scalar conservation law or\n", "Adding a non-conservative equation." ], "metadata": {} }, { "cell_type": "markdown", "source": [ "The linear scalar advection equation in two spatial dimensions\n", "$$\n", "\\frac{\\partial}{\\partial t}u + \\frac{\\partial}{\\partial x} (a_1 u) + \\frac{\\partial}{\\partial y} (a_2 u) = 0\n", "$$\n", "is already implemented in Trixi.jl as\n", "`LinearScalarAdvectionEquation2D`, for which we need to define a two-dimensional parameter\n", "`advection_velocity` describing the parameters $a_1$ and $a_2$. Appropriate for our problem is `(0.2, -0.7)`." ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "advection_velocity = (0.2, -0.7)\n", "equations = LinearScalarAdvectionEquation2D(advection_velocity)" ], "metadata": {}, "execution_count": null }, { "cell_type": "markdown", "source": [ "To solve our problem numerically using Trixi.jl, we have to discretize the spatial\n", "domain, for which we set up a mesh. One of the most used meshes in Trixi.jl is the\n", "`TreeMesh`. The spatial domain used is $[-1, 1] \\times [-1, 1]$. We set an initial number\n", "of elements in the mesh using `initial_refinement_level`, which describes the initial number of\n", "hierarchical refinements. In this simple case, the total number of elements is `2^initial_refinement_level`\n", "throughout the simulation. The variable `n_cells_max` is used to limit the number of elements in the mesh,\n", "which cannot be exceeded when using adaptive mesh refinement." ], "metadata": {} }, { "cell_type": "markdown", "source": [ "All minimum and all maximum coordinates must be combined into `Tuples`." ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "coordinates_min = (-1.0, -1.0)\n", "coordinates_max = ( 1.0, 1.0)\n", "mesh = TreeMesh(coordinates_min, coordinates_max,\n", " initial_refinement_level = 4,\n", " n_cells_max = 30_000)" ], "metadata": {}, "execution_count": null }, { "cell_type": "markdown", "source": [ "To approximate the solution of the defined model, we create a `DGSEM` solver.\n", "The solution in each of the recently defined mesh elements will be approximated by a polynomial\n", "of degree `polydeg`. For more information about discontinuous Galerkin methods,\n", "check out the Introduction to DG methods tutorial. By default,\n", "in the weak formulation `DGSEM` initializes the surface flux as `flux_central` and uses the physical flux for\n", "the volume integral." ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "solver = DGSEM(polydeg=3)" ], "metadata": {}, "execution_count": null }, { "cell_type": "markdown", "source": [ "Now we need to define an initial condition for our problem. All the already implemented\n", "initial conditions for `LinearScalarAdvectionEquation2D` can be found in\n", "[`src/equations/linear_scalar_advection_2d.jl`](https://github.com/trixi-framework/Trixi.jl/blob/main/src/equations/linear_scalar_advection_2d.jl).\n", "If you want to use, for example, a Gaussian pulse, it can be used as follows:\n", "```julia\n", "initial_conditions = initial_condition_gauss\n", "```\n", "But to show you how an arbitrary initial condition can be implemented in a way suitable for\n", "Trixi.jl, we define our own initial conditions.\n", "$$\n", "u(x, y, 0) = \\sin\\bigl(\\pi x \\bigr) \\sin\\bigl(\\pi y \\bigr).\n", "$$\n", "The initial conditions function must take spatial coordinates, time and equation as arguments\n", "and returns an initial condition as a statically sized vector `SVector`. Following the same structure, you\n", "can define your own initial conditions. The time variable `t` can be unused in the initial\n", "condition, but might also be used to describe an analytical solution if known. If you use the\n", "initial condition as analytical solution, you can analyze your numerical solution by computing\n", "the error, see also the\n", "[section about analyzing the solution](https://trixi-framework.github.io/Trixi.jl/stable/callbacks/#Analyzing-the-numerical-solution)." ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "function initial_condition_sinpi(x, t, equations::LinearScalarAdvectionEquation2D)\n", " u = sinpi(x[1]) * sinpi(x[2])\n", " return SVector(u)\n", "end\n", "initial_condition = initial_condition_sinpi" ], "metadata": {}, "execution_count": null }, { "cell_type": "markdown", "source": [ "The next step is to define a function of the source term corresponding to our problem.\n", "$$\n", "f(u, x, y, t) = - 2 e^{-t} \\sin\\bigl(2 \\pi (x - t) \\bigr) \\sin\\bigl(2 \\pi (y - t) \\bigr)\n", "$$\n", "This function must take the state variable, the spatial coordinates, the time and the\n", "equation itself as arguments and returns the source term as a static vector `SVector`." ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "function source_term_exp_sinpi(u, x, t, equations::LinearScalarAdvectionEquation2D)\n", " u = - 2 * exp(-t) * sinpi(2*(x[1] - t)) * sinpi(2*(x[2] - t))\n", " return SVector(u)\n", "end" ], "metadata": {}, "execution_count": null }, { "cell_type": "markdown", "source": [ "Now we collect all the information that is necessary to define a spatial discretization," ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "semi = SemidiscretizationHyperbolic(mesh, equations, initial_condition, solver;\n", " source_terms = source_term_exp_sinpi)" ], "metadata": {}, "execution_count": null }, { "cell_type": "markdown", "source": [ "which leaves us with an ODE problem in time with a span from `0.0` to `1.0`.\n", "This approach is commonly referred to as the method of lines." ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "tspan = (0.0, 1.0)\n", "ode = semidiscretize(semi, tspan)" ], "metadata": {}, "execution_count": null }, { "cell_type": "markdown", "source": [ "At this point, our problem is defined. We will use the `solve` function defined in\n", "[OrdinaryDiffEq.jl](https://github.com/SciML/OrdinaryDiffEq.jl) to get the solution.\n", "OrdinaryDiffEq.jl gives us the ability to customize the solver\n", "using callbacks without actually modifying it. Trixi.jl already has some implemented\n", "Callbacks. The most widely used callbacks in Trixi.jl are\n", "[step control callbacks](https://docs.sciml.ai/DiffEqCallbacks/stable/step_control/) that are\n", "activated at the end of each time step to perform some actions, e.g. to print statistics.\n", "We will show you how to use some of the common callbacks." ], "metadata": {} }, { "cell_type": "markdown", "source": [ "To print a summary of the simulation setup at the beginning\n", "and to reset timers to zero, we use the `SummaryCallback`." ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "summary_callback = SummaryCallback()" ], "metadata": {}, "execution_count": null }, { "cell_type": "markdown", "source": [ "We also want to analyze the current state of the solution in regular intervals.\n", "The `AnalysisCallback` outputs some useful statistical information during the simulation\n", "every `interval` time steps." ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "analysis_callback = AnalysisCallback(semi, interval = 20)" ], "metadata": {}, "execution_count": null }, { "cell_type": "markdown", "source": [ "To indicate that a simulation is still running, we utilize the inexpensive `AliveCallback`\n", "to periodically print information to the screen, such as the\n", "current time, every `alive_interval` time steps." ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "alive_callback = AliveCallback(alive_interval = 10)" ], "metadata": {}, "execution_count": null }, { "cell_type": "markdown", "source": [ "It is also possible to control the time step size using the `StepsizeCallback` if the time\n", "integration method isn't adaptive itself. To get more details, look at\n", "CFL based step size control." ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "stepsize_callback = StepsizeCallback(cfl = 0.9)" ], "metadata": {}, "execution_count": null }, { "cell_type": "markdown", "source": [ "To save the current solution in regular intervals we use the `SaveSolutionCallback`.\n", "We would like to save the initial and final solutions as well. The data\n", "will be saved as HDF5 files located in the `out` folder. Afterwards it is possible to visualize\n", "a solution from saved files using Trixi2Vtk.jl and ParaView, which is described below in the\n", "section Visualize the solution." ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "save_solution = SaveSolutionCallback(interval = 20,\n", " save_initial_solution = true,\n", " save_final_solution = true)" ], "metadata": {}, "execution_count": null }, { "cell_type": "markdown", "source": [ "Alternatively, we have the option to print solution files at fixed time intervals.\n", "```julua\n", "save_solution = SaveSolutionCallback(dt = 0.1,\n", " save_initial_solution = true,\n", " save_final_solution = true)\n", "```" ], "metadata": {} }, { "cell_type": "markdown", "source": [ "Another useful callback is the `SaveRestartCallback`. It saves information for restarting\n", "in regular intervals. We are interested in saving a restart file for the final solution as\n", "well. To perform a restart, you need to configure the restart setup in a special way, which is\n", "described in the section Restart simulation." ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "save_restart = SaveRestartCallback(interval = 100, save_final_restart = true)" ], "metadata": {}, "execution_count": null }, { "cell_type": "markdown", "source": [ "Create a `CallbackSet` to collect all callbacks so that they can be passed to the `solve`\n", "function." ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "callbacks = CallbackSet(summary_callback, analysis_callback, alive_callback, stepsize_callback,\n", " save_solution, save_restart);" ], "metadata": {}, "execution_count": null }, { "cell_type": "markdown", "source": [ "The last step is to choose the time integration method. OrdinaryDiffEq.jl defines a wide range of\n", "[ODE solvers](https://docs.sciml.ai/DiffEqDocs/latest/solvers/ode_solve/), including the\n", "three-stage, third-order strong stability preserving Runge-Kutta method `SSPRK33`. We will pass\n", "the ODE problem, the ODE solver and the callbacks to the `solve` function. Also, to use\n", "`StepsizeCallback`, we must explicitly specify the initial trial time step `dt`, the selected\n", "value is not important, because it will be overwritten by the `StepsizeCallback`. And there is no\n", "need to save every step of the solution, as we are only interested the output provided by\n", "our callback `SaveSolutionCallback`." ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "sol = solve(ode, SSPRK33(); dt = 1.0, save_everystep = false, callback = callbacks);" ], "metadata": {}, "execution_count": null }, { "cell_type": "markdown", "source": [ "Finally, we print the timer summary." ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "summary_callback()" ], "metadata": {}, "execution_count": null }, { "cell_type": "markdown", "source": [ "Now you can plot the solution as shown below, analyze it and improve the stability, accuracy or\n", "efficiency of your setup." ], "metadata": {} }, { "cell_type": "markdown", "source": [ "## Visualize the solution" ], "metadata": {} }, { "cell_type": "markdown", "source": [ "In the previous part of the tutorial, we calculated the final solution of the given problem, now we want\n", "to visualize it. A more detailed explanation of visualization methods can be found in the section\n", "Visualization." ], "metadata": {} }, { "cell_type": "markdown", "source": [ "### Using Plots.jl" ], "metadata": {} }, { "cell_type": "markdown", "source": [ "The first option is to use the [Plots.jl](https://github.com/JuliaPlots/Plots.jl) package\n", "directly after calculations, when the solution is saved in the `sol` variable." ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "using Plots" ], "metadata": {}, "execution_count": null }, { "cell_type": "markdown", "source": [ "As was shown in the Getting started section, you can plot all\n", "variables from the system of equations by executing the following.\n", "```julia\n", "plot(sol)\n", "```\n", "Alternatively, you can configure the plot more precisely. Trixi.jl provides a special data type,\n", "`PlotData2D`, to extract the visualization data from the solution." ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "pd = PlotData2D(sol);" ], "metadata": {}, "execution_count": null }, { "cell_type": "markdown", "source": [ "You can plot specific variables from the system of equations by referring to their names.\n", "To obtain the names of all variables, execute the following." ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "@show pd.variable_names;" ], "metadata": {}, "execution_count": null }, { "cell_type": "markdown", "source": [ "Plot the variable named \"scalar\" (which is the name of the variable for the\n", "linear advection equation in Trixi.jl)." ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "plot(pd[\"scalar\"])" ], "metadata": {}, "execution_count": null }, { "cell_type": "markdown", "source": [ "Mesh extraction is possible using the `getmesh` function.\n", "Plots.jl has the `plot!` function that allows you to modify an already built graph." ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "plot!(getmesh(pd))" ], "metadata": {}, "execution_count": null }, { "cell_type": "markdown", "source": [ "### Using Trixi2Vtk.jl" ], "metadata": {} }, { "cell_type": "markdown", "source": [ "Another way to visualize a solution is to extract it from a saved HDF5 file. After we used the\n", "`solve` function with `SaveSolutionCallback` there is a file with the final solution.\n", "It is located in the `out` folder and is named as follows: `solution_index.h5`. The `index`\n", "is the final time step of the solution that is padded to 6 digits with zeros from the beginning.\n", "With Trixi2Vtk you can convert the HDF5 output file generated by Trixi.jl into a VTK/VTU\n", "files. VTK/VTU are specialized formats designed to store structured data required for\n", "visualization purposes. This can be used in visualization tools such as\n", "[ParaView](https://www.paraview.org) or [VisIt](https://visit.llnl.gov) to plot the solution." ], "metadata": {} }, { "cell_type": "markdown", "source": [ "If you haven't added Trixi2Vtk.jl to your project yet, you can add it as follows.\n", "```julia\n", "import Pkg\n", "Pkg.add([\"Trixi2Vtk\"])\n", "```\n", "Now we load the Trixi2Vtk.jl package and convert the file `out/solution_000032.h5` with\n", "the final solution using the `trixi2vtk` function saving the resulting file in the\n", "`out` folder." ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "using Trixi2Vtk\n", "trixi2vtk(joinpath(\"out\", \"solution_000032.h5\"), output_directory=\"out\")" ], "metadata": {}, "execution_count": null }, { "cell_type": "markdown", "source": [ "Now two files `solution_000032.vtu` and `solution_000032_celldata.vtu` have been generated in the\n", "`out` folder. The first one contains all the information for visualizing the solution, the\n", "second one contains all the cell-based or discretization-based information." ], "metadata": {} }, { "cell_type": "markdown", "source": [ "Now let's visualize the solution from the generated files in ParaView. Follow this short\n", "instruction to get the visualization.\n", "- Download, install and open [ParaView](https://www.paraview.org/download/).\n", "- Press `Ctrl+O` and select the generated files `solution_000032.vtu` and\n", " `solution_000032_celldata.vtu` from the `out` folder.\n", "- In the upper-left corner in the Pipeline Browser window, left-click on the eye-icon near\n", " `solution_000032.vtu`.\n", "- In the lower-left corner in the Properties window, change the Coloring from Solid Color to\n", " scalar. This already generates the visualization of the final solution.\n", "- Now let's add the mesh to the visualization. In the upper-left corner in the\n", " Pipeline Browser window, left-click on the eye-icon near `solution_000032_celldata.vtu`.\n", "- In the lower-left corner in the Properties window, change the Representation from Surface\n", " to Wireframe. Then a white grid should appear on the visualization.\n", "Now, if you followed the instructions exactly, you should get a similar image as shown in the\n", "section Using Plots.jl:" ], "metadata": {} }, { "cell_type": "markdown", "source": [ "![paraview_trixi2vtk_example](https://github.com/trixi-framework/Trixi.jl/assets/119304909/0c29139b-6c5d-4d5c-86e1-f4ebc95aca7e)" ], "metadata": {} }, { "cell_type": "markdown", "source": [ "After completing this tutorial you are able to set up your own simulations with\n", "Trixi.jl. If you have an interest in contributing to Trixi.jl as a developer, refer to the third\n", "part of the introduction titled Changing Trixi.jl itself." ], "metadata": {} } ], "nbformat_minor": 3, "metadata": { "language_info": { "file_extension": ".jl", "mimetype": "application/julia", "name": "julia", "version": "1.10.3" }, "kernelspec": { "name": "julia-1.10", "display_name": "Julia 1.10.3", "language": "julia" } }, "nbformat": 4 }