Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
trixi-framework
GitHub Repository: trixi-framework/Trixi.jl
Path: blob/main/docs/literate/src/files/first_steps/getting_started.jl
5593 views
1
#src # Getting started
2
3
# Trixi.jl is a numerical simulation framework for conservation laws and
4
# is written in the [Julia programming language](https://julialang.org/).
5
# This tutorial is intended for beginners in Julia and Trixi.jl.
6
# After reading it, you will know how to install Julia and Trixi.jl on your computer,
7
# and you will be able to download setup files from our GitHub repository, modify them,
8
# and run simulations.
9
10
# The contents of this tutorial:
11
# - [Julia installation](@ref Julia-installation)
12
# - [Trixi.jl installation](@ref Trixi.jl-installation)
13
# - [Running a simulation](@ref Running-a-simulation)
14
# - [Getting an existing setup file](@ref Getting-an-existing-setup-file)
15
# - [Modifying an existing setup](@ref Modifying-an-existing-setup)
16
17
# ## Julia installation
18
19
# Trixi.jl is compatible with the latest stable release of Julia. Additional details regarding Julia
20
# support can be found in the [`README.md`](https://github.com/trixi-framework/Trixi.jl#installation)
21
# file. After installation, the current default Julia version can be managed through the command
22
# line tool `juliaup`. You may follow our concise installation guidelines for Windows, Linux, and
23
# MacOS provided below. In the event of any issues during the installation process, please consult
24
# the official [Julia installation instruction](https://julialang.org/downloads/).
25
26
# ### Windows
27
28
# - Open a terminal by pressing `Win+r` and entering `cmd` in the opened window.
29
# - To install Julia, execute the following command in the terminal:
30
# ```shell
31
# winget install julia -s msstore
32
# ```
33
# Note: For this installation an MS Store account is necessary to proceed.
34
# - Verify the successful installation of Julia by executing the following command in the terminal:
35
# ```shell
36
# julia
37
# ```
38
# To exit Julia, execute `exit()` or press `Ctrl+d`.
39
40
# ### Linux and MacOS
41
42
# - To install Julia, run the following command in a terminal:
43
# ```shell
44
# curl -fsSL https://install.julialang.org | sh
45
# ```
46
# Follow the instructions displayed in the terminal during the installation process.
47
# - If an error occurs during the execution of the previous command, you may need to install
48
# `curl`. On Ubuntu-type systems, you can use the following command:
49
# ```shell
50
# sudo apt install curl
51
# ```
52
# After installing `curl`, repeat the first step once more to proceed with Julia installation.
53
# - Verify the successful installation of Julia by executing the following command in the terminal:
54
# ```shell
55
# julia
56
# ```
57
# To exit Julia, execute `exit()` or press `Ctrl+d`.
58
59
# ## Trixi.jl installation
60
61
# Trixi.jl and its related tools are registered Julia packages, thus their installation
62
# happens inside Julia.
63
# For a smooth workflow experience with Trixi.jl, you need to install
64
# [Trixi.jl](https://github.com/trixi-framework/Trixi.jl),
65
# [Plots.jl](https://github.com/JuliaPlots/Plots.jl), and sub-packages of
66
# [OrdinaryDiffEq.jl](https://github.com/SciML/OrdinaryDiffEq.jl).
67
68
# - Open a terminal and start Julia.
69
# - Execute the following commands to install all mentioned packages. Please note that the
70
# installation process involves downloading and precompiling the source code, which may take
71
# some time depending on your machine.
72
# ````julia
73
# import Pkg
74
# Pkg.add(["OrdinaryDiffEqLowStorageRK", "OrdinaryDiffEqSSPRK", "Plots", "Trixi"])
75
# ````
76
# - On Windows, the firewall may request permission to install packages.
77
78
# Besides Trixi.jl you have now installed additional packages:
79
# sub-packages of [OrdinaryDiffEq.jl](https://github.com/SciML/OrdinaryDiffEq.jl) provide time
80
# integration schemes used by Trixi.jl and [Plots.jl](https://github.com/JuliaPlots/Plots.jl)
81
# can be used to directly visualize Trixi.jl results from the Julia REPL.
82
83
# ## Usage
84
85
# ### Running a simulation
86
87
# To get you started, Trixi.jl has a large set
88
# of [example setups](https://github.com/trixi-framework/Trixi.jl/tree/main/examples), that can be
89
# taken as a basis for your future investigations. In Trixi.jl, we call these setup files
90
# "elixirs", since they contain Julia code that takes parts of Trixi.jl and combines them into
91
# something new.
92
93
# Any of the examples can be executed using the [`trixi_include`](@ref)
94
# function. `trixi_include(...)` expects
95
# a single string argument with a path to a file containing Julia code.
96
# For convenience, the [`examples_dir`](@ref) function returns a path to the
97
# [`examples`](https://github.com/trixi-framework/Trixi.jl/tree/main/examples)
98
# folder, which has been locally downloaded while installing Trixi.jl.
99
# `joinpath(...)` can be used to join path components into a full path.
100
101
# Let's execute a short two-dimensional problem setup. It approximates the solution of
102
# the compressible Euler equations in 2D for an ideal gas ([`CompressibleEulerEquations2D`](@ref))
103
# with a weak blast wave as the initial condition and periodic boundary conditions.
104
105
# The compressible Euler equations in two spatial dimensions are given by
106
# ```math
107
# \frac{\partial}{\partial t}
108
# \begin{pmatrix}
109
# \rho \\ \rho v_1 \\ \rho v_2 \\ \rho e_{\text{total}}
110
# \end{pmatrix}
111
# +
112
# \frac{\partial}{\partial x}
113
# \begin{pmatrix}
114
# \rho v_1 \\ \rho v_1^2 + p \\ \rho v_1 v_2 \\ (\rho e_{\text{total}} + p) v_1
115
# \end{pmatrix}
116
# +
117
# \frac{\partial}{\partial y}
118
# \begin{pmatrix}
119
# \rho v_2 \\ \rho v_1 v_2 \\ \rho v_2^2 + p \\ (\rho e_{\text{total}} + p) v_2
120
# \end{pmatrix}
121
# =
122
# \begin{pmatrix}
123
# 0 \\ 0 \\ 0 \\ 0
124
# \end{pmatrix},
125
# ```
126
# for an ideal gas with the specific heat ratio ``\gamma``.
127
# Here, ``\rho`` is the density, ``v_1`` and ``v_2`` are the velocities, ``e_{\text{total}}`` is the specific
128
# total energy, and
129
# ```math
130
# p = (\gamma - 1) \left( \rho e_{\text{total}} - \frac{1}{2} \rho (v_1^2 + v_2^2) \right)
131
# ```
132
# is the pressure.
133
134
# The [`initial_condition_weak_blast_wave`](@ref) is specified in
135
# [`compressible_euler_2d.jl`](https://github.com/trixi-framework/Trixi.jl/blob/main/src/equations/compressible_euler_2d.jl)
136
137
# Start Julia in a terminal and execute the following code:
138
139
# ````julia
140
# using Trixi, OrdinaryDiffEq
141
# trixi_include(joinpath(examples_dir(), "tree_2d_dgsem", "elixir_euler_ec.jl"))
142
# ````
143
using Trixi, OrdinaryDiffEqLowStorageRK #hide #md
144
trixi_include(@__MODULE__, joinpath(examples_dir(), "tree_2d_dgsem", "elixir_euler_ec.jl")); #hide #md
145
146
# The output contains a recap of the setup and various information about the course of the simulation.
147
# For instance, the solution was approximated over the [`TreeMesh`](@ref) with 1024 effective cells using
148
# the `CarpenterKennedy2N54` ODE
149
# solver. Further details about the ODE solver can be found in the
150
# [documentation of OrdinaryDiffEq.jl](https://docs.sciml.ai/DiffEqDocs/stable/solvers/ode_solve/#Low-Storage-Methods)
151
152
# To analyze the result of the computation, we can use the Plots.jl package and the function
153
# `plot(...)`, which creates a graphical representation of the solution. `sol` is a variable
154
# defined in the executed example and it contains the solution after the simulation
155
# finishes. `sol.u` holds the vector of values at each saved timestep, while `sol.t` holds the
156
# corresponding times for each saved timestep. In this instance, only two timesteps were saved: the
157
# initial and final ones. The plot depicts the distribution of the weak blast wave at the final moment
158
# of time, showing the density, velocities, and pressure of the ideal gas across a 2D domain.
159
160
using Plots
161
plot(sol)
162
163
# ### Getting an existing setup file
164
165
# To obtain a list of all Trixi.jl elixirs execute
166
# [`get_examples`](@ref). It returns the paths to all example setups.
167
168
get_examples()
169
170
# Editing an existing elixir is the best way to start your first own investigation using Trixi.jl.
171
172
# To edit an existing elixir, you first have to find a suitable one and then copy it to a local
173
# folder. Let's have a look at how to download the `elixir_euler_ec.jl` elixir used in the previous
174
# section from the [Trixi.jl GitHub repository](https://github.com/trixi-framework/Trixi.jl).
175
176
# - All examples are located inside
177
# the [`examples`](https://github.com/trixi-framework/Trixi.jl/tree/main/examples) folder.
178
# - Navigate to the
179
# file [`elixir_euler_ec.jl`](https://github.com/trixi-framework/Trixi.jl/blob/main/examples/tree_2d_dgsem/elixir_euler_ec.jl).
180
# - Right-click the `Raw` button on the right side of the webpage and choose `Save as...`
181
# (or `Save Link As...`).
182
# - Choose a folder and save the file.
183
184
# ### Modifying an existing setup
185
186
# As an example, we will change the initial condition for calculations that occur in
187
# `elixir_euler_ec.jl`. Initial conditions for [`CompressibleEulerEquations2D`](@ref) consist of
188
# initial values for ``\rho``, ``\rho v_1``, ``\rho v_2`` and ``\rho e_{\text{total}}``. One of the common initial
189
# conditions for the compressible Euler equations is a simple density wave. Let's implement it.
190
191
# - Open the downloaded file `elixir_euler_ec.jl` with a text editor.
192
# - Go to the line with the following code:
193
# ````julia
194
# initial_condition = initial_condition_weak_blast_wave
195
# ````
196
# Here, [`initial_condition_weak_blast_wave`](@ref) is used as the initial condition.
197
# - Comment out the line using the `#` symbol:
198
# ````julia
199
# # initial_condition = initial_condition_weak_blast_wave
200
# ````
201
# - Now you can create your own initial conditions. Add the following code after the
202
# commented line:
203
204
function initial_condition_density_waves(x, t, equations::CompressibleEulerEquations2D)
205
v1 = 0.1 # velocity along x-axis
206
v2 = 0.2 # velocity along y-axis
207
rho = 1.0 + 0.98 * sinpi(sum(x) - t * (v1 + v2)) # density wave profile
208
p = 20 # pressure
209
rho_e_total = p / (equations.gamma - 1) + 1 / 2 * rho * (v1^2 + v2^2)
210
return SVector(rho, rho * v1, rho * v2, rho_e_total)
211
end
212
initial_condition = initial_condition_density_waves
213
nothing; #hide #md
214
215
# - Execute the following code one more time, but instead of `path/to/file` paste the path to the
216
# `elixir_euler_ec.jl` file that you just edited.
217
# ````julia
218
# using Trixi
219
# trixi_include(path/to/file);
220
# using Plots
221
# plot(sol)
222
# ````
223
# Then you will obtain a new solution from running the simulation with a different initial
224
# condition.
225
226
trixi_include(@__MODULE__, joinpath(examples_dir(), "tree_2d_dgsem", "elixir_euler_ec.jl"), #hide #md
227
initial_condition = initial_condition); #hide #md
228
pd = PlotData2D(sol) #hide #md
229
p1 = plot(pd["rho"]) #hide #md
230
p2 = plot(pd["v1"], clim = (0.05, 0.15)) #hide #md
231
p3 = plot(pd["v2"], clim = (0.15, 0.25)) #hide #md
232
p4 = plot(pd["p"], clim = (10, 30)) #hide #md
233
plot(p1, p2, p3, p4) #hide #md
234
235
# To get exactly the same picture execute the following.
236
# ````julia
237
# pd = PlotData2D(sol)
238
# p1 = plot(pd["rho"])
239
# p2 = plot(pd["v1"], clim=(0.05, 0.15))
240
# p3 = plot(pd["v2"], clim=(0.15, 0.25))
241
# p4 = plot(pd["p"], clim=(10, 30))
242
# plot(p1, p2, p3, p4)
243
# ````
244
245
# Feel free to make further changes to the initial condition to observe different solutions.
246
247
# Now you are able to download, modify and execute simulation setups for Trixi.jl. To explore
248
# further details on setting up a new simulation with Trixi.jl, refer to the second part of
249
# the introduction titled [Create your first setup](@ref create_first_setup).
250
251
Sys.rm("out"; recursive = true, force = true) #hide #md
252
253