Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
trixi-framework
GitHub Repository: trixi-framework/Trixi.jl
Path: blob/main/docs/literate/src/files/hohqmesh_tutorial.jl
5591 views
1
#src # Unstructured meshes with HOHQMesh.jl
2
3
# Trixi.jl supports numerical approximations on unstructured quadrilateral meshes
4
# with the [`UnstructuredMesh2D`](@ref) mesh type.
5
6
# The purpose of this tutorial is to demonstrate how to use the `UnstructuredMesh2D`
7
# functionality of Trixi.jl. This begins by running and visualizing an available unstructured
8
# quadrilateral mesh example. Then, the tutorial will demonstrate how to
9
# conceptualize a problem with curved boundaries, generate
10
# a curvilinear mesh using the available software in the Trixi.jl ecosystem,
11
# and then run a simulation using Trixi.jl on said mesh.
12
13
# Unstructured quadrilateral meshes can be made
14
# with the [High-Order Hex-Quad Mesh (HOHQMesh) generator](https://github.com/trixi-framework/HOHQMesh)
15
# created and developed by David Kopriva.
16
# HOHQMesh is a mesh generator specifically designed for spectral element methods.
17
# It provides high-order boundary curve information (needed to accurately set boundary conditions)
18
# and elements can be larger (due to the high accuracy of the spatial approximation)
19
# compared to traditional finite element mesh generators.
20
# For more information about the design and features of HOHQMesh one can refer to its
21
# [official documentation](https://trixi-framework.github.io/HOHQMesh/).
22
23
# HOHQMesh is incorporated into the Trixi.jl framework via the registered Julia package
24
# [HOHQMesh.jl](https://github.com/trixi-framework/HOHQMesh.jl).
25
# This package provides a Julia wrapper for the HOHQMesh generator that allows users to easily create mesh
26
# files without the need to build HOHQMesh from source. To install the HOHQMesh package execute
27
# ````julia
28
# import Pkg; Pkg.add("HOHQMesh")
29
# ````
30
# Now we are ready to generate an unstructured quadrilateral mesh that can be used by Trixi.jl.
31
32
# ## Running and visualizing an unstructured simulation
33
34
# Trixi.jl supports solving hyperbolic problems on several mesh types.
35
# There is a default example for this mesh type that can be executed by
36
37
using Trixi
38
rm("out", force = true, recursive = true) #hide #md
39
trixi_include(default_example_unstructured())
40
41
# This will compute a smooth, manufactured solution test case for the 2D compressible Euler equations
42
# on the curved quadrilateral mesh described in the
43
# [Trixi.jl documentation](https://trixi-framework.github.io/TrixiDocumentation/stable/meshes/unstructured_quad_mesh/).
44
# For references on the method of manufactured solutions (MMS) see the following publications:
45
# - Kambiz Salari and Patrick Knupp (2000)
46
# Code Verification by the Method of Manufactured Solutions
47
# [DOI: 10.2172/759450](https://doi.org/10.2172/759450)
48
# - Patrick J. Roache (2002)
49
# Code Verification by the Method of Manufactured Solutions
50
# [DOI: 10.1115/1.1436090](https://doi.org/10.1115/1.1436090)
51
52
# Apart from the usual error and timing output provided by the Trixi.jl run, it is useful to visualize and inspect
53
# the solution. One option available in the Trixi.jl framework to visualize the solution on
54
# unstructured quadrilateral meshes is post-processing the
55
# Trixi.jl output file(s) with the [`Trixi2Vtk`](https://github.com/trixi-framework/Trixi2Vtk.jl) tool
56
# and plotting them with [ParaView](https://www.paraview.org/download/).
57
58
# To convert the HDF5-formatted `.h5` output file(s) from Trixi.jl into VTK format execute the following
59
60
using Trixi2Vtk
61
trixi2vtk("out/solution_000000180.h5", output_directory = "out")
62
63
# Note this step takes about 15-30 seconds as the package `Trixi2Vtk` must be precompiled and executed for the first time
64
# in your REPL session. The `trixi2vtk` command above will convert the solution file at the final time into a `.vtu` file
65
# which can be read in and visualized with ParaView. Optional arguments for `trixi2vtk` are: (1) Pointing to the `output_directory`
66
# where the new files will be saved; it defaults to the current directory. (2) Specifying a higher number of
67
# visualization nodes. For instance, if we want to use 12 uniformly spaced nodes for visualization we can execute
68
69
trixi2vtk("out/solution_000000180.h5", output_directory = "out", nvisnodes = 12)
70
71
# By default `trixi2vtk` sets `nvisnodes` to be the same as the number of nodes specified in
72
# the `elixir` file used to run the simulation.
73
74
# Finally, if you want to convert all the solution files to VTK execute
75
76
trixi2vtk("out/solution_000*.h5", output_directory = "out", nvisnodes = 12)
77
78
# then it is possible to open the `.pvd` file with ParaView and create a video of the simulation.
79
80
# ## Creating a mesh using HOHQMesh
81
82
# The creation of an unstructured quadrilateral mesh using HOHQMesh.jl is driven by a **control file**. In this file the user dictates
83
# the domain to be meshed, prescribes any desired boundary curvature, the polynomial order of said boundaries, etc.
84
# In this tutorial we cover several basic features of the possible control inputs. For a complete discussion
85
# on this topic see the [HOHQMesh control file documentation](https://trixi-framework.github.io/HOHQMesh/the-control-file/).
86
87
# To begin, we provide a complete control file in this tutorial. After this we give a breakdown
88
# of the control file components to explain the chosen parameters.
89
90
# Suppose we want to create a mesh of a domain with straight sided
91
# outer boundaries and a curvilinear "ice cream cone" shaped object at its center.
92
93
# ![mesh_boundary_cartoon](https://user-images.githubusercontent.com/25242486/129603954-9788500d-bba8-49be-8e6f-7555099dbf7c.png)
94
95
# The associated `ice_cream_straight_sides.control` file is created below.
96
open("out/ice_cream_straight_sides.control", "w") do io
97
println(io, raw"""
98
\begin{CONTROL_INPUT}
99
\begin{RUN_PARAMETERS}
100
mesh file name = ice_cream_straight_sides.mesh
101
plot file name = ice_cream_straight_sides.tec
102
stats file name = none
103
mesh file format = ISM-v2
104
polynomial order = 4
105
plot file format = skeleton
106
\end{RUN_PARAMETERS}
107
108
\begin{BACKGROUND_GRID}
109
x0 = [-8.0, -8.0, 0.0]
110
dx = [1.0, 1.0, 0.0]
111
N = [16,16,1]
112
\end{BACKGROUND_GRID}
113
114
\begin{SPRING_SMOOTHER}
115
smoothing = ON
116
smoothing type = LinearAndCrossBarSpring
117
number of iterations = 25
118
\end{SPRING_SMOOTHER}
119
120
\end{CONTROL_INPUT}
121
122
\begin{MODEL}
123
124
\begin{INNER_BOUNDARIES}
125
126
\begin{CHAIN}
127
name = IceCreamCone
128
\begin{END_POINTS_LINE}
129
name = LeftSlant
130
xStart = [-2.0, 1.0, 0.0]
131
xEnd = [ 0.0, -3.0, 0.0]
132
\end{END_POINTS_LINE}
133
134
\begin{END_POINTS_LINE}
135
name = RightSlant
136
xStart = [ 0.0, -3.0, 0.0]
137
xEnd = [ 2.0, 1.0, 0.0]
138
\end{END_POINTS_LINE}
139
140
\begin{CIRCULAR_ARC}
141
name = IceCream
142
units = degrees
143
center = [ 0.0, 1.0, 0.0]
144
radius = 2.0
145
start angle = 0.0
146
end angle = 180.0
147
\end{CIRCULAR_ARC}
148
\end{CHAIN}
149
150
\end{INNER_BOUNDARIES}
151
152
\end{MODEL}
153
\end{FILE}
154
""")
155
return nothing
156
end
157
158
# The first three blocks of information are wrapped within a `CONTROL_INPUT` environment block as they define the
159
# core components of the quadrilateral mesh that will be generated.
160
161
# The first block of information in `RUN_PARAMETERS` is
162
# ```
163
# \begin{RUN_PARAMETERS}
164
# mesh file name = ice_cream_straight_sides.mesh
165
# plot file name = ice_cream_straight_sides.tec
166
# stats file name = none
167
# mesh file format = ISM-v2
168
# polynomial order = 4
169
# plot file format = skeleton
170
# \end{RUN_PARAMETERS}
171
# ```
172
173
# The mesh and plot file names will be the files created by HOHQMesh once successfully executed. The stats file name is
174
# available if you wish to also save a collection of mesh statistics. For this example it is deactivated.
175
# These file names given within `RUN_PARAMETERS` **should match** that of the control file, and although this is not required by
176
# HOHQMesh, it is a useful style convention.
177
# The mesh file format `ISM-v2` in the format currently required by Trixi.jl. The `polynomial order` prescribes the order
178
# of an interpolant constructed on the Chebyshev-Gauss-Lobatto nodes that is used to represent any curved boundaries on a particular element.
179
# The plot file format of `skeleton` means that visualizing the plot file will only draw the element boundaries (and no internal nodes).
180
# Alternatively, the format can be set to `sem` to visualize the interior nodes of the approximation as well.
181
182
# The second block of information in `BACKGROUND_GRID` is
183
# ```
184
# \begin{BACKGROUND_GRID}
185
# x0 = [-8.0, -8.0, 0.0]
186
# dx = [1.0, 1.0, 0.0]
187
# N = [16,16,1]
188
# \end{BACKGROUND_GRID}
189
# ```
190
191
# This lays a grid of Cartesian elements for the domain beginning at the point `x0` as its bottom-left corner.
192
# The value of `dx`, which could differ in each direction if desired, controls the step size taken in each Cartesian direction.
193
# The values in `N` set how many Cartesian box elements are set in each coordinate direction.
194
# The above parameters define a $16\times 16$ element square mesh on $[-8,8]^2$.
195
# Further, this sets up four outer boundaries of the domain that are given the default names: `Top, Left, Bottom, Right`.
196
197
# The third block of information in `SPRING_SMOOTHER` is
198
# ```
199
# \begin{SPRING_SMOOTHER}
200
# smoothing = ON
201
# smoothing type = LinearAndCrossBarSpring
202
# number of iterations = 25
203
# \end{SPRING_SMOOTHER}
204
# ```
205
206
# Once HOHQMesh generates the mesh, a spring-mass-dashpot model is created to smooth the mesh and create "nicer" quadrilateral elements.
207
# The [default parameters of Hooke's law](https://trixi-framework.github.io/HOHQMesh/the-control-input/#the-smoother)
208
# for the spring-mass-dashpot model have been selected after a fair amount of experimentation across many meshes.
209
# If you wish to deactivate this feature you can set `smoothing = OFF` (or remove this block from the control file).
210
211
# After the `CONTROL_INPUT` environment block comes the `MODEL` environment block. It is here where the user prescribes curved boundary information with either:
212
# * An `OUTER_BOUNDARY` (covered in the next section of this tutorial).
213
# * One or more `INNER_BOUNDARIES`.
214
215
# There are several options to describe the boundary curve data to HOHQMesh like splines or parametric curves.
216
217
# For the example `ice_cream_straight_sides.control` we define three internal boundaries; two straight-sided and
218
# one as a circular arc.
219
# Within the HOHQMesh control input each curve must be assigned to a `CHAIN` as shown below in the complete
220
# `INNER_BOUNDARIES` block.
221
# ```
222
# \begin{INNER_BOUNDARIES}
223
#
224
# \begin{CHAIN}
225
# name = IceCreamCone
226
# \begin{END_POINTS_LINE}
227
# name = LeftSlant
228
# xStart = [-2.0, 1.0, 0.0]
229
# xEnd = [ 0.0, -3.0, 0.0]
230
# \end{END_POINTS_LINE}
231
#
232
# \begin{END_POINTS_LINE}
233
# name = RightSlant
234
# xStart = [ 0.0, -3.0, 0.0]
235
# xEnd = [ 2.0, 1.0, 0.0]
236
# \end{END_POINTS_LINE}
237
#
238
# \begin{CIRCULAR_ARC}
239
# name = IceCream
240
# units = degrees
241
# center = [ 0.0, 1.0, 0.0]
242
# radius = 2.0
243
# start angle = 0.0
244
# end angle = 180.0
245
# \end{CIRCULAR_ARC}
246
# \end{CHAIN}
247
#
248
# \end{INNER_BOUNDARIES}
249
# ```
250
251
# It is important to note there are two `name` quantities one for the `CHAIN` and one for the `PARAMETRIC_EQUATION_CURVE`.
252
# The name for the `CHAIN` is used internally by HOHQMesh, so if you have multiple `CHAIN`s they **must be given a unique name**.
253
# The name for the `PARAMETRIC_EQUATION_CURVE` will be printed to the appropriate boundaries within the `.mesh` file produced by
254
# HOHQMesh.
255
256
# We create the mesh file `ice_cream_straight_sides.mesh` and its associated file for plotting
257
# `ice_cream_straight_sides.tec` by using HOHQMesh.jl's function `generate_mesh`.
258
using HOHQMesh
259
control_file = joinpath("out", "ice_cream_straight_sides.control")
260
output = generate_mesh(control_file);
261
262
# The mesh file `ice_cream_straight_sides.mesh` and its associated file for plotting
263
# `ice_cream_straight_sides.tec` are placed in the `out` folder.
264
# The resulting mesh generated by HOHQMesh.jl is given in the following figure.
265
266
# ![mesh_straight_sides](https://user-images.githubusercontent.com/25242486/129603958-08e4b874-53d5-4511-9a54-6daf4c21edca.png)
267
268
# We note that Trixi.jl uses the boundary name information from the control file
269
# to assign boundary conditions in an elixir file.
270
# Therefore, the name should start with a letter and consist only of alphanumeric characters and underscores. Please note that the name will be treated as case sensitive.
271
272
# ## Example simulation on `ice_cream_straight_sides.mesh`
273
274
# With this newly generated mesh we are ready to run a Trixi.jl simulation on an unstructured quadrilateral mesh.
275
# For this we must create a new elixir file.
276
277
# The elixir file given below creates an initial condition for a
278
# uniform background flow state with a free stream Mach number of 0.3.
279
# A focus for this part of the tutorial is to specify the boundary conditions and to construct the new mesh from the
280
# file that was generated in the previous exercise.
281
282
# It is straightforward to set the different boundary
283
# condition types in an elixir by assigning a particular function to a boundary name inside a
284
# Julia dictionary, `Dict`, variable. Observe that the names of these boundaries match those provided by HOHQMesh
285
# either by default, e.g. `Bottom`, or user assigned, e.g. `IceCream`. For this problem setup use
286
# * Freestream boundary conditions on the four box edges.
287
# * Free slip wall boundary condition on the interior curved boundaries.
288
289
# To construct the unstructured quadrilateral mesh from the HOHQMesh file we point to the appropriate location
290
# with the variable `mesh_file` and then feed this into the constructor for the [`UnstructuredMesh2D`](@ref) type in Trixi.jl
291
292
# ````julia
293
# # create the unstructured mesh from your mesh file
294
# using Trixi
295
# mesh_file = joinpath("out", "ice_cream_straight_sides.mesh")
296
# mesh = UnstructuredMesh2D(mesh_file);
297
# ````
298
299
# The complete elixir file for this simulation example is given below.
300
301
using OrdinaryDiffEqLowStorageRK
302
using Trixi
303
304
equations = CompressibleEulerEquations2D(1.4) # set gas gamma = 1.4
305
306
## freestream flow state with Ma_inf = 0.3
307
@inline function uniform_flow_state(x, t, equations::CompressibleEulerEquations2D)
308
309
## set the freestream flow parameters
310
rho_freestream = 1.0
311
u_freestream = 0.3
312
p_freestream = inv(equations.gamma)
313
314
theta = 0.0 # zero angle of attack
315
si, co = sincos(theta)
316
v1 = u_freestream * co
317
v2 = u_freestream * si
318
319
prim = SVector(rho_freestream, v1, v2, p_freestream)
320
return prim2cons(prim, equations)
321
end
322
323
## initial condition
324
initial_condition = uniform_flow_state
325
326
## boundary condition types
327
boundary_condition_uniform_flow = BoundaryConditionDirichlet(uniform_flow_state)
328
329
## boundary conditions (NamedTuple)
330
boundary_conditions = (; Bottom = boundary_condition_uniform_flow,
331
Top = boundary_condition_uniform_flow,
332
Right = boundary_condition_uniform_flow,
333
Left = boundary_condition_uniform_flow,
334
LeftSlant = boundary_condition_slip_wall,
335
RightSlant = boundary_condition_slip_wall,
336
IceCream = boundary_condition_slip_wall);
337
338
## DGSEM solver.
339
## 1) polydeg must be >= the polynomial order set in the HOHQMesh control file to guarantee
340
## freestream preservation. As a extra task try setting polydeg=3
341
## 2) VolumeIntegralFluxDifferencing with central volume flux is activated
342
## for dealiasing
343
volume_flux = flux_ranocha
344
solver = DGSEM(polydeg = 4, surface_flux = flux_hll,
345
volume_integral = VolumeIntegralFluxDifferencing(volume_flux))
346
347
## create the unstructured mesh from your mesh file
348
mesh_file = joinpath("out", "ice_cream_straight_sides.mesh")
349
mesh = UnstructuredMesh2D(mesh_file)
350
351
## Create semidiscretization with all spatial discretization-related components
352
semi = SemidiscretizationHyperbolic(mesh, equations, initial_condition, solver,
353
boundary_conditions = boundary_conditions)
354
355
## Create ODE problem from semidiscretization with time span from 0.0 to 2.0
356
tspan = (0.0, 2.0)
357
ode = semidiscretize(semi, tspan)
358
359
## Create the callbacks to output solution files and adapt the time step
360
summary_callback = SummaryCallback()
361
save_solution = SaveSolutionCallback(interval = 10,
362
save_initial_solution = true,
363
save_final_solution = true)
364
stepsize_callback = StepsizeCallback(cfl = 1.0)
365
366
callbacks = CallbackSet(summary_callback, save_solution, stepsize_callback)
367
368
## Evolve ODE problem in time using `solve` from OrdinaryDiffEq
369
sol = solve(ode, CarpenterKennedy2N54(williamson_condition = false);
370
dt = 1.0, # solve needs some value here but it will be overwritten by the stepsize_callback
371
ode_default_options()..., callback = callbacks)
372
373
# Visualization of the solution is carried out in a similar way as above. That is, one converts the `.h5`
374
# output files with `trixi2vtk` and then plot the solution in ParaView. An example plot of the pressure
375
# at the final time is shown below.
376
377
# ![simulation_straight_sides](https://user-images.githubusercontent.com/25242486/129733926-6ef80676-779b-4f1e-9826-3ebf750cf382.png)
378
379
# ## Making a mesh with a curved outer boundary
380
381
# Let us modify the mesh from the previous task and place a circular outer boundary instead
382
# of straight-sided outer boundaries.
383
# Note, the "ice cream cone" shape is still placed at the center of the domain.
384
385
# We create the new control file `ice_cream_curved_sides.control` file below and will then highlight the
386
# major differences compared to `ice_cream_straight_sides.control`.
387
open("out/ice_cream_curved_sides.control", "w") do io
388
println(io, raw"""
389
\begin{CONTROL_INPUT}
390
\begin{RUN_PARAMETERS}
391
mesh file name = ice_cream_curved_sides.mesh
392
plot file name = ice_cream_curved_sides.tec
393
stats file name = none
394
mesh file format = ISM-v2
395
polynomial order = 4
396
plot file format = skeleton
397
\end{RUN_PARAMETERS}
398
399
\begin{BACKGROUND_GRID}
400
background grid size = [1.0, 1.0, 0.0]
401
\end{BACKGROUND_GRID}
402
403
\begin{SPRING_SMOOTHER}
404
smoothing = ON
405
smoothing type = LinearAndCrossBarSpring
406
number of iterations = 25
407
\end{SPRING_SMOOTHER}
408
409
\end{CONTROL_INPUT}
410
411
\begin{MODEL}
412
413
\begin{OUTER_BOUNDARY}
414
\begin{PARAMETRIC_EQUATION_CURVE}
415
name = OuterCircle
416
xEqn = x(t) = 8.0*sin(2.0*pi*t)
417
yEqn = y(t) = 8.0*cos(2.0*pi*t)
418
zEqn = z(t) = 0.0
419
\end{PARAMETRIC_EQUATION_CURVE}
420
421
\end{OUTER_BOUNDARY}
422
423
\begin{INNER_BOUNDARIES}
424
425
\begin{CHAIN}
426
name = IceCreamCone
427
\begin{END_POINTS_LINE}
428
name = LeftSlant
429
xStart = [-2.0, 1.0, 0.0]
430
xEnd = [ 0.0, -3.0, 0.0]
431
\end{END_POINTS_LINE}
432
433
\begin{END_POINTS_LINE}
434
name = RightSlant
435
xStart = [ 0.0, -3.0, 0.0]
436
xEnd = [ 2.0, 1.0, 0.0]
437
\end{END_POINTS_LINE}
438
439
\begin{CIRCULAR_ARC}
440
name = IceCream
441
units = degrees
442
center = [ 0.0, 1.0, 0.0]
443
radius = 2.0
444
start angle = 0.0
445
end angle = 180.0
446
\end{CIRCULAR_ARC}
447
\end{CHAIN}
448
449
\end{INNER_BOUNDARIES}
450
451
\end{MODEL}
452
\end{FILE}
453
""")
454
return nothing
455
end
456
457
# The first alteration is that we have altered the second block of information
458
# `BACKGROUND_GRID` within the `CONTROL_INPUT` to be
459
# ```
460
# \begin{BACKGROUND_GRID}
461
# background grid size = [1.0, 1.0, 0.0]
462
# \end{BACKGROUND_GRID}
463
# ```
464
# This mesh control file has an outer boundary that determines the extent of the domain to be meshed.
465
# Therefore, we only need to supply the `background grid size` to the `BACKGROUND_GRID` control input.
466
467
# The second alteration is that the `MODEL` now contains information for an `OUTER_BOUNDARY`.
468
# In this case it is a circle of radius `8` centered at `[0.0, 0.0, 0.0]` written as a set of
469
# `PARAMETRIC_EQUATION_CURVE`s.
470
# ```
471
# \begin{OUTER_BOUNDARY}
472
#
473
# \begin{PARAMETRIC_EQUATION_CURVE}
474
# name = OuterCircle
475
# xEqn = x(t) = 8.0*sin(2.0*pi*t)
476
# yEqn = y(t) = 8.0*cos(2.0*pi*t)
477
# zEqn = z(t) = 0.0
478
# \end{PARAMETRIC_EQUATION_CURVE}
479
#
480
# \end{OUTER_BOUNDARY}
481
# ```
482
# Just as with the inner boundary curves, we must assign a name to the `OUTER_BOUNDARY`. It will be included
483
# in the generated `.mesh` file and is used within the Trixi.jl elixir file to set boundary conditions.
484
485
# Again, we create the `.mesh` and `.tec` files with HOHQMesh.jl's function `generate_mesh`
486
control_file = joinpath("out", "ice_cream_curved_sides.control")
487
output = generate_mesh(control_file);
488
# The files are placed in the `out` folder.
489
490
# The resulting mesh generated by HOHQMesh.jl is given in the following figure.
491
492
# ![mesh_curved_sides](https://user-images.githubusercontent.com/25242486/129603957-6a92618f-9ed8-4072-b6ab-05533bea746a.png)
493
494
# ## Running Trixi.jl on `ice_cream_curved_sides.mesh`
495
496
# We can reuse much of the elixir file to setup the uniform flow over an ice cream cone from the
497
# previous part of this tutorial. The only component of the elixir file that must be changed is the boundary condition
498
# `NamedTuple` because we now have a boundary named `OuterCircle` instead of four edges of a bounding box.
499
500
## boundary conditions (NamedTuple)
501
boundary_conditions = (; OuterCircle = boundary_condition_uniform_flow,
502
LeftSlant = boundary_condition_slip_wall,
503
RightSlant = boundary_condition_slip_wall,
504
IceCream = boundary_condition_slip_wall);
505
506
# Also, we must update the construction of the mesh from our new mesh file `ice_cream_curved_sides.mesh` that
507
# is located in the `out` folder.
508
509
## create the unstructured mesh from your mesh file
510
mesh_file = joinpath("out", "ice_cream_curved_sides.mesh")
511
mesh = UnstructuredMesh2D(mesh_file);
512
513
# We can then post-process the solution file at the final time on the new mesh with `Trixi2Vtk` and visualize with ParaView.
514
515
# ![simulation_curved_sides](https://user-images.githubusercontent.com/25242486/129733924-778795c1-9119-419a-8b89-bcbe13e33cd7.png)
516
517
# ## Setting up a simulation with AMR via `P4estMesh`
518
# The above explained mesh file format of `ISM-V2` only works with `UnstructuredMesh2D` and so does
519
# not support AMR. On the other hand, the mesh type [`P4estMesh`](@ref) allows AMR. The mesh
520
# constructor for the `P4estMesh` imports an unstructured, conforming mesh from an Abaqus mesh file
521
# (`.inp`).
522
523
# As described above, the first block of the HOHQMesh control file contains the parameter
524
# `mesh file format`. If you set `mesh file format = ABAQUS` instead of `ISM-V2`,
525
# HOHQMesh.jl's function `generate_mesh` creates an Abaqus mesh file `.inp`.
526
# ````julia
527
# using HOHQMesh
528
# control_file = joinpath("out", "ice_cream_straight_sides.control")
529
# output = generate_mesh(control_file);
530
# ````
531
532
# Now, you can create a `P4estMesh` from your mesh file. It is described in detail in the
533
# [P4est-based mesh](https://trixi-framework.github.io/TrixiDocumentation/stable/meshes/p4est_mesh/#P4est-based-mesh)
534
# part of the Trixi.jl docs.
535
# ````julia
536
# using Trixi
537
# mesh_file = joinpath("out", "ice_cream_straight_sides.inp")
538
# mesh = P4estMesh{2}(mesh_file)
539
# ````
540
541
# Since `P4estMesh` supports AMR, we just have to extend the setup from the first example by the
542
# standard AMR procedure. For more information about AMR in Trixi.jl, see the [matching tutorial](@ref adaptive_mesh_refinement).
543
544
# ````julia
545
# amr_indicator = IndicatorLöhner(semi, variable=density)
546
547
# amr_controller = ControllerThreeLevel(semi, amr_indicator,
548
# base_level=0,
549
# med_level =1, med_threshold=0.05,
550
# max_level =3, max_threshold=0.1)
551
552
# amr_callback = AMRCallback(semi, amr_controller,
553
# interval=5,
554
# adapt_initial_condition=true,
555
# adapt_initial_condition_only_refine=true)
556
557
# callbacks = CallbackSet(..., amr_callback)
558
# ````
559
560
# We can then post-process the solution file at the final time on the new mesh with `Trixi2Vtk` and visualize
561
# with ParaView, see the appropriate [visualization section](https://trixi-framework.github.io/TrixiDocumentation/stable/visualization/#Trixi2Vtk)
562
# for details.
563
564
# ![simulation_straight_sides_p4est_amr](https://user-images.githubusercontent.com/74359358/168049930-8abce6ac-cd47-4d04-b40b-0fa459bbd98d.png)
565
566
# ## Package versions
567
568
# These results were obtained using the following versions.
569
570
using InteractiveUtils
571
versioninfo()
572
573
using Pkg
574
Pkg.status(["Trixi", "OrdinaryDiffEqLowStorageRK", "Plots", "Trixi2Vtk", "HOHQMesh"],
575
mode = PKGMODE_MANIFEST)
576
577