Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
trixi-framework
GitHub Repository: trixi-framework/Trixi.jl
Path: blob/main/src/meshes/p4est_mesh.jl
5586 views
1
# By default, Julia/LLVM does not use fused multiply-add operations (FMAs).
2
# Since these FMAs can increase the performance of many numerical algorithms,
3
# we need to opt-in explicitly.
4
# See https://ranocha.de/blog/Optimizing_EC_Trixi for further details.
5
@muladd begin
6
#! format: noindent
7
8
"""
9
P4estMesh{NDIMS, NDIMS_AMBIENT} <: AbstractMesh{NDIMS}
10
11
An unstructured curved mesh based on trees that uses the C library `p4est`
12
to manage trees and mesh refinement.
13
14
The parameter `NDIMS` denotes the dimension of the spatial domain or manifold represented
15
by the mesh itself, while `NDIMS_AMBIENT` denotes the dimension of the ambient space in
16
which the mesh is embedded. For example, the type `P4estMesh{3, 3}` corresponds to a
17
standard mesh for a three-dimensional volume, whereas `P4estMesh{2, 3}` corresponds to a
18
mesh for a two-dimensional surface or shell in three-dimensional space.
19
20
!!! warning "Experimental implementation"
21
The use of `NDIMS != NDIMS_AMBIENT` is an experimental feature and may change in future
22
releases.
23
"""
24
mutable struct P4estMesh{NDIMS, NDIMS_AMBIENT, RealT <: Real, IsParallel, P, Ghost,
25
NDIMSP2, NNODES} <:
26
AbstractMesh{NDIMS}
27
const p4est::P # Either PointerWrapper{p4est_t} or PointerWrapper{p8est_t}
28
const is_parallel::IsParallel
29
# Either `PointerWrapper{p4est_ghost_t}` or `PointerWrapper{p8est_ghost_t}`.
30
# Required for ghost/halo layers in parallel runs, thus mutable.
31
ghost::Ghost
32
# Coordinates at the nodes specified by the tensor product of `nodes` (NDIMS times).
33
# This specifies the geometry interpolation for each tree.
34
const tree_node_coordinates::Array{RealT, NDIMSP2} # [dimension, i, j, k, tree]
35
const nodes::SVector{NNODES, RealT}
36
const boundary_names::Array{Symbol, 2} # [face direction, tree]
37
current_filename::String
38
unsaved_changes::Bool
39
const p4est_partition_allow_for_coarsening::Bool
40
41
function P4estMesh{NDIMS}(p4est, tree_node_coordinates, nodes, boundary_names,
42
current_filename, unsaved_changes,
43
p4est_partition_allow_for_coarsening) where {NDIMS}
44
if NDIMS == 2
45
@assert p4est isa Ptr{p4est_t}
46
elseif NDIMS == 3
47
@assert p4est isa Ptr{p8est_t}
48
end
49
50
if mpi_isparallel()
51
if !P4est.uses_mpi()
52
error("p4est library does not support MPI")
53
end
54
is_parallel = True()
55
else
56
is_parallel = False()
57
end
58
59
p4est_pw = PointerWrapper(p4est)
60
61
ghost = ghost_new_p4est(p4est)
62
ghost_pw = PointerWrapper(ghost)
63
64
# To enable the treatment of a manifold of dimension NDIMS embedded within an
65
# ambient space of dimension NDIMS_AMBIENT, we store both as type parameters and
66
# allow them to differ in the general case. This functionality is used for
67
# constructing discretizations on spherical shell domains for applications in
68
# global atmospheric modelling. The ambient dimension NDIMS_AMBIENT is therefore
69
# set here in the inner constructor to size(tree_node_coordinates, 1).
70
mesh = new{NDIMS, size(tree_node_coordinates, 1),
71
eltype(tree_node_coordinates), typeof(is_parallel),
72
typeof(p4est_pw), typeof(ghost_pw), NDIMS + 2, length(nodes)}(p4est_pw,
73
is_parallel,
74
ghost_pw,
75
tree_node_coordinates,
76
nodes,
77
boundary_names,
78
current_filename,
79
unsaved_changes,
80
p4est_partition_allow_for_coarsening)
81
82
# Destroy `p4est` structs when the mesh is garbage collected
83
finalizer(destroy_mesh, mesh)
84
85
return mesh
86
end
87
end
88
89
const P4estMeshSerial{NDIMS} = P4estMesh{NDIMS, <:Any, <:Real, <:False}
90
const P4estMeshParallel{NDIMS} = P4estMesh{NDIMS, <:Any, <:Real, <:True}
91
92
@inline mpi_parallel(mesh::P4estMeshSerial) = False()
93
@inline mpi_parallel(mesh::P4estMeshParallel) = True()
94
95
function destroy_mesh(mesh::P4estMesh{2})
96
connectivity = mesh.p4est.connectivity
97
p4est_ghost_destroy(mesh.ghost)
98
p4est_destroy(mesh.p4est)
99
return p4est_connectivity_destroy(connectivity)
100
end
101
102
function destroy_mesh(mesh::P4estMesh{3})
103
connectivity = mesh.p4est.connectivity
104
p8est_ghost_destroy(mesh.ghost)
105
p8est_destroy(mesh.p4est)
106
return p8est_connectivity_destroy(connectivity)
107
end
108
109
@inline Base.ndims(::P4estMesh{NDIMS}) where {NDIMS} = NDIMS
110
@inline Base.real(::P4estMesh{NDIMS, NDIMS_AMBIENT, RealT}) where {NDIMS, NDIMS_AMBIENT, RealT} = RealT
111
@inline ndims_ambient(::P4estMesh{NDIMS, NDIMS_AMBIENT}) where {NDIMS, NDIMS_AMBIENT} = NDIMS_AMBIENT
112
113
@inline function ntrees(mesh::P4estMesh)
114
return mesh.p4est.trees.elem_count[]
115
end
116
# returns Int32 by default which causes a weird method error when creating the cache
117
@inline ncells(mesh::P4estMesh) = Int(mesh.p4est.local_num_quadrants[])
118
@inline ncellsglobal(mesh::P4estMesh) = Int(mesh.p4est.global_num_quadrants[])
119
120
function Base.show(io::IO, mesh::P4estMesh)
121
print(io, "P4estMesh{", ndims(mesh), ", ", ndims_ambient(mesh), ", ", real(mesh),
122
"}")
123
return nothing
124
end
125
126
function Base.show(io::IO, ::MIME"text/plain", mesh::P4estMesh)
127
if get(io, :compact, false)
128
show(io, mesh)
129
else
130
setup = [
131
"#trees" => ntrees(mesh),
132
"current #cells" => ncellsglobal(mesh),
133
"polydeg" => length(mesh.nodes) - 1
134
]
135
summary_box(io,
136
"P4estMesh{" * string(ndims(mesh)) * ", " *
137
string(ndims_ambient(mesh)) *
138
", " * string(real(mesh)) * "}", setup)
139
end
140
end
141
142
"""
143
P4estMesh(trees_per_dimension; polydeg,
144
mapping=nothing, faces=nothing, coordinates_min=nothing, coordinates_max=nothing,
145
RealT=Float64, initial_refinement_level=0, periodicity=false, unsaved_changes=true,
146
p4est_partition_allow_for_coarsening=true)
147
148
Create a structured curved/higher-order `P4estMesh` of the specified size.
149
150
There are three ways to map the mesh to the physical domain.
151
1. Define a `mapping` that maps the hypercube `[-1, 1]^n`.
152
2. Specify a `Tuple` `faces` of functions that parametrize each face.
153
3. Create a rectangular mesh by specifying `coordinates_min` and `coordinates_max`.
154
155
Non-periodic boundaries will be called `:x_neg`, `:x_pos`, `:y_neg`, `:y_pos`, `:z_neg`, `:z_pos`.
156
157
# Arguments
158
- `trees_per_dimension::NTupleE{NDIMS, Int}`: the number of trees in each dimension.
159
- `polydeg::Integer`: polynomial degree used to store the geometry of the mesh.
160
The mapping will be approximated by an interpolation polynomial
161
of the specified degree for each tree.
162
- `mapping`: a function of `NDIMS` variables to describe the mapping that transforms
163
the reference mesh (`[-1, 1]^n`) to the physical domain.
164
Use only one of `mapping`, `faces` and `coordinates_min`/`coordinates_max`.
165
- `faces::NTuple{2*NDIMS}`: a tuple of `2 * NDIMS` functions that describe the faces of the domain.
166
Each function must take `NDIMS-1` arguments.
167
`faces[1]` describes the face onto which the face in negative x-direction
168
of the unit hypercube is mapped. The face in positive x-direction of
169
the unit hypercube will be mapped onto the face described by `faces[2]`.
170
`faces[3:4]` describe the faces in positive and negative y-direction respectively
171
(in 2D and 3D).
172
`faces[5:6]` describe the faces in positive and negative z-direction respectively (in 3D).
173
Use only one of `mapping`, `faces` and `coordinates_min`/`coordinates_max`.
174
- `coordinates_min`: vector or tuple of the coordinates of the corner in the negative direction of each dimension
175
to create a rectangular mesh.
176
Use only one of `mapping`, `faces` and `coordinates_min`/`coordinates_max`.
177
- `coordinates_max`: vector or tuple of the coordinates of the corner in the positive direction of each dimension
178
to create a rectangular mesh.
179
Use only one of `mapping`, `faces` and `coordinates_min`/`coordinates_max`.
180
- `RealT::Type`: the type that should be used for coordinates.
181
- `initial_refinement_level::Integer`: refine the mesh uniformly to this level before the simulation starts.
182
- `periodicity`: either a `Bool` deciding if all of the boundaries are periodic or an `NTuple{NDIMS, Bool}`
183
deciding for each dimension if the boundaries in this dimension are periodic.
184
- `unsaved_changes::Bool`: if set to `true`, the mesh will be saved to a mesh file.
185
- `p4est_partition_allow_for_coarsening::Bool`: Must be `true` when using AMR to make mesh adaptivity
186
independent of domain partitioning. Should be `false` for static meshes
187
to permit more fine-grained partitioning.
188
"""
189
function P4estMesh(trees_per_dimension; polydeg,
190
mapping = nothing, faces = nothing, coordinates_min = nothing,
191
coordinates_max = nothing,
192
RealT = Float64, initial_refinement_level = 0, periodicity = false,
193
unsaved_changes = true,
194
p4est_partition_allow_for_coarsening = true)
195
@assert ((coordinates_min === nothing)===(coordinates_max === nothing)) "Either both or none of coordinates_min and coordinates_max must be specified"
196
197
coordinates_min_max_check(coordinates_min, coordinates_max)
198
199
@assert count(i -> i !== nothing,
200
(mapping, faces, coordinates_min))==1 "Exactly one of mapping, faces and coordinates_min/max must be specified"
201
202
# Extract mapping
203
if faces !== nothing
204
validate_faces(faces)
205
mapping = transfinite_mapping(faces)
206
elseif coordinates_min !== nothing
207
mapping = coordinates2mapping(coordinates_min, coordinates_max)
208
end
209
210
NDIMS = length(trees_per_dimension)
211
212
# Convert periodicity to a Tuple of a Bool for every dimension
213
if all(periodicity)
214
# Also catches case where periodicity = true
215
periodicity = ntuple(_ -> true, NDIMS)
216
elseif !any(periodicity)
217
# Also catches case where periodicity = false
218
periodicity = ntuple(_ -> false, NDIMS)
219
else
220
# Default case if periodicity is an iterable
221
periodicity = Tuple(periodicity)
222
end
223
224
basis = LobattoLegendreBasis(RealT, polydeg)
225
nodes = basis.nodes
226
tree_node_coordinates = Array{RealT, NDIMS + 2}(undef, NDIMS,
227
ntuple(_ -> length(nodes),
228
NDIMS)...,
229
prod(trees_per_dimension))
230
calc_tree_node_coordinates!(tree_node_coordinates, nodes, mapping,
231
trees_per_dimension)
232
233
# p4est_connectivity_new_brick has trees in Z-order, so use our own function for this
234
connectivity = connectivity_structured(trees_per_dimension..., periodicity)
235
236
p4est = new_p4est(connectivity, initial_refinement_level)
237
238
# Non-periodic boundaries
239
boundary_names = fill(Symbol("---"), 2 * NDIMS, prod(trees_per_dimension))
240
241
structured_boundary_names!(boundary_names, trees_per_dimension, periodicity)
242
243
return P4estMesh{NDIMS}(p4est, tree_node_coordinates, nodes,
244
boundary_names, "", unsaved_changes,
245
p4est_partition_allow_for_coarsening)
246
end
247
248
# 2D version
249
function structured_boundary_names!(boundary_names, trees_per_dimension::NTuple{2},
250
periodicity)
251
linear_indices = LinearIndices(trees_per_dimension)
252
253
# Boundaries in x-direction
254
if !periodicity[1]
255
for cell_y in 1:trees_per_dimension[2]
256
tree = linear_indices[1, cell_y]
257
boundary_names[1, tree] = :x_neg
258
259
tree = linear_indices[end, cell_y]
260
boundary_names[2, tree] = :x_pos
261
end
262
end
263
264
# Boundaries in y-direction
265
if !periodicity[2]
266
for cell_x in 1:trees_per_dimension[1]
267
tree = linear_indices[cell_x, 1]
268
boundary_names[3, tree] = :y_neg
269
270
tree = linear_indices[cell_x, end]
271
boundary_names[4, tree] = :y_pos
272
end
273
end
274
end
275
276
# 3D version
277
function structured_boundary_names!(boundary_names, trees_per_dimension::NTuple{3},
278
periodicity)
279
linear_indices = LinearIndices(trees_per_dimension)
280
281
# Boundaries in x-direction
282
if !periodicity[1]
283
for cell_z in 1:trees_per_dimension[3], cell_y in 1:trees_per_dimension[2]
284
tree = linear_indices[1, cell_y, cell_z]
285
boundary_names[1, tree] = :x_neg
286
287
tree = linear_indices[end, cell_y, cell_z]
288
boundary_names[2, tree] = :x_pos
289
end
290
end
291
292
# Boundaries in y-direction
293
if !periodicity[2]
294
for cell_z in 1:trees_per_dimension[3], cell_x in 1:trees_per_dimension[1]
295
tree = linear_indices[cell_x, 1, cell_z]
296
boundary_names[3, tree] = :y_neg
297
298
tree = linear_indices[cell_x, end, cell_z]
299
boundary_names[4, tree] = :y_pos
300
end
301
end
302
303
# Boundaries in z-direction
304
if !periodicity[3]
305
for cell_y in 1:trees_per_dimension[2], cell_x in 1:trees_per_dimension[1]
306
tree = linear_indices[cell_x, cell_y, 1]
307
boundary_names[5, tree] = :z_neg
308
309
tree = linear_indices[cell_x, cell_y, end]
310
boundary_names[6, tree] = :z_pos
311
end
312
end
313
end
314
315
"""
316
P4estMesh{NDIMS}(meshfile::String;
317
mapping=nothing, polydeg=1, RealT=Float64,
318
initial_refinement_level=0, unsaved_changes=true,
319
p4est_partition_allow_for_coarsening=true,
320
boundary_symbols = nothing)
321
322
Main mesh constructor for the `P4estMesh` that imports an unstructured, conforming
323
mesh from an Abaqus mesh file (`.inp`). Each element of the conforming mesh parsed
324
from the `meshfile` is created as a [`p4est`](https://github.com/cburstedde/p4est)
325
tree datatype.
326
327
To create a curved/higher-order unstructured mesh `P4estMesh` two strategies are available:
328
329
- `p4est_mesh_from_hohqmesh_abaqus`: High-order (curved) boundary information created by
330
[`HOHQMesh.jl`](https://github.com/trixi-framework/HOHQMesh.jl) is
331
available in the `meshfile`. The mesh polynomial degree `polydeg`
332
of the boundaries is provided from the `meshfile`. The computation of
333
the mapped tree coordinates is done with transfinite interpolation
334
with linear blending similar to [`UnstructuredMesh2D`](@ref). Boundary name
335
information is also parsed from the `meshfile` such that different boundary
336
conditions can be set at each named boundary on a given tree.
337
- `p4est_mesh_from_standard_abaqus`: By default, with `mapping=nothing` and `polydeg=1`, this creates a
338
straight-sided from the information parsed from the `meshfile`. If a mapping
339
function is specified then it computes the mapped tree coordinates via polynomial
340
interpolants with degree `polydeg`. The mesh created by this function will only
341
have one boundary `:all` if `boundary_symbols` is not specified.
342
If `boundary_symbols` is specified the `meshfile` will be parsed for nodesets defining
343
the boundary nodes from which boundary edges (2D) and faces (3D) will be assigned.
344
345
Note that the `mapping` and `polydeg` keyword arguments are only used by the `p4est_mesh_from_standard_abaqus`
346
function. The `p4est_mesh_from_hohqmesh_abaqus` function obtains the mesh `polydeg` directly from the `meshfile`
347
and constructs the transfinite mapping internally.
348
349
The particular strategy is selected according to the header present in the `meshfile` where
350
the constructor checks whether or not the `meshfile` was created with
351
[HOHQMesh.jl](https://github.com/trixi-framework/HOHQMesh.jl).
352
If the Abaqus file header is not present in the `meshfile` then the `P4estMesh` is created
353
with the function `p4est_mesh_from_standard_abaqus`.
354
355
The default keyword argument `initial_refinement_level=0` corresponds to a forest
356
where the number of trees is the same as the number of elements in the original `meshfile`.
357
Increasing the `initial_refinement_level` allows one to uniformly refine the base mesh given
358
in the `meshfile` to create a forest with more trees before the simulation begins.
359
For example, if a two-dimensional base mesh contains 25 elements then setting
360
`initial_refinement_level=1` creates an initial forest of `2^2 * 25 = 100` trees.
361
362
# Arguments
363
- `meshfile::String`: an uncurved Abaqus mesh file that can be imported by `p4est`.
364
- `mapping`: a function of `NDIMS` variables to describe the mapping that transforms
365
the imported mesh to the physical domain. Use `nothing` for the identity map.
366
- `polydeg::Integer`: polynomial degree used to store the geometry of the mesh.
367
The mapping will be approximated by an interpolation polynomial
368
of the specified degree for each tree.
369
The default of `1` creates an uncurved geometry. Use a higher value if the mapping
370
will curve the imported uncurved mesh.
371
- `RealT::Type`: the type that should be used for coordinates.
372
- `initial_refinement_level::Integer`: refine the mesh uniformly to this level before the simulation starts.
373
- `unsaved_changes::Bool`: if set to `true`, the mesh will be saved to a mesh file.
374
- `p4est_partition_allow_for_coarsening::Bool`: Must be `true` when using AMR to make mesh adaptivity
375
independent of domain partitioning. Should be `false` for static meshes
376
to permit more fine-grained partitioning.
377
- `boundary_symbols::Vector{Symbol}`: A vector of symbols that correspond to the boundary names in the `meshfile`.
378
If `nothing` is passed then all boundaries are named `:all`.
379
"""
380
function P4estMesh{NDIMS}(meshfile::String;
381
mapping = nothing, polydeg = 1, RealT = Float64,
382
initial_refinement_level = 0, unsaved_changes = true,
383
p4est_partition_allow_for_coarsening = true,
384
boundary_symbols = nothing) where {NDIMS}
385
# Prevent `p4est` from crashing Julia if the file doesn't exist
386
@assert isfile(meshfile)
387
388
# Read in the Header of the meshfile to determine which constructor is appropriate
389
header = open(meshfile, "r") do io
390
readline(io) # *Header of the Abaqus file; discarded
391
return readline(io) # Readin the actual header information
392
end
393
394
# Check if the meshfile was generated using HOHQMesh
395
if header == " File created by HOHQMesh"
396
# Mesh curvature and boundary naming is handled with additional information available in meshfile
397
p4est, tree_node_coordinates, nodes, boundary_names = p4est_mesh_from_hohqmesh_abaqus(meshfile,
398
initial_refinement_level,
399
NDIMS,
400
RealT)
401
else
402
# Mesh curvature is handled directly by applying the mapping keyword argument
403
p4est, tree_node_coordinates, nodes, boundary_names = p4est_mesh_from_standard_abaqus(meshfile,
404
mapping,
405
polydeg,
406
initial_refinement_level,
407
NDIMS,
408
RealT,
409
boundary_symbols)
410
end
411
412
return P4estMesh{NDIMS}(p4est, tree_node_coordinates, nodes,
413
boundary_names, "", unsaved_changes,
414
p4est_partition_allow_for_coarsening)
415
end
416
417
# Wrapper for `p4est_connectivity_from_hohqmesh_abaqus`. The latter is used
418
# by `T8codeMesh`, too.
419
function p4est_mesh_from_hohqmesh_abaqus(meshfile, initial_refinement_level,
420
n_dimensions, RealT)
421
connectivity, tree_node_coordinates, nodes, boundary_names = p4est_connectivity_from_hohqmesh_abaqus(meshfile,
422
n_dimensions,
423
RealT)
424
425
p4est = new_p4est(connectivity, initial_refinement_level)
426
427
return p4est, tree_node_coordinates, nodes, boundary_names
428
end
429
430
# Wrapper for `p4est_connectivity_from_standard_abaqus`. The latter is used
431
# by `T8codeMesh`, too.
432
function p4est_mesh_from_standard_abaqus(meshfile, mapping, polydeg,
433
initial_refinement_level, n_dimensions, RealT,
434
boundary_symbols)
435
connectivity, tree_node_coordinates, nodes, boundary_names = p4est_connectivity_from_standard_abaqus(meshfile,
436
mapping,
437
polydeg,
438
n_dimensions,
439
RealT,
440
boundary_symbols)
441
442
p4est = new_p4est(connectivity, initial_refinement_level)
443
444
return p4est, tree_node_coordinates, nodes, boundary_names
445
end
446
447
# Create the mesh connectivity, mapped node coordinates within each tree, reference nodes in [-1,1]
448
# and a list of boundary names for the `P4estMesh`. High-order boundary curve information as well as
449
# the boundary names on each tree are provided by the `meshfile` created by
450
# [`HOHQMesh.jl`](https://github.com/trixi-framework/HOHQMesh.jl).
451
function p4est_connectivity_from_hohqmesh_abaqus(meshfile,
452
n_dimensions, RealT)
453
# Create the mesh connectivity using `p4est`
454
connectivity = read_inp_p4est(meshfile, Val(n_dimensions))
455
connectivity_pw = PointerWrapper(connectivity)
456
457
# These need to be of the type Int for unsafe_wrap below to work
458
n_trees::Int = connectivity_pw.num_trees[] # = number elements
459
n_vertices::Int = connectivity_pw.num_vertices[]
460
461
# Extract a copy of the element vertices to compute the tree node coordinates
462
# `vertices` stores coordinates of all three dimensions (even for the 2D case)
463
# since the Abaqus `.inp` format always stores 3D coordinates.
464
vertices = unsafe_wrap(Array, connectivity_pw.vertices, (3, n_vertices))
465
466
# Readin all the information from the mesh file into a string array
467
file_lines = readlines(open(meshfile))
468
469
# Get the file index where the mesh polynomial degree is given in the meshfile
470
file_idx = findfirst(contains("** mesh polynomial degree"), file_lines)
471
472
# Get the polynomial order of the mesh boundary information
473
current_line = split(file_lines[file_idx])
474
mesh_polydeg = parse(Int, current_line[6])
475
mesh_nnodes = mesh_polydeg + 1
476
477
# Create the Chebyshev-Gauss-Lobatto nodes used by HOHQMesh to represent the boundaries
478
cheby_nodes, _ = chebyshev_gauss_lobatto_nodes_weights(mesh_nnodes)
479
nodes = SVector{mesh_nnodes}(cheby_nodes)
480
481
# Allocate the memory for the tree node coordinates
482
tree_node_coordinates = Array{RealT, n_dimensions + 2}(undef, n_dimensions,
483
ntuple(_ -> length(nodes),
484
n_dimensions)...,
485
n_trees)
486
487
# Compute the tree node coordinates and return the updated file index
488
file_idx = calc_tree_node_coordinates!(tree_node_coordinates, file_lines, nodes,
489
vertices, RealT)
490
491
# Allocate the memory for the boundary labels
492
boundary_names = Array{Symbol}(undef, (2 * n_dimensions, n_trees))
493
494
# Read in the boundary names from the last portion of the meshfile
495
# Note here the boundary names where "---" means an internal connection
496
for tree in 1:n_trees
497
current_line = split(file_lines[file_idx])
498
boundary_names[:, tree] = map(Symbol, current_line[2:end])
499
file_idx += 1
500
end
501
502
return connectivity, tree_node_coordinates, nodes, boundary_names
503
end
504
505
# This removes irrelevant elements from the meshfile, i.e., trusses/beams which are 1D elements.
506
# Thus, if `n_dimensions == 2`, only quads are kept, and if `n_dimensions == 3`,
507
# only hexes are kept, i.e., in that case also quads are removed.
508
function preprocess_standard_abaqus(meshfile,
509
linear_quads, linear_hexes,
510
quadratic_quads, quadratic_hexes,
511
n_dimensions)
512
meshfile_preproc = replace(meshfile, ".inp" => "_preproc.inp")
513
514
# Line number where the element section begins
515
elements_begin_idx = 0
516
# Line number where the node and element sets begin (if present)
517
sets_begin_idx = 0
518
519
open(meshfile, "r") do infile
520
# Copy header and node data to pre-processed file
521
open(meshfile_preproc, "w") do outfile
522
for (line_index, line) in enumerate(eachline(infile))
523
println(outfile, line)
524
if occursin("******* E L E M E N T S *************", line)
525
elements_begin_idx = line_index + 1
526
break
527
end
528
end
529
end
530
531
# Find the line number where the node and element sets begin
532
for (line_index, line) in enumerate(eachline(infile))
533
if startswith(line, "*ELSET") || startswith(line, "*NSET")
534
sets_begin_idx = line_index
535
break
536
end
537
end
538
end
539
# Catch the case where there are no node or element sets
540
if sets_begin_idx == 0
541
sets_begin_idx = length(readlines(meshfile)) + 1
542
else
543
# Need to add `elements_begin_idx - 1` since the loop above starts at `elements_begin_idx`
544
sets_begin_idx += elements_begin_idx - 1
545
end
546
547
element_index = 0
548
preproc_element_section_lines = 0
549
open(meshfile, "r") do infile
550
open(meshfile_preproc, "a") do outfile
551
print_following_lines = false
552
553
for (line_index, line) in enumerate(eachline(infile))
554
# Act only in the element section
555
if elements_begin_idx <= line_index < sets_begin_idx
556
# Check if a new element type/element set is defined
557
if startswith(line, "*ELEMENT")
558
# Retrieve element type
559
current_element_type = match(r"\*ELEMENT, type=([^,]+)", line).captures[1]
560
561
# Keep only quads (2D) or hexes (3D), i.e., eliminate all other elements
562
# like trusses (2D/3D) or quads (3D).
563
if n_dimensions == 2 &&
564
(occursin(linear_quads, current_element_type) ||
565
occursin(quadratic_quads, current_element_type))
566
print_following_lines = true
567
println(outfile, line)
568
preproc_element_section_lines += 1
569
elseif n_dimensions == 3 &&
570
(occursin(linear_hexes, current_element_type) ||
571
occursin(quadratic_hexes, current_element_type))
572
print_following_lines = true
573
println(outfile, line)
574
preproc_element_section_lines += 1
575
else
576
print_following_lines = false
577
end
578
else # Element data in line
579
if print_following_lines
580
element_index += 1
581
parts = split(line, ',')
582
# Exchange element index
583
parts[1] = string(element_index)
584
println(outfile, join(parts, ','))
585
preproc_element_section_lines += 1
586
end
587
end
588
end
589
590
# Print node and element sets as they are
591
if line_index >= sets_begin_idx
592
println(outfile, line)
593
end
594
end
595
end
596
end
597
# Adjust `sets_begin_idx` to correct line number after removing unnecessary elements
598
sets_begin_idx = elements_begin_idx + preproc_element_section_lines
599
600
return elements_begin_idx, sets_begin_idx
601
end
602
603
# p4est can handle only linear elements. This function checks the `meshfile_pre_proc`
604
# for quadratic elements (highest order supported by standard Abaqus) and
605
# replaces them with linear elements. The higher-order (quadratic) boundaries are handled
606
# "internally" by Trixi as for the HOHQMesh-Abaqus case.
607
function preprocess_standard_abaqus_for_p4est(meshfile_pre_proc,
608
linear_quads, linear_hexes,
609
quadratic_quads, quadratic_hexes,
610
elements_begin_idx,
611
sets_begin_idx)
612
meshfile_p4est_rdy = replace(meshfile_pre_proc,
613
"_preproc.inp" => "_p4est_ready.inp")
614
order = 1 # Assume linear elements by default
615
616
# Some useful function-wide variables
617
current_element_type = ""
618
new_line = ""
619
620
open(meshfile_pre_proc, "r") do infile
621
open(meshfile_p4est_rdy, "w") do outfile
622
for (line_index, line) in enumerate(eachline(infile))
623
# Copy header and node data
624
if line_index < elements_begin_idx || line_index >= sets_begin_idx
625
println(outfile, line)
626
else # Element data
627
# Check if a new element type/element set is defined
628
if startswith(line, "*ELEMENT")
629
# Retrieve element type
630
current_element_type = match(r"\*ELEMENT, type=([^,]+)", line).captures[1]
631
632
# Check for linear elements - then we just copy the line
633
if occursin(linear_quads, current_element_type) ||
634
occursin(linear_hexes, current_element_type)
635
println(outfile, line)
636
else # Quadratic element - replace with linear
637
order = 2
638
if occursin(quadratic_quads, current_element_type)
639
linear_quad_type = "CPS4"
640
new_line = replace(line,
641
current_element_type => linear_quad_type)
642
elseif occursin(quadratic_hexes, current_element_type)
643
linear_hex_type = "C3D8"
644
new_line = replace(line,
645
current_element_type => linear_hex_type)
646
end
647
println(outfile, new_line)
648
end
649
else # Element data in line
650
# Check for linear elements - then we just copy the line
651
if occursin(linear_quads, current_element_type) ||
652
occursin(linear_hexes, current_element_type)
653
println(outfile, line)
654
else
655
parts = split(line, ',')
656
if occursin(quadratic_quads, current_element_type)
657
# Print the first (element), second to fifth (vertices 1-4) indices to file.
658
# For node order of quadratic (second-order) quads, see
659
# http://130.149.89.49:2080/v2016/books/usb/default.htm?startat=pt06ch28s01ael02.html
660
new_line = join(parts[1:5], ',')
661
elseif occursin(quadratic_hexes, current_element_type)
662
# Print the first (element), second to ninth (vertices 1-8) indices to file.
663
# The node order is fortunately the same for hexes/bricks of type "C3D27", see
664
# http://130.149.89.49:2080/v2016/books/usb/default.htm?startat=pt06ch28s01ael03.html
665
new_line = join(parts[1:9], ',')
666
end
667
println(outfile, new_line)
668
end
669
end
670
end
671
end
672
end
673
end
674
675
return order
676
end
677
678
# Read all nodes (not only vertices, i.e., endpoints of elements) into a dict.
679
# Those are required to enable higher-order boundaries for quadratic elements.
680
function read_nodes_standard_abaqus(meshfile, n_dimensions,
681
elements_begin_idx, RealT)
682
mesh_nodes = Dict{Int, SVector{n_dimensions, RealT}}()
683
684
node_section = false
685
open(meshfile, "r") do infile
686
for (line_index, line) in enumerate(eachline(infile))
687
if line_index < elements_begin_idx
688
if startswith(line, "*NODE")
689
node_section = true
690
end
691
if node_section
692
parts = split(line, ',')
693
if length(parts) >= n_dimensions + 1
694
node_id = parse(Int, parts[1])
695
coordinates = SVector{n_dimensions, RealT}([parse(RealT,
696
parts[i])
697
for i in 2:(n_dimensions + 1)])
698
mesh_nodes[node_id] = coordinates
699
end
700
end
701
end
702
end
703
end
704
705
return mesh_nodes
706
end
707
708
# Create the mesh connectivity, mapped node coordinates within each tree, reference nodes in [-1,1]
709
# and a list of boundary names for the `P4estMesh`. For linear meshes, the tree node coordinates are
710
# computed according to the `mapping` passed to this function using polynomial interpolants of degree `polydeg`.
711
# For quadratic (second-order) meshes, the tree node coordinates are read from the meshfile, similar as for
712
# `p4est_connectivity_from_hohqmesh_abaqus`.
713
function p4est_connectivity_from_standard_abaqus(meshfile, mapping, polydeg,
714
n_dimensions,
715
RealT,
716
boundary_symbols)
717
718
### Regular expressions for Abaqus element types ###
719
720
# These are the standard Abaqus linear quads.
721
# Note that there are many(!) more variants designed for specific purposes in the Abaqus solver, see
722
# http://130.149.89.49:2080/v2016/books/usb/default.htm?startat=pt10eli01.html
723
# To keep it simple we support only basic quads, membranes, and shells here.
724
linear_quads = r"^(CPE4|CPEG4|CPS4|M3D4|S4|SFM3D4).*$"
725
726
# Same logic as for linear quads:
727
# Support only basic quads, membranes, and shells here.
728
quadratic_quads = r"^(CPE8|CPS8|CAX8|S8|M3D8|M3D9).*$"
729
730
# In 3D only standard hexahedra/bricks are supported.
731
linear_hexes = r"^(C3D8).*$"
732
quadratic_hexes = r"^(C3D27).*$"
733
734
meshfile_preproc = replace(meshfile, ".inp" => "_preproc.inp")
735
736
# Define variables that are retrieved in the MPI-parallel case on root and then bcasted
737
elements_begin_idx = -1
738
sets_begin_idx = -1
739
mesh_polydeg = 1
740
741
if !mpi_isparallel() || (mpi_isparallel() && mpi_isroot())
742
# Preprocess the meshfile to remove lower-dimensional elements
743
elements_begin_idx, sets_begin_idx = preprocess_standard_abaqus(meshfile,
744
linear_quads,
745
linear_hexes,
746
quadratic_quads,
747
quadratic_hexes,
748
n_dimensions)
749
750
# Copy of mesh for p4est with linear elements only
751
mesh_polydeg = preprocess_standard_abaqus_for_p4est(meshfile_preproc,
752
linear_quads,
753
linear_hexes,
754
quadratic_quads,
755
quadratic_hexes,
756
elements_begin_idx,
757
sets_begin_idx)
758
end
759
760
# Broadcast from meshfile retrieved variables across all MPI ranks
761
if mpi_isparallel()
762
if mpi_isroot()
763
MPI.Bcast!(Ref(elements_begin_idx), mpi_root(), mpi_comm())
764
MPI.Bcast!(Ref(sets_begin_idx), mpi_root(), mpi_comm())
765
MPI.Bcast!(Ref(mesh_polydeg), mpi_root(), mpi_comm())
766
else
767
elements_begin_idx = MPI.Bcast!(Ref(0), mpi_root(), mpi_comm())[]
768
sets_begin_idx = MPI.Bcast!(Ref(0), mpi_root(), mpi_comm())[]
769
mesh_polydeg = MPI.Bcast!(Ref(0), mpi_root(), mpi_comm())[]
770
end
771
end
772
773
# Create the mesh connectivity using `p4est`
774
meshfile_p4est_rdy = replace(meshfile, ".inp" => "_p4est_ready.inp")
775
connectivity = read_inp_p4est(meshfile_p4est_rdy, Val(n_dimensions))
776
connectivity_pw = PointerWrapper(connectivity)
777
778
# These need to be of the type Int for unsafe_wrap below to work
779
n_trees::Int = connectivity_pw.num_trees[] # = number elements
780
n_vertices::Int = connectivity_pw.num_vertices[]
781
782
# Extract a copy of the element vertices to compute the tree node coordinates
783
# `vertices` store coordinates of all three dimensions (even for the 2D case)
784
# since the Abaqus `.inp` format always stores 3D coordinates.
785
vertices = unsafe_wrap(Array, connectivity_pw.vertices, (3, n_vertices))
786
787
tree_to_vertex = unsafe_wrap(Array, connectivity_pw.tree_to_vertex,
788
(2^n_dimensions, n_trees))
789
790
basis = LobattoLegendreBasis(RealT, polydeg)
791
nodes = basis.nodes
792
793
# The highest supported element order is quadratic (second-order) in the standard Abaqus format.
794
# Thus, this check is equivalent to checking for higher-order boundary information.
795
if mesh_polydeg == 2
796
mesh_nnodes = mesh_polydeg + 1 # = 3
797
# Note: We ASSUME that the additional node between the end-vertices lies
798
# on the center on that line, such that we can use Chebyshev-Gauss-Lobatto nodes!
799
# For polydeg = 2, we have the 3 nodes [-1, 0, 1] (within the reference element).
800
# Note that these coincide for polydeg = 2 with the Legendre-Gauss-Lobatto nodes.
801
cheby_nodes, _ = chebyshev_gauss_lobatto_nodes_weights(mesh_nnodes)
802
nodes = SVector{mesh_nnodes}(cheby_nodes)
803
end
804
805
tree_node_coordinates = Array{RealT, n_dimensions + 2}(undef, n_dimensions,
806
ntuple(_ -> length(nodes),
807
n_dimensions)...,
808
n_trees)
809
810
# No nonlinearity in the mesh. Rely on the mapping function to realize the curvature.
811
if mesh_polydeg == 1
812
calc_tree_node_coordinates!(tree_node_coordinates, nodes, mapping,
813
vertices, tree_to_vertex)
814
else # mesh_polydeg = 2 => Nonlinearity in supplied mesh
815
mesh_nodes = read_nodes_standard_abaqus(meshfile_preproc, n_dimensions,
816
elements_begin_idx, RealT)
817
818
# Extract element section from pre-processed Abaqus meshfile
819
element_lines = readlines(open(meshfile_preproc))[elements_begin_idx:(sets_begin_idx - 1)]
820
821
if n_dimensions == 2
822
calc_tree_node_coordinates!(tree_node_coordinates, element_lines,
823
nodes, vertices, RealT,
824
linear_quads, mesh_nodes)
825
else # n_dimensions == 3
826
calc_tree_node_coordinates!(tree_node_coordinates, element_lines,
827
nodes, vertices, RealT,
828
linear_hexes, mesh_nodes)
829
end
830
end
831
832
if boundary_symbols === nothing
833
# There's no simple and generic way to distinguish boundaries without any information given.
834
# Name all of them :all.
835
boundary_names = fill(:all, 2 * n_dimensions, n_trees)
836
else # Boundary information given
837
# Read in nodes belonging to boundaries
838
node_set_dict = parse_node_sets(meshfile_p4est_rdy, boundary_symbols)
839
# Read in all elements with associated nodes to specify the boundaries
840
element_node_matrix = parse_elements(meshfile_p4est_rdy, n_trees, n_dimensions,
841
elements_begin_idx, sets_begin_idx)
842
843
# Initialize boundary information matrix with symbol for no boundary / internal connection
844
boundary_names = fill(Symbol("---"), 2 * n_dimensions, n_trees)
845
846
# Fill `boundary_names` such that it can be processed by p4est
847
assign_boundaries_standard_abaqus!(boundary_names, n_trees,
848
element_node_matrix, node_set_dict,
849
Val(n_dimensions))
850
end
851
852
return connectivity, tree_node_coordinates, nodes, boundary_names
853
end
854
855
function parse_elements(meshfile, n_trees, n_dims,
856
elements_begin_idx,
857
sets_begin_idx)
858
# 2D quads: 4 nodes + element index, 3D hexes: 8 nodes + element index
859
expected_content_length = n_dims == 2 ? 5 : 9
860
861
element_node_matrix = Matrix{Int64}(undef, n_trees, expected_content_length - 1)
862
el_list_follows = false
863
tree_id = 1
864
865
open(meshfile, "r") do infile
866
for (line_index, line) in enumerate(eachline(infile))
867
if elements_begin_idx <= line_index < sets_begin_idx
868
if startswith(line, "*ELEMENT")
869
el_list_follows = true
870
elseif el_list_follows
871
content = split(line, ",")
872
if length(content) == expected_content_length # Check that we still read in connectivity data
873
content_int = parse.(Int64, content)
874
# Add constituent nodes to the element_node_matrix.
875
# Important: Do not use index from the Abaqus file, but the one from p4est.
876
element_node_matrix[tree_id, :] = content_int[2:end] # First entry is element id
877
tree_id += 1
878
else # Processed all elements for this `ELSET`
879
el_list_follows = false
880
end
881
end
882
end
883
end
884
end
885
886
return element_node_matrix
887
end
888
889
function parse_node_sets(meshfile, boundary_symbols)
890
nodes_dict = Dict{Symbol, Set{Int64}}()
891
current_symbol = nothing
892
current_nodes = Int64[]
893
894
open(meshfile, "r") do file
895
for line in eachline(file)
896
# Check if the line contains nodes assembled in a special set, i.e., a physical boundary
897
if startswith(line, "*NSET,NSET=")
898
# Safe the previous nodeset
899
if current_symbol !== nothing
900
nodes_dict[current_symbol] = current_nodes
901
end
902
903
current_symbol = Symbol(split(line, "=")[2])
904
if current_symbol in boundary_symbols
905
# New nodeset
906
current_nodes = Set{Int64}()
907
else # Read only boundary node sets
908
current_symbol = nothing
909
end
910
elseif current_symbol !== nothing # Read only if there was already a nodeset specified
911
try # Check if line contains nodes
912
# There is always a trailing comma, remove the corresponding empty string
913
union!(current_nodes, parse.(Int64, split(line, ",")[1:(end - 1)]))
914
catch # Something different, stop reading in nodes
915
# If parsing fails, set current_symbol to nothing
916
nodes_dict[current_symbol] = current_nodes
917
current_symbol = nothing
918
end
919
end
920
end
921
# Safe the previous nodeset
922
if current_symbol !== nothing
923
nodes_dict[current_symbol] = current_nodes
924
end
925
end
926
927
for symbol in boundary_symbols
928
if !haskey(nodes_dict, symbol)
929
@warn "No nodes found for nodeset :" * "$symbol" * " !"
930
end
931
end
932
933
return nodes_dict
934
end
935
936
# This function assigns the edges of elements to boundaries by
937
# checking if the nodes that define the edges are part of nodesets which correspond to boundaries.
938
function assign_boundaries_standard_abaqus!(boundary_names, n_trees,
939
element_node_matrix, node_set_dict,
940
::Val{2}) # 2D version
941
@threaded for tree in 1:n_trees
942
tree_nodes = element_node_matrix[tree, :]
943
# For node labeling, see
944
# https://docs.software.vt.edu/abaqusv2022/English/SIMACAEELMRefMap/simaelm-r-2delem.htm#simaelm-r-2delem-t-nodedef1
945
# and search for "Node ordering and face numbering on elements"
946
for boundary in keys(node_set_dict) # Loop over specified boundaries
947
# Check bottom edge
948
if tree_nodes[1] in node_set_dict[boundary] &&
949
tree_nodes[2] in node_set_dict[boundary]
950
# Bottom boundary is position 3 in p4est indexing
951
boundary_names[3, tree] = boundary
952
end
953
# Check right edge
954
if tree_nodes[2] in node_set_dict[boundary] &&
955
tree_nodes[3] in node_set_dict[boundary]
956
# Right boundary is position 2 in p4est indexing
957
boundary_names[2, tree] = boundary
958
end
959
# Check top edge
960
if tree_nodes[3] in node_set_dict[boundary] &&
961
tree_nodes[4] in node_set_dict[boundary]
962
# Top boundary is position 4 in p4est indexing
963
boundary_names[4, tree] = boundary
964
end
965
# Check left edge
966
if tree_nodes[4] in node_set_dict[boundary] &&
967
tree_nodes[1] in node_set_dict[boundary]
968
# Left boundary is position 1 in p4est indexing
969
boundary_names[1, tree] = boundary
970
end
971
end
972
end
973
974
return boundary_names
975
end
976
977
# This function assigns the edges of elements to boundaries by
978
# checking if the nodes that define the faces are part of nodesets which correspond to boundaries.
979
function assign_boundaries_standard_abaqus!(boundary_names, n_trees,
980
element_node_matrix, node_set_dict,
981
::Val{3}) # 3D version
982
@threaded for tree in 1:n_trees
983
tree_nodes = element_node_matrix[tree, :]
984
# For node labeling, see
985
# https://web.mit.edu/calculix_v2.7/CalculiX/ccx_2.7/doc/ccx/node26.html
986
for boundary in keys(node_set_dict) # Loop over specified boundaries
987
# Check "front face" (y_min)
988
if tree_nodes[1] in node_set_dict[boundary] &&
989
tree_nodes[2] in node_set_dict[boundary] &&
990
tree_nodes[5] in node_set_dict[boundary] &&
991
tree_nodes[6] in node_set_dict[boundary]
992
# Front face is position 3 in p4est indexing
993
boundary_names[3, tree] = boundary
994
end
995
# Check "back face" (y_max)
996
if tree_nodes[3] in node_set_dict[boundary] &&
997
tree_nodes[4] in node_set_dict[boundary] &&
998
tree_nodes[7] in node_set_dict[boundary] &&
999
tree_nodes[8] in node_set_dict[boundary]
1000
# Front face is position 4 in p4est indexing
1001
boundary_names[4, tree] = boundary
1002
end
1003
# Check "left face" (x_min)
1004
if tree_nodes[1] in node_set_dict[boundary] &&
1005
tree_nodes[4] in node_set_dict[boundary] &&
1006
tree_nodes[5] in node_set_dict[boundary] &&
1007
tree_nodes[8] in node_set_dict[boundary]
1008
# Left face is position 1 in p4est indexing
1009
boundary_names[1, tree] = boundary
1010
end
1011
# Check "right face" (x_max)
1012
if tree_nodes[2] in node_set_dict[boundary] &&
1013
tree_nodes[3] in node_set_dict[boundary] &&
1014
tree_nodes[6] in node_set_dict[boundary] &&
1015
tree_nodes[7] in node_set_dict[boundary]
1016
# Right face is position 2 in p4est indexing
1017
boundary_names[2, tree] = boundary
1018
end
1019
# Check "bottom face" (z_min)
1020
if tree_nodes[1] in node_set_dict[boundary] &&
1021
tree_nodes[2] in node_set_dict[boundary] &&
1022
tree_nodes[3] in node_set_dict[boundary] &&
1023
tree_nodes[4] in node_set_dict[boundary]
1024
# Bottom face is position 5 in p4est indexing
1025
boundary_names[5, tree] = boundary
1026
end
1027
# Check "top face" (z_max)
1028
if tree_nodes[5] in node_set_dict[boundary] &&
1029
tree_nodes[6] in node_set_dict[boundary] &&
1030
tree_nodes[7] in node_set_dict[boundary] &&
1031
tree_nodes[8] in node_set_dict[boundary]
1032
# Top face is position 6 in p4est indexing
1033
boundary_names[6, tree] = boundary
1034
end
1035
end
1036
end
1037
1038
return boundary_names
1039
end
1040
1041
"""
1042
P4estMeshCubedSphere(trees_per_face_dimension, layers, inner_radius, thickness;
1043
polydeg, RealT=Float64,
1044
initial_refinement_level=0, unsaved_changes=true,
1045
p4est_partition_allow_for_coarsening=true)
1046
1047
Build a "Cubed Sphere" mesh as `P4estMesh` with
1048
`6 * trees_per_face_dimension^2 * layers` trees.
1049
1050
The mesh will have two boundaries, `:inside` and `:outside`.
1051
1052
# Arguments
1053
- `trees_per_face_dimension::Integer`: the number of trees in the first two local dimensions of
1054
each face.
1055
- `layers::Integer`: the number of trees in the third local dimension of each face, i.e., the number
1056
of layers of the sphere.
1057
- `inner_radius::RealT`: the inner radius of the sphere.
1058
- `thickness::RealT`: the thickness of the spherical shell. The outer radius will be `inner_radius + thickness`.
1059
- `polydeg::Integer`: polynomial degree used to store the geometry of the mesh.
1060
The mapping will be approximated by an interpolation polynomial
1061
of the specified degree for each tree.
1062
- `RealT::Type`: the type that should be used for coordinates.
1063
- `initial_refinement_level::Integer`: refine the mesh uniformly to this level before the simulation starts.
1064
- `unsaved_changes::Bool`: if set to `true`, the mesh will be saved to a mesh file.
1065
- `p4est_partition_allow_for_coarsening::Bool`: Must be `true` when using AMR to make mesh adaptivity
1066
independent of domain partitioning. Should be `false` for static meshes
1067
to permit more fine-grained partitioning.
1068
"""
1069
function P4estMeshCubedSphere(trees_per_face_dimension, layers, inner_radius, thickness;
1070
polydeg, RealT = Float64,
1071
initial_refinement_level = 0, unsaved_changes = true,
1072
p4est_partition_allow_for_coarsening = true)
1073
connectivity = connectivity_cubed_sphere(trees_per_face_dimension, layers)
1074
1075
n_trees = 6 * trees_per_face_dimension^2 * layers
1076
1077
basis = LobattoLegendreBasis(RealT, polydeg)
1078
nodes = basis.nodes
1079
1080
tree_node_coordinates = Array{RealT, 5}(undef, 3,
1081
ntuple(_ -> length(nodes), 3)...,
1082
n_trees)
1083
calc_tree_node_coordinates!(tree_node_coordinates, nodes,
1084
trees_per_face_dimension, layers,
1085
inner_radius, thickness)
1086
1087
p4est = new_p4est(connectivity, initial_refinement_level)
1088
1089
boundary_names = fill(Symbol("---"), 2 * 3, n_trees)
1090
boundary_names[5, :] .= Symbol("inside")
1091
boundary_names[6, :] .= Symbol("outside")
1092
1093
return P4estMesh{3}(p4est, tree_node_coordinates, nodes,
1094
boundary_names, "", unsaved_changes,
1095
p4est_partition_allow_for_coarsening)
1096
end
1097
1098
# Create a new p4est_connectivity that represents a structured rectangle.
1099
# Similar to p4est_connectivity_new_brick, but doesn't use Morton order.
1100
# This order makes `calc_tree_node_coordinates!` below and the calculation
1101
# of `boundary_names` above easier but is irrelevant otherwise.
1102
# 2D version
1103
function connectivity_structured(n_cells_x, n_cells_y, periodicity)
1104
linear_indices = LinearIndices((n_cells_x, n_cells_y))
1105
1106
# Vertices represent the coordinates of the forest. This is used by `p4est`
1107
# to write VTK files.
1108
# Trixi.jl doesn't use the coordinates from `p4est`, so the vertices can be empty.
1109
n_vertices = 0
1110
n_trees = n_cells_x * n_cells_y
1111
# No corner connectivity is needed
1112
n_corners = 0
1113
vertices = C_NULL
1114
tree_to_vertex = C_NULL
1115
1116
tree_to_tree = Array{p4est_topidx_t, 2}(undef, 4, n_trees)
1117
tree_to_face = Array{Int8, 2}(undef, 4, n_trees)
1118
1119
for cell_y in 1:n_cells_y, cell_x in 1:n_cells_x
1120
tree = linear_indices[cell_x, cell_y]
1121
1122
# Subtract 1 because `p4est` uses zero-based indexing
1123
# Negative x-direction
1124
if cell_x > 1
1125
tree_to_tree[1, tree] = linear_indices[cell_x - 1, cell_y] - 1
1126
tree_to_face[1, tree] = 1
1127
elseif periodicity[1]
1128
tree_to_tree[1, tree] = linear_indices[n_cells_x, cell_y] - 1
1129
tree_to_face[1, tree] = 1
1130
else # Non-periodic boundary, tree and face point to themselves (zero-based indexing)
1131
tree_to_tree[1, tree] = tree - 1
1132
tree_to_face[1, tree] = 0
1133
end
1134
1135
# Positive x-direction
1136
if cell_x < n_cells_x
1137
tree_to_tree[2, tree] = linear_indices[cell_x + 1, cell_y] - 1
1138
tree_to_face[2, tree] = 0
1139
elseif periodicity[1]
1140
tree_to_tree[2, tree] = linear_indices[1, cell_y] - 1
1141
tree_to_face[2, tree] = 0
1142
else # Non-periodic boundary, tree and face point to themselves (zero-based indexing)
1143
tree_to_tree[2, tree] = tree - 1
1144
tree_to_face[2, tree] = 1
1145
end
1146
1147
# Negative y-direction
1148
if cell_y > 1
1149
tree_to_tree[3, tree] = linear_indices[cell_x, cell_y - 1] - 1
1150
tree_to_face[3, tree] = 3
1151
elseif periodicity[2]
1152
tree_to_tree[3, tree] = linear_indices[cell_x, n_cells_y] - 1
1153
tree_to_face[3, tree] = 3
1154
else # Non-periodic boundary, tree and face point to themselves (zero-based indexing)
1155
tree_to_tree[3, tree] = tree - 1
1156
tree_to_face[3, tree] = 2
1157
end
1158
1159
# Positive y-direction
1160
if cell_y < n_cells_y
1161
tree_to_tree[4, tree] = linear_indices[cell_x, cell_y + 1] - 1
1162
tree_to_face[4, tree] = 2
1163
elseif periodicity[2]
1164
tree_to_tree[4, tree] = linear_indices[cell_x, 1] - 1
1165
tree_to_face[4, tree] = 2
1166
else # Non-periodic boundary, tree and face point to themselves (zero-based indexing)
1167
tree_to_tree[4, tree] = tree - 1
1168
tree_to_face[4, tree] = 3
1169
end
1170
end
1171
1172
tree_to_corner = C_NULL
1173
# `p4est` docs: "in trivial cases it is just a pointer to a p4est_topix value of 0."
1174
# We don't need corner connectivity, so this is a trivial case.
1175
ctt_offset = zeros(p4est_topidx_t, 1)
1176
1177
corner_to_tree = C_NULL
1178
corner_to_corner = C_NULL
1179
1180
connectivity = p4est_connectivity_new_copy(n_vertices, n_trees, n_corners,
1181
vertices, tree_to_vertex,
1182
tree_to_tree, tree_to_face,
1183
tree_to_corner, ctt_offset,
1184
corner_to_tree, corner_to_corner)
1185
1186
@assert p4est_connectivity_is_valid(connectivity) == 1
1187
1188
return connectivity
1189
end
1190
1191
# 3D version
1192
function connectivity_structured(n_cells_x, n_cells_y, n_cells_z, periodicity)
1193
linear_indices = LinearIndices((n_cells_x, n_cells_y, n_cells_z))
1194
1195
# Vertices represent the coordinates of the forest. This is used by `p4est`
1196
# to write VTK files.
1197
# Trixi.jl doesn't use the coordinates from `p4est`, so the vertices can be empty.
1198
n_vertices = 0
1199
n_trees = n_cells_x * n_cells_y * n_cells_z
1200
# No edge connectivity is needed
1201
n_edges = 0
1202
# No corner connectivity is needed
1203
n_corners = 0
1204
vertices = C_NULL
1205
tree_to_vertex = C_NULL
1206
1207
tree_to_tree = Array{p4est_topidx_t, 2}(undef, 6, n_trees)
1208
tree_to_face = Array{Int8, 2}(undef, 6, n_trees)
1209
1210
for cell_z in 1:n_cells_z, cell_y in 1:n_cells_y, cell_x in 1:n_cells_x
1211
tree = linear_indices[cell_x, cell_y, cell_z]
1212
1213
# Subtract 1 because `p4est` uses zero-based indexing
1214
# Negative x-direction
1215
if cell_x > 1
1216
tree_to_tree[1, tree] = linear_indices[cell_x - 1, cell_y, cell_z] - 1
1217
tree_to_face[1, tree] = 1
1218
elseif periodicity[1]
1219
tree_to_tree[1, tree] = linear_indices[n_cells_x, cell_y, cell_z] - 1
1220
tree_to_face[1, tree] = 1
1221
else # Non-periodic boundary, tree and face point to themselves (zero-based indexing)
1222
tree_to_tree[1, tree] = tree - 1
1223
tree_to_face[1, tree] = 0
1224
end
1225
1226
# Positive x-direction
1227
if cell_x < n_cells_x
1228
tree_to_tree[2, tree] = linear_indices[cell_x + 1, cell_y, cell_z] - 1
1229
tree_to_face[2, tree] = 0
1230
elseif periodicity[1]
1231
tree_to_tree[2, tree] = linear_indices[1, cell_y, cell_z] - 1
1232
tree_to_face[2, tree] = 0
1233
else # Non-periodic boundary, tree and face point to themselves (zero-based indexing)
1234
tree_to_tree[2, tree] = tree - 1
1235
tree_to_face[2, tree] = 1
1236
end
1237
1238
# Negative y-direction
1239
if cell_y > 1
1240
tree_to_tree[3, tree] = linear_indices[cell_x, cell_y - 1, cell_z] - 1
1241
tree_to_face[3, tree] = 3
1242
elseif periodicity[2]
1243
tree_to_tree[3, tree] = linear_indices[cell_x, n_cells_y, cell_z] - 1
1244
tree_to_face[3, tree] = 3
1245
else # Non-periodic boundary, tree and face point to themselves (zero-based indexing)
1246
tree_to_tree[3, tree] = tree - 1
1247
tree_to_face[3, tree] = 2
1248
end
1249
1250
# Positive y-direction
1251
if cell_y < n_cells_y
1252
tree_to_tree[4, tree] = linear_indices[cell_x, cell_y + 1, cell_z] - 1
1253
tree_to_face[4, tree] = 2
1254
elseif periodicity[2]
1255
tree_to_tree[4, tree] = linear_indices[cell_x, 1, cell_z] - 1
1256
tree_to_face[4, tree] = 2
1257
else # Non-periodic boundary, tree and face point to themselves (zero-based indexing)
1258
tree_to_tree[4, tree] = tree - 1
1259
tree_to_face[4, tree] = 3
1260
end
1261
1262
# Negative z-direction
1263
if cell_z > 1
1264
tree_to_tree[5, tree] = linear_indices[cell_x, cell_y, cell_z - 1] - 1
1265
tree_to_face[5, tree] = 5
1266
elseif periodicity[3]
1267
tree_to_tree[5, tree] = linear_indices[cell_x, cell_y, n_cells_z] - 1
1268
tree_to_face[5, tree] = 5
1269
else # Non-periodic boundary, tree and face point to themselves (zero-based indexing)
1270
tree_to_tree[5, tree] = tree - 1
1271
tree_to_face[5, tree] = 4
1272
end
1273
1274
# Positive z-direction
1275
if cell_z < n_cells_z
1276
tree_to_tree[6, tree] = linear_indices[cell_x, cell_y, cell_z + 1] - 1
1277
tree_to_face[6, tree] = 4
1278
elseif periodicity[3]
1279
tree_to_tree[6, tree] = linear_indices[cell_x, cell_y, 1] - 1
1280
tree_to_face[6, tree] = 4
1281
else # Non-periodic boundary, tree and face point to themselves (zero-based indexing)
1282
tree_to_tree[6, tree] = tree - 1
1283
tree_to_face[6, tree] = 5
1284
end
1285
end
1286
1287
tree_to_edge = C_NULL
1288
# `p4est` docs: "in trivial cases it is just a pointer to a p4est_topix value of 0."
1289
# We don't need edge connectivity, so this is a trivial case.
1290
ett_offset = zeros(p4est_topidx_t, 1)
1291
edge_to_tree = C_NULL
1292
edge_to_edge = C_NULL
1293
1294
tree_to_corner = C_NULL
1295
# `p4est` docs: "in trivial cases it is just a pointer to a p4est_topix value of 0."
1296
# We don't need corner connectivity, so this is a trivial case.
1297
ctt_offset = zeros(p4est_topidx_t, 1)
1298
1299
corner_to_tree = C_NULL
1300
corner_to_corner = C_NULL
1301
1302
connectivity = p8est_connectivity_new_copy(n_vertices, n_trees, n_corners, n_edges,
1303
vertices, tree_to_vertex,
1304
tree_to_tree, tree_to_face,
1305
tree_to_edge, ett_offset,
1306
edge_to_tree, edge_to_edge,
1307
tree_to_corner, ctt_offset,
1308
corner_to_tree, corner_to_corner)
1309
1310
@assert p8est_connectivity_is_valid(connectivity) == 1
1311
1312
return connectivity
1313
end
1314
1315
function connectivity_cubed_sphere(trees_per_face_dimension, layers)
1316
n_cells_x = n_cells_y = trees_per_face_dimension
1317
n_cells_z = layers
1318
1319
linear_indices = LinearIndices((trees_per_face_dimension, trees_per_face_dimension,
1320
layers, 6))
1321
1322
# Vertices represent the coordinates of the forest. This is used by `p4est`
1323
# to write VTK files.
1324
# Trixi.jl doesn't use the coordinates from `p4est`, so the vertices can be empty.
1325
n_vertices = 0
1326
n_trees = 6 * n_cells_x * n_cells_y * n_cells_z
1327
# No edge connectivity is needed
1328
n_edges = 0
1329
# No corner connectivity is needed
1330
n_corners = 0
1331
vertices = C_NULL
1332
tree_to_vertex = C_NULL
1333
1334
tree_to_tree = Array{p4est_topidx_t, 2}(undef, 6, n_trees)
1335
tree_to_face = Array{Int8, 2}(undef, 6, n_trees)
1336
1337
# Illustration of the local coordinates of each face. ξ and η are the first
1338
# local coordinates of each face. The third local coordinate ζ is always
1339
# pointing outwards, which yields a right-handed coordinate system for each face.
1340
# ┌────────────────────────────────────────────────────┐
1341
# ╱│ ╱│
1342
# ╱ │ ξ <───┐ ╱ │
1343
# ╱ │ ╱ ╱ │
1344
# ╱ │ 4 (+y) V ╱ │
1345
# ╱ │ η ╱ │
1346
# ╱ │ ╱ │
1347
# ╱ │ ╱ │
1348
# ╱ │ ╱ │
1349
# ╱ │ ╱ │
1350
# ╱ │ 5 (-z) η ╱ │
1351
# ╱ │ ↑ ╱ │
1352
# ╱ │ │ ╱ │
1353
# ╱ │ ξ <───┘ ╱ │
1354
# ┌────────────────────────────────────────────────────┐ 2 (+x) │
1355
# │ │ │ │
1356
# │ │ │ ξ │
1357
# │ │ │ ↑ │
1358
# │ 1 (-x) │ │ │ │
1359
# │ │ │ │ │
1360
# │ ╱│ │ │ ╱ │
1361
# │ V │ │ │ V │
1362
# │ η ↓ │ │ η │
1363
# │ ξ └──────────────────────────────────────│─────────────┘
1364
# │ ╱ η 6 (+z) │ ╱
1365
# │ ╱ ↑ │ ╱
1366
# │ ╱ │ │ ╱
1367
# │ ╱ └───> ξ │ ╱
1368
# │ ╱ │ ╱
1369
# │ ╱ │ ╱ Global coordinates:
1370
# │ ╱ │ ╱ y
1371
# │ ╱ ┌───> ξ │ ╱ ↑
1372
# │ ╱ ╱ │ ╱ │
1373
# │ ╱ V 3 (-y) │ ╱ │
1374
# │ ╱ η │ ╱ └─────> x
1375
# │ ╱ │ ╱ ╱
1376
# │╱ │╱ V
1377
# └────────────────────────────────────────────────────┘ z
1378
for direction in 1:6
1379
for cell_z in 1:n_cells_z, cell_y in 1:n_cells_y, cell_x in 1:n_cells_x
1380
tree = linear_indices[cell_x, cell_y, cell_z, direction]
1381
1382
# Subtract 1 because `p4est` uses zero-based indexing
1383
# Negative x-direction
1384
if cell_x > 1 # Connect to tree at the same face
1385
tree_to_tree[1, tree] = linear_indices[cell_x - 1, cell_y, cell_z,
1386
direction] - 1
1387
tree_to_face[1, tree] = 1
1388
elseif direction == 1 # This is the -x face
1389
target = 4
1390
tree_to_tree[1, tree] = linear_indices[end, cell_y, cell_z, target] - 1
1391
tree_to_face[1, tree] = 1
1392
elseif direction == 2 # This is the +x face
1393
target = 3
1394
tree_to_tree[1, tree] = linear_indices[end, cell_y, cell_z, target] - 1
1395
tree_to_face[1, tree] = 1
1396
elseif direction == 3 # This is the -y face
1397
target = 1
1398
tree_to_tree[1, tree] = linear_indices[end, cell_y, cell_z, target] - 1
1399
tree_to_face[1, tree] = 1
1400
elseif direction == 4 # This is the +y face
1401
target = 2
1402
tree_to_tree[1, tree] = linear_indices[end, cell_y, cell_z, target] - 1
1403
tree_to_face[1, tree] = 1
1404
elseif direction == 5 # This is the -z face
1405
target = 2
1406
tree_to_tree[1, tree] = linear_indices[cell_y, 1, cell_z, target] - 1
1407
tree_to_face[1, tree] = 2
1408
else # direction == 6, this is the +z face
1409
target = 1
1410
tree_to_tree[1, tree] = linear_indices[end - cell_y + 1, end, cell_z,
1411
target] - 1
1412
tree_to_face[1, tree] = 9 # first face dimensions are oppositely oriented, add 6
1413
end
1414
1415
# Positive x-direction
1416
if cell_x < n_cells_x # Connect to tree at the same face
1417
tree_to_tree[2, tree] = linear_indices[cell_x + 1, cell_y, cell_z,
1418
direction] - 1
1419
tree_to_face[2, tree] = 0
1420
elseif direction == 1 # This is the -x face
1421
target = 3
1422
tree_to_tree[2, tree] = linear_indices[1, cell_y, cell_z, target] - 1
1423
tree_to_face[2, tree] = 0
1424
elseif direction == 2 # This is the +x face
1425
target = 4
1426
tree_to_tree[2, tree] = linear_indices[1, cell_y, cell_z, target] - 1
1427
tree_to_face[2, tree] = 0
1428
elseif direction == 3 # This is the -y face
1429
target = 2
1430
tree_to_tree[2, tree] = linear_indices[1, cell_y, cell_z, target] - 1
1431
tree_to_face[2, tree] = 0
1432
elseif direction == 4 # This is the +y face
1433
target = 1
1434
tree_to_tree[2, tree] = linear_indices[1, cell_y, cell_z, target] - 1
1435
tree_to_face[2, tree] = 0
1436
elseif direction == 5 # This is the -z face
1437
target = 1
1438
tree_to_tree[2, tree] = linear_indices[end - cell_y + 1, 1, cell_z,
1439
target] - 1
1440
tree_to_face[2, tree] = 8 # first face dimensions are oppositely oriented, add 6
1441
else # direction == 6, this is the +z face
1442
target = 2
1443
tree_to_tree[2, tree] = linear_indices[cell_y, end, cell_z, target] - 1
1444
tree_to_face[2, tree] = 3
1445
end
1446
1447
# Negative y-direction
1448
if cell_y > 1 # Connect to tree at the same face
1449
tree_to_tree[3, tree] = linear_indices[cell_x, cell_y - 1, cell_z,
1450
direction] - 1
1451
tree_to_face[3, tree] = 3
1452
elseif direction == 1
1453
target = 5
1454
tree_to_tree[3, tree] = linear_indices[end, end - cell_x + 1, cell_z,
1455
target] - 1
1456
tree_to_face[3, tree] = 7 # first face dimensions are oppositely oriented, add 6
1457
elseif direction == 2
1458
target = 5
1459
tree_to_tree[3, tree] = linear_indices[1, cell_x, cell_z, target] - 1
1460
tree_to_face[3, tree] = 0
1461
elseif direction == 3
1462
target = 5
1463
tree_to_tree[3, tree] = linear_indices[end - cell_x + 1, 1, cell_z,
1464
target] - 1
1465
tree_to_face[3, tree] = 8 # first face dimensions are oppositely oriented, add 6
1466
elseif direction == 4
1467
target = 5
1468
tree_to_tree[3, tree] = linear_indices[cell_x, end, cell_z, target] - 1
1469
tree_to_face[3, tree] = 3
1470
elseif direction == 5
1471
target = 3
1472
tree_to_tree[3, tree] = linear_indices[end - cell_x + 1, 1, cell_z,
1473
target] - 1
1474
tree_to_face[3, tree] = 8 # first face dimensions are oppositely oriented, add 6
1475
else # direction == 6
1476
target = 3
1477
tree_to_tree[3, tree] = linear_indices[cell_x, end, cell_z, target] - 1
1478
tree_to_face[3, tree] = 3
1479
end
1480
1481
# Positive y-direction
1482
if cell_y < n_cells_y # Connect to tree at the same face
1483
tree_to_tree[4, tree] = linear_indices[cell_x, cell_y + 1, cell_z,
1484
direction] - 1
1485
tree_to_face[4, tree] = 2
1486
elseif direction == 1
1487
target = 6
1488
tree_to_tree[4, tree] = linear_indices[1, end - cell_x + 1, cell_z,
1489
target] - 1
1490
tree_to_face[4, tree] = 6 # first face dimensions are oppositely oriented, add 6
1491
elseif direction == 2
1492
target = 6
1493
tree_to_tree[4, tree] = linear_indices[end, cell_x, cell_z, target] - 1
1494
tree_to_face[4, tree] = 1
1495
elseif direction == 3
1496
target = 6
1497
tree_to_tree[4, tree] = linear_indices[cell_x, 1, cell_z, target] - 1
1498
tree_to_face[4, tree] = 2
1499
elseif direction == 4
1500
target = 6
1501
tree_to_tree[4, tree] = linear_indices[end - cell_x + 1, end, cell_z,
1502
target] - 1
1503
tree_to_face[4, tree] = 9 # first face dimensions are oppositely oriented, add 6
1504
elseif direction == 5
1505
target = 4
1506
tree_to_tree[4, tree] = linear_indices[cell_x, 1, cell_z, target] - 1
1507
tree_to_face[4, tree] = 2
1508
else # direction == 6
1509
target = 4
1510
tree_to_tree[4, tree] = linear_indices[end - cell_x + 1, end, cell_z,
1511
target] - 1
1512
tree_to_face[4, tree] = 9 # first face dimensions are oppositely oriented, add 6
1513
end
1514
1515
# Negative z-direction
1516
if cell_z > 1
1517
tree_to_tree[5, tree] = linear_indices[cell_x, cell_y, cell_z - 1,
1518
direction] - 1
1519
tree_to_face[5, tree] = 5
1520
else # Non-periodic boundary, tree and face point to themselves (zero-based indexing)
1521
tree_to_tree[5, tree] = tree - 1
1522
tree_to_face[5, tree] = 4
1523
end
1524
1525
# Positive z-direction
1526
if cell_z < n_cells_z
1527
tree_to_tree[6, tree] = linear_indices[cell_x, cell_y, cell_z + 1,
1528
direction] - 1
1529
tree_to_face[6, tree] = 4
1530
else # Non-periodic boundary, tree and face point to themselves (zero-based indexing)
1531
tree_to_tree[6, tree] = tree - 1
1532
tree_to_face[6, tree] = 5
1533
end
1534
end
1535
end
1536
1537
tree_to_edge = C_NULL
1538
# `p4est` docs: "in trivial cases it is just a pointer to a p4est_topix value of 0."
1539
# We don't need edge connectivity, so this is a trivial case.
1540
ett_offset = zeros(p4est_topidx_t, 1)
1541
edge_to_tree = C_NULL
1542
edge_to_edge = C_NULL
1543
1544
tree_to_corner = C_NULL
1545
# `p4est` docs: "in trivial cases it is just a pointer to a p4est_topix value of 0."
1546
# We don't need corner connectivity, so this is a trivial case.
1547
ctt_offset = zeros(p4est_topidx_t, 1)
1548
1549
corner_to_tree = C_NULL
1550
corner_to_corner = C_NULL
1551
1552
connectivity = p8est_connectivity_new_copy(n_vertices, n_trees, n_corners, n_edges,
1553
vertices, tree_to_vertex,
1554
tree_to_tree, tree_to_face,
1555
tree_to_edge, ett_offset,
1556
edge_to_tree, edge_to_edge,
1557
tree_to_corner, ctt_offset,
1558
corner_to_tree, corner_to_corner)
1559
1560
@assert p8est_connectivity_is_valid(connectivity) == 1
1561
1562
return connectivity
1563
end
1564
1565
# Calculate physical coordinates of each node of a structured mesh.
1566
# This function assumes a structured mesh with trees in row order.
1567
# 2D version
1568
function calc_tree_node_coordinates!(node_coordinates::AbstractArray{<:Any, 4},
1569
nodes, mapping, trees_per_dimension)
1570
linear_indices = LinearIndices(trees_per_dimension)
1571
1572
# Get cell length in reference mesh
1573
dx = 2 / trees_per_dimension[1]
1574
dy = 2 / trees_per_dimension[2]
1575
1576
for cell_y in 1:trees_per_dimension[2], cell_x in 1:trees_per_dimension[1]
1577
tree_id = linear_indices[cell_x, cell_y]
1578
1579
# Calculate node coordinates of reference mesh
1580
cell_x_offset = -1 + (cell_x - 1) * dx + dx / 2
1581
cell_y_offset = -1 + (cell_y - 1) * dy + dy / 2
1582
1583
for j in eachindex(nodes), i in eachindex(nodes)
1584
# node_coordinates are the mapped reference node coordinates
1585
node_coordinates[:, i, j, tree_id] .= mapping(cell_x_offset +
1586
dx / 2 * nodes[i],
1587
cell_y_offset +
1588
dy / 2 * nodes[j])
1589
end
1590
end
1591
1592
return nothing
1593
end
1594
1595
# 3D version
1596
function calc_tree_node_coordinates!(node_coordinates::AbstractArray{<:Any, 5},
1597
nodes, mapping, trees_per_dimension)
1598
linear_indices = LinearIndices(trees_per_dimension)
1599
1600
# Get cell length in reference mesh
1601
dx = 2 / trees_per_dimension[1]
1602
dy = 2 / trees_per_dimension[2]
1603
dz = 2 / trees_per_dimension[3]
1604
1605
for cell_z in 1:trees_per_dimension[3],
1606
cell_y in 1:trees_per_dimension[2],
1607
cell_x in 1:trees_per_dimension[1]
1608
1609
tree_id = linear_indices[cell_x, cell_y, cell_z]
1610
1611
# Calculate node coordinates of reference mesh
1612
cell_x_offset = -1 + (cell_x - 1) * dx + dx / 2
1613
cell_y_offset = -1 + (cell_y - 1) * dy + dy / 2
1614
cell_z_offset = -1 + (cell_z - 1) * dz + dz / 2
1615
1616
for k in eachindex(nodes), j in eachindex(nodes), i in eachindex(nodes)
1617
# node_coordinates are the mapped reference node coordinates
1618
node_coordinates[:, i, j, k, tree_id] .= mapping(cell_x_offset +
1619
dx / 2 * nodes[i],
1620
cell_y_offset +
1621
dy / 2 * nodes[j],
1622
cell_z_offset +
1623
dz / 2 * nodes[k])
1624
end
1625
end
1626
1627
return nothing
1628
end
1629
1630
# Calculate physical coordinates of each node of an unstructured mesh.
1631
# Extract corners of each tree from the connectivity,
1632
# interpolate to requested interpolation nodes,
1633
# map the resulting coordinates with the specified mapping.
1634
# 2D version
1635
function calc_tree_node_coordinates!(node_coordinates::AbstractArray{RealT, 4},
1636
nodes, mapping,
1637
vertices, tree_to_vertex) where {RealT}
1638
nodes_in = [-1.0, 1.0]
1639
matrix = polynomial_interpolation_matrix(nodes_in, nodes)
1640
data_in = Array{RealT, 3}(undef, 2, 2, 2)
1641
tmp1 = zeros(RealT, 2, length(nodes), length(nodes_in))
1642
1643
for tree in 1:size(tree_to_vertex, 2)
1644
# Tree vertices are stored in Z-order, ignore z-coordinate in 2D, zero-based indexing
1645
@views data_in[:, 1, 1] .= vertices[1:2, tree_to_vertex[1, tree] + 1]
1646
@views data_in[:, 2, 1] .= vertices[1:2, tree_to_vertex[2, tree] + 1]
1647
@views data_in[:, 1, 2] .= vertices[1:2, tree_to_vertex[3, tree] + 1]
1648
@views data_in[:, 2, 2] .= vertices[1:2, tree_to_vertex[4, tree] + 1]
1649
1650
# Interpolate corner coordinates to specified nodes
1651
multiply_dimensionwise!(view(node_coordinates, :, :, :, tree),
1652
matrix, matrix,
1653
data_in,
1654
tmp1)
1655
end
1656
1657
return map_node_coordinates!(node_coordinates, mapping)
1658
end
1659
1660
function map_node_coordinates!(node_coordinates::AbstractArray{<:Any, 4}, mapping)
1661
for tree in axes(node_coordinates, 4),
1662
j in axes(node_coordinates, 3),
1663
i in axes(node_coordinates, 2)
1664
1665
node_coordinates[:, i, j, tree] .= mapping(node_coordinates[1, i, j, tree],
1666
node_coordinates[2, i, j, tree])
1667
end
1668
1669
return node_coordinates
1670
end
1671
1672
function map_node_coordinates!(node_coordinates::AbstractArray{<:Any, 4},
1673
mapping::Nothing)
1674
return node_coordinates
1675
end
1676
1677
# 3D version
1678
function calc_tree_node_coordinates!(node_coordinates::AbstractArray{RealT, 5},
1679
nodes, mapping,
1680
vertices, tree_to_vertex) where {RealT}
1681
nodes_in = [-1.0, 1.0]
1682
matrix = polynomial_interpolation_matrix(nodes_in, nodes)
1683
data_in = Array{RealT, 4}(undef, 3, 2, 2, 2)
1684
1685
for tree in 1:size(tree_to_vertex, 2)
1686
# Tree vertices are stored in Z-order, zero-based indexing
1687
@views data_in[:, 1, 1, 1] .= vertices[:, tree_to_vertex[1, tree] + 1]
1688
@views data_in[:, 2, 1, 1] .= vertices[:, tree_to_vertex[2, tree] + 1]
1689
@views data_in[:, 1, 2, 1] .= vertices[:, tree_to_vertex[3, tree] + 1]
1690
@views data_in[:, 2, 2, 1] .= vertices[:, tree_to_vertex[4, tree] + 1]
1691
@views data_in[:, 1, 1, 2] .= vertices[:, tree_to_vertex[5, tree] + 1]
1692
@views data_in[:, 2, 1, 2] .= vertices[:, tree_to_vertex[6, tree] + 1]
1693
@views data_in[:, 1, 2, 2] .= vertices[:, tree_to_vertex[7, tree] + 1]
1694
@views data_in[:, 2, 2, 2] .= vertices[:, tree_to_vertex[8, tree] + 1]
1695
1696
# Interpolate corner coordinates to specified nodes
1697
multiply_dimensionwise!(view(node_coordinates, :, :, :, :, tree),
1698
matrix, matrix, matrix,
1699
data_in)
1700
end
1701
1702
return map_node_coordinates!(node_coordinates, mapping)
1703
end
1704
1705
function map_node_coordinates!(node_coordinates::AbstractArray{<:Any, 5}, mapping)
1706
for tree in axes(node_coordinates, 5),
1707
k in axes(node_coordinates, 4),
1708
j in axes(node_coordinates, 3),
1709
i in axes(node_coordinates, 2)
1710
1711
node_coordinates[:, i, j, k, tree] .= mapping(node_coordinates[1, i, j, k,
1712
tree],
1713
node_coordinates[2, i, j, k,
1714
tree],
1715
node_coordinates[3, i, j, k,
1716
tree])
1717
end
1718
1719
return node_coordinates
1720
end
1721
1722
function map_node_coordinates!(node_coordinates::AbstractArray{<:Any, 5},
1723
mapping::Nothing)
1724
return node_coordinates
1725
end
1726
1727
# Calculate physical coordinates of each node of a cubed sphere mesh.
1728
function calc_tree_node_coordinates!(node_coordinates::AbstractArray{<:Any, 5},
1729
nodes, trees_per_face_dimension, layers,
1730
inner_radius, thickness)
1731
n_cells_x = n_cells_y = trees_per_face_dimension
1732
n_cells_z = layers
1733
1734
linear_indices = LinearIndices((n_cells_x, n_cells_y, n_cells_z, 6))
1735
1736
# Get cell length in reference mesh
1737
dx = 2 / n_cells_x
1738
dy = 2 / n_cells_y
1739
dz = 2 / n_cells_z
1740
1741
for direction in 1:6
1742
for cell_z in 1:n_cells_z, cell_y in 1:n_cells_y, cell_x in 1:n_cells_x
1743
tree = linear_indices[cell_x, cell_y, cell_z, direction]
1744
1745
x_offset = -1 + (cell_x - 1) * dx + dx / 2
1746
y_offset = -1 + (cell_y - 1) * dy + dy / 2
1747
z_offset = -1 + (cell_z - 1) * dz + dz / 2
1748
1749
for k in eachindex(nodes), j in eachindex(nodes), i in eachindex(nodes)
1750
# node_coordinates are the mapped reference node coordinates
1751
node_coordinates[:, i, j, k, tree] .= cubed_sphere_mapping(x_offset +
1752
dx / 2 *
1753
nodes[i],
1754
y_offset +
1755
dy / 2 *
1756
nodes[j],
1757
z_offset +
1758
dz / 2 *
1759
nodes[k],
1760
inner_radius,
1761
thickness,
1762
direction)
1763
end
1764
end
1765
end
1766
1767
return nothing
1768
end
1769
1770
# Map the computational coordinates xi, eta, zeta to the specified side of a cubed sphere
1771
# with the specified inner radius and thickness.
1772
function cubed_sphere_mapping(xi, eta, zeta, inner_radius, thickness, direction)
1773
alpha = xi * pi / 4
1774
beta = eta * pi / 4
1775
1776
# Equiangular projection
1777
x = tan(alpha)
1778
y = tan(beta)
1779
1780
# Coordinates on unit cube per direction, see illustration above in the function connectivity_cubed_sphere
1781
cube_coordinates = (SVector(-1, -x, y),
1782
SVector(1, x, y),
1783
SVector(x, -1, y),
1784
SVector(-x, 1, y),
1785
SVector(-x, y, -1),
1786
SVector(x, y, 1))
1787
1788
# Radius on cube surface
1789
r = sqrt(1 + x^2 + y^2)
1790
1791
# Radius of the sphere
1792
R = inner_radius + thickness * (0.5f0 * (zeta + 1))
1793
1794
# Projection onto the sphere
1795
return R / r * cube_coordinates[direction]
1796
end
1797
1798
# Calculate physical coordinates of each element of an unstructured mesh read
1799
# in from a HOHQMesh file. This calculation is done with the transfinite interpolation
1800
# routines found in `mappings_geometry_curved_2d.jl` or `mappings_geometry_straight_2d.jl`
1801
function calc_tree_node_coordinates!(node_coordinates::AbstractArray{<:Any, 4},
1802
file_lines::Vector{String}, nodes, vertices, RealT)
1803
# Get the number of trees (elements) and the number of interpolation nodes
1804
n_trees = last(size(node_coordinates))
1805
nnodes = length(nodes)
1806
1807
# Setup the starting file index to read in element indices and the additional
1808
# higher-order (curved) boundary information provided by HOHQMesh.
1809
file_idx = findfirst(contains("** mesh polynomial degree"), file_lines) + 1
1810
1811
# Create a work set of Gamma curves to create the node coordinates
1812
CurvedSurfaceT = CurvedSurface{RealT}
1813
surface_curves = Array{CurvedSurfaceT}(undef, 4)
1814
1815
# Create other work arrays to perform the mesh construction
1816
element_node_ids = Array{Int}(undef, 4)
1817
curved_check = Vector{Int}(undef, 4)
1818
quad_vertices = Array{RealT}(undef, (4, 2))
1819
quad_vertices_flipped = Array{RealT}(undef, (4, 2))
1820
curve_values = Array{RealT}(undef, (nnodes, 2))
1821
1822
# Create the barycentric weights used for the surface interpolations
1823
bary_weights_ = barycentric_weights(nodes)
1824
bary_weights = SVector{nnodes}(bary_weights_)
1825
1826
# Loop through all the trees, i.e., the elements generated by HOHQMesh and create the node coordinates.
1827
# When we extract information from the `current_line` we start at index 2 in order to
1828
# avoid the Abaqus comment character "** "
1829
for tree in 1:n_trees
1830
# Pull the vertex node IDs
1831
current_line = split(file_lines[file_idx])
1832
element_node_ids[1] = parse(Int, current_line[2])
1833
element_node_ids[2] = parse(Int, current_line[3])
1834
element_node_ids[3] = parse(Int, current_line[4])
1835
element_node_ids[4] = parse(Int, current_line[5])
1836
1837
# Pull the (x,y) values of the four vertices of the current tree out of the global vertices array
1838
for i in 1:4
1839
quad_vertices[i, :] .= vertices[1:2, element_node_ids[i]] # 2D => 1:2
1840
end
1841
# Pull the information to check if boundary is curved/higher-order in order to read in additional data
1842
file_idx += 1
1843
current_line = split(file_lines[file_idx])
1844
# Note: This strategy is HOHQMesh-Abaqus.inp specific!
1845
curved_check[1] = parse(Int, current_line[2])
1846
curved_check[2] = parse(Int, current_line[3])
1847
curved_check[3] = parse(Int, current_line[4])
1848
curved_check[4] = parse(Int, current_line[5])
1849
if sum(curved_check) == 0
1850
# Create the node coordinates on this particular element
1851
calc_node_coordinates!(node_coordinates, tree, nodes, quad_vertices)
1852
else
1853
# Quadrilateral element has at least one curved/higher-order side
1854
# Flip node ordering to make sure the element is right-handed for the interpolations
1855
m1 = 1
1856
m2 = 2
1857
@views quad_vertices_flipped[1, :] .= quad_vertices[4, :]
1858
@views quad_vertices_flipped[2, :] .= quad_vertices[2, :]
1859
@views quad_vertices_flipped[3, :] .= quad_vertices[3, :]
1860
@views quad_vertices_flipped[4, :] .= quad_vertices[1, :]
1861
for i in 1:4
1862
if curved_check[i] == 0
1863
# When curved_check[i] is 0 then the "curve" from vertex `i` to vertex `i+1` is a straight line.
1864
# Evaluate a linear interpolant between the two points at each of the nodes.
1865
for k in 1:nnodes
1866
curve_values[k, 1] = linear_interpolate(nodes[k],
1867
quad_vertices_flipped[m1,
1868
1],
1869
quad_vertices_flipped[m2,
1870
1])
1871
curve_values[k, 2] = linear_interpolate(nodes[k],
1872
quad_vertices_flipped[m1,
1873
2],
1874
quad_vertices_flipped[m2,
1875
2])
1876
end
1877
else
1878
# When curved_check[i] is 1 this curved/higher-order boundary information is supplied by the mesh
1879
# generator. So we just read it into a work array
1880
for k in 1:nnodes
1881
file_idx += 1
1882
current_line = split(file_lines[file_idx])
1883
curve_values[k, 1] = parse(RealT, current_line[2])
1884
curve_values[k, 2] = parse(RealT, current_line[3])
1885
end
1886
end
1887
# Construct the curve interpolant for the current side
1888
surface_curves[i] = CurvedSurfaceT(nodes, bary_weights,
1889
copy(curve_values))
1890
# Indexing update that contains a "flip" to ensure correct element orientation.
1891
# If we need to construct the straight line "curves" when curved_check[i] == 0
1892
m1 += 1
1893
if i == 3
1894
m2 = 1
1895
else
1896
m2 += 1
1897
end
1898
end
1899
# Create the node coordinates on this particular element
1900
calc_node_coordinates!(node_coordinates, tree, nodes, surface_curves)
1901
end
1902
# Move file index to the next tree
1903
file_idx += 1
1904
end
1905
1906
return file_idx
1907
end
1908
1909
# Version for quadratic 2D elements, i.e., second-order quads in
1910
# standard Abaqus .inp format as exported from e.g. Gmsh.
1911
function calc_tree_node_coordinates!(node_coordinates::AbstractArray{<:Any, 4},
1912
element_lines, nodes, vertices, RealT,
1913
linear_quads, mesh_nodes)
1914
nnodes = length(nodes)
1915
1916
# Create a work set of Gamma curves to create the node coordinates
1917
CurvedSurfaceT = CurvedSurface{RealT}
1918
surface_curves = Array{CurvedSurfaceT}(undef, 4)
1919
1920
# Create other work arrays to perform the mesh construction
1921
quad_vertices = Array{RealT}(undef, (4, 2))
1922
curve_values = Array{RealT}(undef, (nnodes, 2))
1923
1924
# Create the barycentric weights used for the surface interpolations
1925
bary_weights_ = barycentric_weights(nodes)
1926
bary_weights = SVector{nnodes}(bary_weights_)
1927
1928
element_set_order = 1
1929
tree = 0
1930
for line_idx in 1:length(element_lines)
1931
line = element_lines[line_idx]
1932
1933
# Check if a new element type/element set is defined
1934
if startswith(line, "*ELEMENT")
1935
# Retrieve element type
1936
current_element_type = match(r"\*ELEMENT, type=([^,]+)", line).captures[1]
1937
1938
# Check if these are linear elements
1939
if occursin(linear_quads, current_element_type)
1940
element_set_order = 1
1941
else
1942
element_set_order = 2
1943
end
1944
else # Element data
1945
tree += 1
1946
1947
# Pull the vertex node IDs
1948
line_split = split(line, r",\s+")
1949
# Store the indices of the element-defining nodes
1950
element_nodes = parse.(Int, line_split)[2:end] # First entry is the element number, can be ignored
1951
1952
if element_set_order == 1
1953
# Create the node coordinates on this particular element
1954
# Pull the (x,y) values of the four vertices of the current tree out of the global vertices array
1955
for i in 1:4
1956
quad_vertices[i, :] .= vertices[1:2, element_nodes[i]] # 2D => 1:2
1957
end
1958
calc_node_coordinates!(node_coordinates, tree, nodes, quad_vertices)
1959
else # element_set_order == 2
1960
for edge in 1:4
1961
node1 = element_nodes[edge] # "Left" node
1962
node2 = element_nodes[edge + 4] # "Middle" node
1963
node3 = edge == 4 ? element_nodes[1] : element_nodes[edge + 1] # "Right" node
1964
1965
node1_coords = mesh_nodes[node1]
1966
node2_coords = mesh_nodes[node2]
1967
node3_coords = mesh_nodes[node3]
1968
1969
# The nodes for an Abaqus element are labeled following a closed path
1970
# around the element:
1971
#
1972
# <----
1973
# *-----*-----*
1974
# | |
1975
# | | | ^
1976
# | * * |
1977
# v | | |
1978
# | |
1979
# *-----*-----*
1980
# ---->
1981
# `curve_values`, however, requires to sort the nodes into a
1982
# valid coordinate system,
1983
#
1984
# *-----*-----*
1985
# | |
1986
# | |
1987
# * *
1988
# η | |
1989
# ↑ | |
1990
# | *-----*-----*
1991
# |----> ξ
1992
# thus we need to flip the node order for the second xi and eta edges met.
1993
1994
if edge in [1, 2]
1995
curve_values[1, 1] = node1_coords[1]
1996
curve_values[1, 2] = node1_coords[2]
1997
1998
curve_values[2, 1] = node2_coords[1]
1999
curve_values[2, 2] = node2_coords[2]
2000
2001
curve_values[3, 1] = node3_coords[1]
2002
curve_values[3, 2] = node3_coords[2]
2003
else # Flip "left" and "right" nodes
2004
curve_values[1, 1] = node3_coords[1]
2005
curve_values[1, 2] = node3_coords[2]
2006
2007
curve_values[2, 1] = node2_coords[1]
2008
curve_values[2, 2] = node2_coords[2]
2009
2010
curve_values[3, 1] = node1_coords[1]
2011
curve_values[3, 2] = node1_coords[2]
2012
end
2013
2014
# Construct the curve interpolant for the current side
2015
surface_curves[edge] = CurvedSurfaceT(nodes, bary_weights,
2016
copy(curve_values))
2017
end
2018
2019
# Create the node coordinates on this particular element
2020
calc_node_coordinates!(node_coordinates, tree, nodes, surface_curves)
2021
end
2022
end
2023
end
2024
2025
return nothing
2026
end
2027
2028
# Calculate physical coordinates of each element of an unstructured mesh read
2029
# in from a HOHQMesh file. This calculation is done with the transfinite interpolation
2030
# routines found in `transfinite_mappings_3d.jl`
2031
function calc_tree_node_coordinates!(node_coordinates::AbstractArray{<:Any, 5},
2032
file_lines::Vector{String}, nodes, vertices, RealT)
2033
# Get the number of trees and the number of interpolation nodes
2034
n_trees = last(size(node_coordinates))
2035
nnodes = length(nodes)
2036
2037
# Setup the starting file index to read in element indices and the additional
2038
# curved/higher-order boundary information provided by HOHQMesh.
2039
file_idx = findfirst(contains("** mesh polynomial degree"), file_lines) + 1
2040
2041
# Create a work set of Gamma curves to create the node coordinates
2042
CurvedFaceT = CurvedFace{RealT}
2043
face_curves = Array{CurvedFaceT}(undef, 6)
2044
2045
# Create other work arrays to perform the mesh construction
2046
element_node_ids = Array{Int}(undef, 8)
2047
curved_check = Vector{Int}(undef, 6)
2048
hex_vertices = Array{RealT}(undef, (3, 8))
2049
face_vertices = Array{RealT}(undef, (3, 4))
2050
curve_values = Array{RealT}(undef, (3, nnodes, nnodes))
2051
2052
# Create the barycentric weights used for the surface interpolations
2053
bary_weights_ = barycentric_weights(nodes)
2054
bary_weights = SVector{nnodes}(bary_weights_)
2055
2056
# Loop through all the trees, i.e., the elements generated by HOHQMesh and create the node coordinates.
2057
# When we extract information from the `current_line` we start at index 2 in order to
2058
# avoid the Abaqus comment character "** "
2059
for tree in 1:n_trees
2060
# pull the vertex node IDs
2061
current_line = split(file_lines[file_idx])
2062
element_node_ids[1] = parse(Int, current_line[2])
2063
element_node_ids[2] = parse(Int, current_line[3])
2064
element_node_ids[3] = parse(Int, current_line[4])
2065
element_node_ids[4] = parse(Int, current_line[5])
2066
element_node_ids[5] = parse(Int, current_line[6])
2067
element_node_ids[6] = parse(Int, current_line[7])
2068
element_node_ids[7] = parse(Int, current_line[8])
2069
element_node_ids[8] = parse(Int, current_line[9])
2070
2071
# Pull the (x, y, z) values of the eight vertices of the current tree out of the global vertices array
2072
for i in 1:8
2073
hex_vertices[:, i] .= vertices[:, element_node_ids[i]]
2074
end
2075
# Pull the information to check if boundary is curved/higher-order in order to read in additional data
2076
file_idx += 1
2077
current_line = split(file_lines[file_idx])
2078
curved_check[1] = parse(Int, current_line[2])
2079
curved_check[2] = parse(Int, current_line[3])
2080
curved_check[3] = parse(Int, current_line[4])
2081
curved_check[4] = parse(Int, current_line[5])
2082
curved_check[5] = parse(Int, current_line[6])
2083
curved_check[6] = parse(Int, current_line[7])
2084
if sum(curved_check) == 0
2085
# Create the node coordinates on this element
2086
calc_node_coordinates!(node_coordinates, tree, nodes, hex_vertices)
2087
else
2088
# Hexahedral element has at least one curved/higher-order side
2089
for face in 1:6
2090
if curved_check[face] == 0
2091
# Face is a flat plane.
2092
# Evaluate a bilinear interpolant between the four vertices
2093
# of the face at each of the nodes.
2094
get_vertices_for_bilinear_interpolant!(face_vertices, face,
2095
hex_vertices)
2096
for q in 1:nnodes, p in 1:nnodes
2097
@views bilinear_interpolation!(curve_values[:, p, q],
2098
face_vertices, nodes[p],
2099
nodes[q])
2100
end
2101
else # curved_check[face] == 1
2102
# Curved/higher-order face boundary information is supplied by
2103
# the mesh file. Just read it into a work array
2104
for q in 1:nnodes, p in 1:nnodes
2105
file_idx += 1
2106
current_line = split(file_lines[file_idx])
2107
curve_values[1, p, q] = parse(RealT, current_line[2])
2108
curve_values[2, p, q] = parse(RealT, current_line[3])
2109
curve_values[3, p, q] = parse(RealT, current_line[4])
2110
end
2111
end
2112
# Construct the curve interpolant for the current side
2113
face_curves[face] = CurvedFaceT(nodes, bary_weights, copy(curve_values))
2114
end
2115
# Create the node coordinates on this particular element
2116
calc_node_coordinates!(node_coordinates, tree, nodes, face_curves)
2117
end
2118
# Move file index to the next tree
2119
file_idx += 1
2120
end
2121
2122
return file_idx
2123
end
2124
2125
function face_curves_quadratic_3d!(face_curves, face_nodes, face,
2126
mesh_nodes, nodes, bary_weights,
2127
RealT, CurvedFaceT)
2128
node1_coords = mesh_nodes[face_nodes[1]]
2129
node2_coords = mesh_nodes[face_nodes[2]]
2130
node3_coords = mesh_nodes[face_nodes[3]]
2131
node4_coords = mesh_nodes[face_nodes[4]]
2132
node5_coords = mesh_nodes[face_nodes[5]]
2133
node6_coords = mesh_nodes[face_nodes[6]]
2134
node7_coords = mesh_nodes[face_nodes[7]]
2135
node8_coords = mesh_nodes[face_nodes[8]]
2136
node9_coords = mesh_nodes[face_nodes[9]]
2137
2138
# The nodes for an Abaqus element are labeled following a closed path
2139
# around the element:
2140
#
2141
# <----
2142
# 7 6 5
2143
# *-----*-----*
2144
# | |
2145
# | | | ^
2146
# | 8* 9* *4 |
2147
# v | | |
2148
# | |
2149
# *-----*-----*
2150
# 1 2 3
2151
# ---->
2152
# `curve_values`, however, requires to sort the nodes into a
2153
# valid coordinate system,
2154
#
2155
# *-----*-----*
2156
# | |
2157
# | |
2158
# * * *
2159
# η | |
2160
# ↑ | |
2161
# | *-----*-----*
2162
# |----> ξ
2163
# thus we need to flip the node order for the second xi and eta edges met.
2164
2165
nnodes = length(nodes)
2166
curve_values = Array{RealT}(undef, (3, nnodes, nnodes))
2167
2168
# Proceed along bottom edge
2169
curve_values[:, 1, 1] = node1_coords
2170
curve_values[:, 2, 1] = node2_coords
2171
curve_values[:, 3, 1] = node3_coords
2172
2173
# Proceed along middle line
2174
curve_values[:, 1, 2] = node8_coords
2175
curve_values[:, 2, 2] = node9_coords
2176
curve_values[:, 3, 2] = node4_coords
2177
2178
# Proceed along top edge
2179
curve_values[:, 1, 3] = node7_coords
2180
curve_values[:, 2, 3] = node6_coords
2181
curve_values[:, 3, 3] = node5_coords
2182
2183
# Construct the curve interpolant for the current side
2184
face_curves[face] = CurvedFaceT(nodes, bary_weights, copy(curve_values))
2185
2186
return nothing
2187
end
2188
2189
# Version for quadratic 3D elements, i.e., second-order hexes in
2190
# standard Abaqus .inp format as exported from e.g. Gmsh.
2191
function calc_tree_node_coordinates!(node_coordinates::AbstractArray{<:Any, 5},
2192
element_lines, nodes, vertices, RealT,
2193
linear_hexes, mesh_nodes)
2194
nnodes = length(nodes)
2195
2196
# Create a work set of Gamma curves to create the node coordinates
2197
CurvedFaceT = CurvedFace{RealT}
2198
face_curves = Array{CurvedFaceT}(undef, 6)
2199
2200
# Create other work arrays to perform the mesh construction
2201
hex_vertices = Array{RealT}(undef, (3, 8))
2202
face_nodes = Array{Int}(undef, 9)
2203
2204
# Create the barycentric weights used for the surface interpolations
2205
bary_weights_ = barycentric_weights(nodes)
2206
bary_weights = SVector{nnodes}(bary_weights_)
2207
2208
element_set_order = 1
2209
tree = 0
2210
for line_idx in 1:length(element_lines)
2211
line = element_lines[line_idx]
2212
2213
# Check if a new element type/element set is defined
2214
if startswith(line, "*ELEMENT")
2215
# Retrieve element type
2216
current_element_type = match(r"\*ELEMENT, type=([^,]+)", line).captures[1]
2217
2218
# Check if these are linear elements
2219
if occursin(linear_hexes, current_element_type)
2220
element_set_order = 1
2221
else
2222
element_set_order = 2
2223
end
2224
else # Element data
2225
tree += 1
2226
2227
# Pull the vertex node IDs
2228
line_split = split(line, r",\s+")
2229
# Store the indices of the element-defining nodes
2230
element_nodes = parse.(Int, line_split)[2:end] # First entry is the element number, can be ignored
2231
2232
if element_set_order == 1
2233
# Create the node coordinates on this particular element
2234
# Pull the (x,y,z) values of the four vertices of the current tree out of the global vertices array
2235
for i in 1:8
2236
hex_vertices[:, i] .= vertices[:, element_nodes[i]] # 3D => : = 1:3
2237
end
2238
calc_node_coordinates!(node_coordinates, tree, nodes, hex_vertices)
2239
else # element_set_order == 2
2240
# The second-order face has the following node distribution:
2241
# NW N NE
2242
# *-----*-----*
2243
# | |
2244
# | C |
2245
# W * * * E
2246
# | |
2247
# | |
2248
# *-----*-----*
2249
# SW S SE
2250
#
2251
# where SW gets placed at the origin of the local face coordinate system:
2252
#
2253
# η
2254
# ↑
2255
# |
2256
# |
2257
# |----> ξ
2258
2259
# We proceed now face by face, with face numbering as sketched out in
2260
# https://trixi-framework.github.io/TrixiDocumentation/stable/meshes/p4est_mesh/#HOHQMesh-Extended-Abaqus-format
2261
# and in `get_vertices_for_bilinear_interpolant!`
2262
#
2263
# The node selection follows from the sketch for the 27-node element presented in the doc
2264
# http://130.149.89.49:2080/v2016/books/usb/default.htm?startat=pt06ch28s01ael03.html
2265
2266
face = 1 # -y, "Front"
2267
face_nodes[1] = element_nodes[1] # "SW" node
2268
face_nodes[2] = element_nodes[9] # "S" node
2269
face_nodes[3] = element_nodes[2] # "SE" node
2270
face_nodes[4] = element_nodes[18] # "E" node
2271
face_nodes[5] = element_nodes[6] # "NE" node
2272
face_nodes[6] = element_nodes[13] # "N" node
2273
face_nodes[7] = element_nodes[5] # "NW" node
2274
face_nodes[8] = element_nodes[17] # "W" node
2275
face_nodes[9] = element_nodes[24] # "C" node
2276
2277
face_curves_quadratic_3d!(face_curves, face_nodes, face,
2278
mesh_nodes, nodes, bary_weights,
2279
RealT, CurvedFaceT)
2280
2281
face = 2 # +y, "Back"
2282
face_nodes[1] = element_nodes[4] # "SW" node
2283
face_nodes[2] = element_nodes[11] # "S" node
2284
face_nodes[3] = element_nodes[3] # "SE" node
2285
face_nodes[4] = element_nodes[19] # "E" node
2286
face_nodes[5] = element_nodes[7] # "NE" node
2287
face_nodes[6] = element_nodes[15] # "N" node
2288
face_nodes[7] = element_nodes[8] # "NW" node
2289
face_nodes[8] = element_nodes[20] # "W" node
2290
face_nodes[9] = element_nodes[26] # "C" node
2291
2292
face_curves_quadratic_3d!(face_curves, face_nodes, face,
2293
mesh_nodes, nodes, bary_weights,
2294
RealT, CurvedFaceT)
2295
2296
face = 3 # -z, "Bottom"
2297
face_nodes[1] = element_nodes[1] # "SW" node
2298
face_nodes[2] = element_nodes[9] # "S" node
2299
face_nodes[3] = element_nodes[2] # "SE" node
2300
face_nodes[4] = element_nodes[10] # "E" node
2301
face_nodes[5] = element_nodes[3] # "NE" node
2302
face_nodes[6] = element_nodes[11] # "N" node
2303
face_nodes[7] = element_nodes[4] # "NW" node
2304
face_nodes[8] = element_nodes[12] # "W" node
2305
face_nodes[9] = element_nodes[22] # "C" node
2306
2307
face_curves_quadratic_3d!(face_curves, face_nodes, face,
2308
mesh_nodes, nodes, bary_weights,
2309
RealT, CurvedFaceT)
2310
2311
# Care: Note the different node ordering compared to `get_vertices_for_bilinear_interpolant!`
2312
# Here, we obtain the node indices in the usual order, the "direction" is defined by the
2313
# corner nodes 2 and 3, which are the same as face 4.
2314
face = 4 # +x, "Right"
2315
face_nodes[1] = element_nodes[2] # "SW" node
2316
face_nodes[2] = element_nodes[10] # "S" node
2317
face_nodes[3] = element_nodes[3] # "SE" node
2318
face_nodes[4] = element_nodes[19] # "E" node
2319
face_nodes[5] = element_nodes[7] # "NE" node
2320
face_nodes[6] = element_nodes[14] # "N" node
2321
face_nodes[7] = element_nodes[6] # "NW" node
2322
face_nodes[8] = element_nodes[18] # "W" node
2323
face_nodes[9] = element_nodes[25] # "C" node
2324
2325
face_curves_quadratic_3d!(face_curves, face_nodes, face,
2326
mesh_nodes, nodes, bary_weights,
2327
RealT, CurvedFaceT)
2328
2329
face = 5 # +z, "Top"
2330
face_nodes[1] = element_nodes[5] # "SW" node
2331
face_nodes[2] = element_nodes[13] # "S" node
2332
face_nodes[3] = element_nodes[6] # "SE" node
2333
face_nodes[4] = element_nodes[14] # "E" node
2334
face_nodes[5] = element_nodes[7] # "NE" node
2335
face_nodes[6] = element_nodes[15] # "N" node
2336
face_nodes[7] = element_nodes[8] # "NW" node
2337
face_nodes[8] = element_nodes[16] # "W" node
2338
face_nodes[9] = element_nodes[23] # "C" node
2339
2340
face_curves_quadratic_3d!(face_curves, face_nodes, face,
2341
mesh_nodes, nodes, bary_weights,
2342
RealT, CurvedFaceT)
2343
2344
face = 6 # -x, "Left"
2345
face_nodes[1] = element_nodes[1] # "SW" node
2346
face_nodes[2] = element_nodes[12] # "S" node
2347
face_nodes[3] = element_nodes[4] # "SE" node
2348
face_nodes[4] = element_nodes[20] # "E" node
2349
face_nodes[5] = element_nodes[8] # "NE" node
2350
face_nodes[6] = element_nodes[16] # "N" node
2351
face_nodes[7] = element_nodes[5] # "NW" node
2352
face_nodes[8] = element_nodes[17] # "W" node
2353
face_nodes[9] = element_nodes[27] # "C" node
2354
2355
face_curves_quadratic_3d!(face_curves, face_nodes, face,
2356
mesh_nodes, nodes, bary_weights,
2357
RealT, CurvedFaceT)
2358
2359
# Note: `element_nodes[21]` remains unused, since it is not part of any face
2360
# (sits in the center of the element)
2361
2362
# Create the node coordinates on this particular element
2363
calc_node_coordinates!(node_coordinates, tree, nodes, face_curves)
2364
end
2365
end
2366
end
2367
2368
return nothing
2369
end
2370
2371
# Given the eight `hex_vertices` for a hexahedral element extract
2372
# the four `face_vertices` for a particular `face_index`.
2373
function get_vertices_for_bilinear_interpolant!(face_vertices, face_index, hex_vertices)
2374
if face_index == 1
2375
@views face_vertices[:, 1] .= hex_vertices[:, 1]
2376
@views face_vertices[:, 2] .= hex_vertices[:, 2]
2377
@views face_vertices[:, 3] .= hex_vertices[:, 6]
2378
@views face_vertices[:, 4] .= hex_vertices[:, 5]
2379
elseif face_index == 2
2380
@views face_vertices[:, 1] .= hex_vertices[:, 4]
2381
@views face_vertices[:, 2] .= hex_vertices[:, 3]
2382
@views face_vertices[:, 3] .= hex_vertices[:, 7]
2383
@views face_vertices[:, 4] .= hex_vertices[:, 8]
2384
elseif face_index == 3
2385
@views face_vertices[:, 1] .= hex_vertices[:, 1]
2386
@views face_vertices[:, 2] .= hex_vertices[:, 2]
2387
@views face_vertices[:, 3] .= hex_vertices[:, 3]
2388
@views face_vertices[:, 4] .= hex_vertices[:, 4]
2389
elseif face_index == 4
2390
@views face_vertices[:, 1] .= hex_vertices[:, 2]
2391
@views face_vertices[:, 2] .= hex_vertices[:, 3]
2392
@views face_vertices[:, 3] .= hex_vertices[:, 6]
2393
@views face_vertices[:, 4] .= hex_vertices[:, 7]
2394
elseif face_index == 5
2395
@views face_vertices[:, 1] .= hex_vertices[:, 5]
2396
@views face_vertices[:, 2] .= hex_vertices[:, 6]
2397
@views face_vertices[:, 3] .= hex_vertices[:, 7]
2398
@views face_vertices[:, 4] .= hex_vertices[:, 8]
2399
else # face_index == 6
2400
@views face_vertices[:, 1] .= hex_vertices[:, 1]
2401
@views face_vertices[:, 2] .= hex_vertices[:, 4]
2402
@views face_vertices[:, 3] .= hex_vertices[:, 8]
2403
@views face_vertices[:, 4] .= hex_vertices[:, 5]
2404
end
2405
2406
return nothing
2407
end
2408
2409
# Evaluate a bilinear interpolant at a point (u,v) given the four vertices where the face is right-handed
2410
# 4 3
2411
# o----------------o
2412
# | |
2413
# | |
2414
# | |
2415
# | |
2416
# | |
2417
# | |
2418
# o----------------o
2419
# 1 2
2420
# and return the 3D coordinate point (x, y, z)
2421
function bilinear_interpolation!(coordinate, face_vertices, u, v)
2422
for j in 1:3
2423
coordinate[j] = 0.25f0 * (face_vertices[j, 1] * (1 - u) * (1 - v)
2424
+ face_vertices[j, 2] * (1 + u) * (1 - v)
2425
+ face_vertices[j, 3] * (1 + u) * (1 + v)
2426
+ face_vertices[j, 4] * (1 - u) * (1 + v))
2427
end
2428
end
2429
2430
function get_global_first_element_ids(mesh::P4estMesh)
2431
return unsafe_wrap(Array, mesh.p4est.global_first_quadrant, mpi_nranks() + 1)
2432
end
2433
2434
function balance!(mesh::P4estMesh{2}, init_fn = C_NULL)
2435
p4est_balance(mesh.p4est, P4EST_CONNECT_FACE, init_fn)
2436
# Due to a bug in `p4est`, the forest needs to be rebalanced twice sometimes
2437
# See https://github.com/cburstedde/p4est/issues/112
2438
return p4est_balance(mesh.p4est, P4EST_CONNECT_FACE, init_fn)
2439
end
2440
2441
function balance!(mesh::P4estMesh{3}, init_fn = C_NULL)
2442
return p8est_balance(mesh.p4est, P8EST_CONNECT_FACE, init_fn)
2443
end
2444
2445
function partition!(mesh::P4estMesh{2}; weight_fn = C_NULL)
2446
return p4est_partition(mesh.p4est, Int(mesh.p4est_partition_allow_for_coarsening),
2447
weight_fn)
2448
end
2449
2450
function partition!(mesh::P4estMesh{3}; weight_fn = C_NULL)
2451
return p8est_partition(mesh.p4est, Int(mesh.p4est_partition_allow_for_coarsening),
2452
weight_fn)
2453
end
2454
2455
function update_ghost_layer!(mesh::P4estMesh)
2456
ghost_destroy_p4est(mesh.ghost)
2457
return mesh.ghost = PointerWrapper(ghost_new_p4est(mesh.p4est))
2458
end
2459
2460
function init_fn(p4est, which_tree, quadrant)
2461
# Unpack quadrant's user data ([global quad ID, controller_value])
2462
# Use `unsafe_load` here since `quadrant.p.user_data isa Ptr{Ptr{Nothing}}`
2463
# and we only need the first (only!) entry
2464
pw = PointerWrapper(Int, unsafe_load(quadrant.p.user_data))
2465
2466
# Initialize quad ID as -1 and controller_value as 0 (don't refine or coarsen)
2467
pw[1] = -1
2468
pw[2] = 0
2469
return nothing
2470
end
2471
2472
# 2D
2473
function cfunction(::typeof(init_fn), ::Val{2})
2474
@cfunction(init_fn, Cvoid,
2475
(Ptr{p4est_t}, Ptr{p4est_topidx_t}, Ptr{p4est_quadrant_t}))
2476
end
2477
# 3D
2478
function cfunction(::typeof(init_fn), ::Val{3})
2479
@cfunction(init_fn, Cvoid,
2480
(Ptr{p8est_t}, Ptr{p4est_topidx_t}, Ptr{p8est_quadrant_t}))
2481
end
2482
2483
function refine_fn(p4est, which_tree, quadrant)
2484
# Controller value has been copied to the quadrant's user data storage before.
2485
# Unpack quadrant's user data ([global quad ID, controller_value]).
2486
# Use `unsafe_load` here since `quadrant.p.user_data isa Ptr{Ptr{Nothing}}`
2487
# and we only need the first (only!) entry
2488
pw = PointerWrapper(Int, unsafe_load(quadrant.p.user_data))
2489
controller_value = pw[2]
2490
2491
if controller_value > 0
2492
# return true (refine)
2493
return Cint(1)
2494
else
2495
# return false (don't refine)
2496
return Cint(0)
2497
end
2498
end
2499
2500
# 2D
2501
function cfunction(::typeof(refine_fn), ::Val{2})
2502
@cfunction(refine_fn, Cint,
2503
(Ptr{p4est_t}, Ptr{p4est_topidx_t}, Ptr{p4est_quadrant_t}))
2504
end
2505
# 3D
2506
function cfunction(::typeof(refine_fn), ::Val{3})
2507
@cfunction(refine_fn, Cint,
2508
(Ptr{p8est_t}, Ptr{p4est_topidx_t}, Ptr{p8est_quadrant_t}))
2509
end
2510
2511
# Refine marked cells and rebalance forest.
2512
# Return a list of all cells that have been refined during refinement or rebalancing.
2513
function refine!(mesh::P4estMesh)
2514
# Copy original element IDs to quad user data storage
2515
original_n_cells = ncells(mesh)
2516
save_original_ids(mesh)
2517
2518
init_fn_c = cfunction(init_fn, Val(ndims(mesh)))
2519
refine_fn_c = cfunction(refine_fn, Val(ndims(mesh)))
2520
2521
# Refine marked cells
2522
@trixi_timeit timer() "refine" refine_p4est!(mesh.p4est, false, refine_fn_c,
2523
init_fn_c)
2524
2525
@trixi_timeit timer() "rebalance" balance!(mesh, init_fn_c)
2526
2527
return collect_changed_cells(mesh, original_n_cells)
2528
end
2529
2530
function coarsen_fn(p4est, which_tree, quadrants_ptr)
2531
quadrants = unsafe_wrap_quadrants(quadrants_ptr, p4est)
2532
2533
# Controller value has been copied to the quadrant's user data storage before.
2534
# Load controller value from quadrant's user data ([global quad ID, controller_value]).
2535
# Use `unsafe_load` here since `quadrant.p.user_data isa Ptr{Ptr{Nothing}}`
2536
# and we only need the first (only!) entry
2537
controller_value(i) = PointerWrapper(Int, unsafe_load(quadrants[i].p.user_data))[2]
2538
2539
# `p4est` calls this function for each 2^ndims quads that could be coarsened to a single one.
2540
# Only coarsen if all these 2^ndims quads have been marked for coarsening.
2541
if all(i -> controller_value(i) < 0, eachindex(quadrants))
2542
# return true (coarsen)
2543
return Cint(1)
2544
else
2545
# return false (don't coarsen)
2546
return Cint(0)
2547
end
2548
end
2549
2550
# 2D
2551
function unsafe_wrap_quadrants(quadrants_ptr, ::Ptr{p4est_t})
2552
return unsafe_wrap(Array, quadrants_ptr, 4)
2553
end
2554
# 3D
2555
function unsafe_wrap_quadrants(quadrants_ptr, ::Ptr{p8est_t})
2556
return unsafe_wrap(Array, quadrants_ptr, 8)
2557
end
2558
2559
# 2D
2560
function cfunction(::typeof(coarsen_fn), ::Val{2})
2561
@cfunction(coarsen_fn, Cint,
2562
(Ptr{p4est_t}, Ptr{p4est_topidx_t}, Ptr{Ptr{p4est_quadrant_t}}))
2563
end
2564
# 3D
2565
function cfunction(::typeof(coarsen_fn), ::Val{3})
2566
@cfunction(coarsen_fn, Cint,
2567
(Ptr{p8est_t}, Ptr{p4est_topidx_t}, Ptr{Ptr{p8est_quadrant_t}}))
2568
end
2569
2570
# Coarsen marked cells if the forest will stay balanced.
2571
# Return a list of all cells that have been coarsened.
2572
function coarsen!(mesh::P4estMesh)
2573
# Copy original element IDs to quad user data storage
2574
original_n_cells = ncells(mesh)
2575
save_original_ids(mesh)
2576
2577
# Coarsen marked cells
2578
coarsen_fn_c = cfunction(coarsen_fn, Val(ndims(mesh)))
2579
init_fn_c = cfunction(init_fn, Val(ndims(mesh)))
2580
2581
@trixi_timeit timer() "coarsen!" coarsen_p4est!(mesh.p4est, false, coarsen_fn_c,
2582
init_fn_c)
2583
2584
# IDs of newly created cells (one-based)
2585
new_cells = collect_new_cells(mesh)
2586
# Old IDs of cells that have been coarsened (one-based)
2587
coarsened_cells_vec = collect_changed_cells(mesh, original_n_cells)
2588
# 2^ndims changed cells should have been coarsened to one new cell.
2589
# This matrix will store the IDs of all cells that have been coarsened to cell new_cells[i]
2590
# in the i-th column.
2591
coarsened_cells = reshape(coarsened_cells_vec, 2^ndims(mesh), length(new_cells))
2592
2593
# Save new original IDs to find out what changed after balancing
2594
intermediate_n_cells = ncells(mesh)
2595
save_original_ids(mesh)
2596
2597
@trixi_timeit timer() "rebalance" balance!(mesh, init_fn_c)
2598
2599
refined_cells = collect_changed_cells(mesh, intermediate_n_cells)
2600
2601
# Some cells may have been coarsened even though they unbalanced the forest.
2602
# These cells have now been refined again by p4est_balance.
2603
# refined_cells contains the intermediate IDs (ID of coarse cell
2604
# between coarsening and balancing) of these cells.
2605
# Find original ID of each cell that has been coarsened and then refined again.
2606
for refined_cell in refined_cells
2607
# i-th cell of the ones that have been created by coarsening has been refined again
2608
i = findfirst(==(refined_cell), new_cells)
2609
2610
# Remove IDs of the 2^ndims cells that have been coarsened to this cell
2611
coarsened_cells[:, i] .= -1
2612
end
2613
2614
# Return all IDs of cells that have been coarsened but not refined again by balancing
2615
return coarsened_cells_vec[coarsened_cells_vec .>= 0]
2616
end
2617
2618
# Copy global quad ID to quad's user data storage, will be called below
2619
function save_original_id_iter_volume(info, user_data)
2620
info_pw = PointerWrapper(info)
2621
2622
# Load tree from global trees array, one-based indexing
2623
tree_pw = load_pointerwrapper_tree(info_pw.p4est, info_pw.treeid[] + 1)
2624
# Quadrant numbering offset of this quadrant
2625
offset = tree_pw.quadrants_offset[]
2626
# Global quad ID
2627
quad_id = offset + info_pw.quadid[]
2628
2629
# Unpack quadrant's user data ([global quad ID, controller_value])
2630
pw = PointerWrapper(Int, info_pw.quad.p.user_data[])
2631
# Save global quad ID
2632
pw[1] = quad_id
2633
return nothing
2634
end
2635
2636
# 2D
2637
function cfunction(::typeof(save_original_id_iter_volume), ::Val{2})
2638
@cfunction(save_original_id_iter_volume, Cvoid,
2639
(Ptr{p4est_iter_volume_info_t}, Ptr{Cvoid}))
2640
end
2641
# 3D
2642
function cfunction(::typeof(save_original_id_iter_volume), ::Val{3})
2643
@cfunction(save_original_id_iter_volume, Cvoid,
2644
(Ptr{p8est_iter_volume_info_t}, Ptr{Cvoid}))
2645
end
2646
2647
# Copy old element IDs to each quad's user data storage
2648
function save_original_ids(mesh::P4estMesh)
2649
iter_volume_c = cfunction(save_original_id_iter_volume, Val(ndims(mesh)))
2650
2651
return iterate_p4est(mesh.p4est, C_NULL; iter_volume_c = iter_volume_c)
2652
end
2653
2654
# Extract information about which cells have been changed
2655
function collect_changed_iter_volume(info, user_data)
2656
info_pw = PointerWrapper(info)
2657
2658
# The original element ID has been saved to user_data before.
2659
# Load original quad ID from quad's user data ([global quad ID, controller_value]).
2660
quad_data_pw = PointerWrapper(Int, info_pw.quad.p.user_data[])
2661
original_id = quad_data_pw[1]
2662
2663
# original_id of cells that have been newly created is -1
2664
if original_id >= 0
2665
# Unpack user_data = original_cells
2666
user_data_pw = PointerWrapper(Int, user_data)
2667
2668
# If quad has an original_id, it existed before refinement/coarsening,
2669
# and therefore wasn't changed.
2670
# Mark original_id as "not changed during refinement/coarsening" in original_cells
2671
user_data_pw[original_id + 1] = 0
2672
end
2673
return nothing
2674
end
2675
2676
# 2D
2677
function cfunction(::typeof(collect_changed_iter_volume), ::Val{2})
2678
@cfunction(collect_changed_iter_volume, Cvoid,
2679
(Ptr{p4est_iter_volume_info_t}, Ptr{Cvoid}))
2680
end
2681
# 3D
2682
function cfunction(::typeof(collect_changed_iter_volume), ::Val{3})
2683
@cfunction(collect_changed_iter_volume, Cvoid,
2684
(Ptr{p8est_iter_volume_info_t}, Ptr{Cvoid}))
2685
end
2686
2687
function collect_changed_cells(mesh::P4estMesh, original_n_cells)
2688
original_cells = collect(1:original_n_cells)
2689
2690
# Iterate over all quads and set original cells that haven't been changed to zero
2691
iter_volume_c = cfunction(collect_changed_iter_volume, Val(ndims(mesh)))
2692
2693
iterate_p4est(mesh.p4est, original_cells; iter_volume_c = iter_volume_c)
2694
2695
# Changed cells are all that haven't been set to zero above
2696
changed_original_cells = original_cells[original_cells .> 0]
2697
2698
return changed_original_cells
2699
end
2700
2701
# Extract newly created cells
2702
function collect_new_iter_volume(info, user_data)
2703
info_pw = PointerWrapper(info)
2704
2705
# The original element ID has been saved to user_data before.
2706
# Unpack quadrant's user data ([global quad ID, controller_value]).
2707
original_id = PointerWrapper(Int, info_pw.quad.p.user_data[])[1]
2708
2709
# original_id of cells that have been newly created is -1
2710
if original_id < 0
2711
# Load tree from global trees array, one-based indexing
2712
tree_pw = load_pointerwrapper_tree(info_pw.p4est, info_pw.treeid[] + 1)
2713
# Quadrant numbering offset of this quadrant
2714
offset = tree_pw.quadrants_offset[]
2715
# Global quad ID
2716
quad_id = offset + info_pw.quadid[]
2717
2718
# Unpack user_data = original_cells
2719
user_data_pw = PointerWrapper(Int, user_data)
2720
2721
# Mark cell as "newly created during refinement/coarsening/balancing"
2722
user_data_pw[quad_id + 1] = 1
2723
end
2724
return nothing
2725
end
2726
2727
# 2D
2728
function cfunction(::typeof(collect_new_iter_volume), ::Val{2})
2729
@cfunction(collect_new_iter_volume, Cvoid,
2730
(Ptr{p4est_iter_volume_info_t}, Ptr{Cvoid}))
2731
end
2732
# 3D
2733
function cfunction(::typeof(collect_new_iter_volume), ::Val{3})
2734
@cfunction(collect_new_iter_volume, Cvoid,
2735
(Ptr{p8est_iter_volume_info_t}, Ptr{Cvoid}))
2736
end
2737
2738
function collect_new_cells(mesh::P4estMesh)
2739
cell_is_new = zeros(Int, ncells(mesh))
2740
2741
# Iterate over all quads and set original cells that have been changed to one
2742
iter_volume_c = cfunction(collect_new_iter_volume, Val(ndims(mesh)))
2743
2744
iterate_p4est(mesh.p4est, cell_is_new; iter_volume_c = iter_volume_c)
2745
2746
# Changed cells are all that haven't been set to zero above
2747
new_cells = findall(==(1), cell_is_new)
2748
2749
return new_cells
2750
end
2751
end # @muladd
2752
2753