Skip to content
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
9 changes: 6 additions & 3 deletions src/auxiliary/precompile.jl
Original file line number Diff line number Diff line change
Expand Up @@ -407,7 +407,8 @@ function _precompile_manual_()
# 2D, serial
@assert Base.precompile(Tuple{typeof(Trixi.init_boundaries), Array{Int, 1},
TreeMesh{2, Trixi.SerialTree{2}, RealT},
Trixi.TreeElementContainer2D{RealT, uEltype}})
Trixi.TreeElementContainer2D{RealT, uEltype},
basis_type_dgsem(RealT, nnodes_)})
@assert Base.precompile(Tuple{typeof(Trixi.init_interfaces), Array{Int, 1},
TreeMesh{2, Trixi.SerialTree{2}, RealT},
Trixi.TreeElementContainer2D{RealT, uEltype}})
Expand All @@ -421,7 +422,8 @@ function _precompile_manual_()
# 2D, parallel
@assert Base.precompile(Tuple{typeof(Trixi.init_boundaries), Array{Int, 1},
TreeMesh{2, Trixi.ParallelTree{2}, RealT},
Trixi.TreeElementContainer2D{RealT, uEltype}})
Trixi.TreeElementContainer2D{RealT, uEltype},
basis_type_dgsem(RealT, nnodes_)})
@assert Base.precompile(Tuple{typeof(Trixi.init_interfaces), Array{Int, 1},
TreeMesh{2, Trixi.ParallelTree{2}, RealT},
Trixi.TreeElementContainer2D{RealT, uEltype}})
Expand All @@ -438,7 +440,8 @@ function _precompile_manual_()
# 3D, serial
@assert Base.precompile(Tuple{typeof(Trixi.init_boundaries), Array{Int, 1},
TreeMesh{3, Trixi.SerialTree{3}, RealT},
Trixi.TreeElementContainer3D{RealT, uEltype}})
Trixi.TreeElementContainer3D{RealT, uEltype},
basis_type_dgsem(RealT, nnodes_)})
@assert Base.precompile(Tuple{typeof(Trixi.init_interfaces), Array{Int, 1},
TreeMesh{3, Trixi.SerialTree{3}, RealT},
Trixi.TreeElementContainer3D{RealT, uEltype}})
Expand Down
33 changes: 19 additions & 14 deletions src/solvers/dgsem_tree/containers_1d.jl
Original file line number Diff line number Diff line change
Expand Up @@ -377,15 +377,16 @@ end
function set_boundary_node_coordinates!(boundaries, element, count, direction,
elements, mesh::TreeMesh1D,
basis::LobattoLegendreBasis)
el_node_coords = elements.node_coordinates
bnd_node_coords = boundaries.node_coordinates

orientation = 1 # always 1 in 1D
if direction == 1
@views boundaries.node_coordinates[:, count] .= elements.node_coordinates[orientation,
1,
element]
bnd_node_coords[orientation, count] = el_node_coords[orientation, 1,
element]
elseif direction == 2
@views boundaries.node_coordinates[:, count] .= elements.node_coordinates[orientation,
end,
element]
bnd_node_coords[orientation, count] = el_node_coords[orientation, end,
element]
else
error("should not happen")
end
Expand All @@ -397,17 +398,21 @@ end
function set_boundary_node_coordinates!(boundaries, element, count, direction,
elements, mesh::TreeMesh1D,
basis::GaussLegendreBasis)
boundary_matrix = basis.boundary_interpolation
el_node_coords = elements.node_coordinates
bnd_node_coords = boundaries.node_coordinates

orientation = 1 # always 1 in 1D
if direction == 1
@views x_interpolated_left = dot(basis.boundary_interpolation[:, 1],
elements.node_coordinates[orientation, :,
element])
boundaries.node_coordinates[:, count] .= x_interpolated_left
@views x_interpolated_left = dot(boundary_matrix[:, 1],
el_node_coords[orientation, :,
element])
bnd_node_coords[orientation, count] = x_interpolated_left
elseif direction == 2
@views x_interpolated_right = dot(basis.boundary_interpolation[:, 2],
elements.node_coordinates[orientation, :,
element])
boundaries.node_coordinates[:, count] .= x_interpolated_right
@views x_interpolated_right = dot(boundary_matrix[:, 2],
el_node_coords[orientation, :,
element])
bnd_node_coords[orientation, count] = x_interpolated_right
else
error("should not happen")
end
Expand Down
107 changes: 89 additions & 18 deletions src/solvers/dgsem_tree/containers_2d.jl
Original file line number Diff line number Diff line change
Expand Up @@ -353,15 +353,15 @@ end

# Create boundaries container and initialize boundary data in `elements`.
function init_boundaries(cell_ids, mesh::TreeMesh2D,
elements::TreeElementContainer2D)
elements::TreeElementContainer2D, basis)
# Initialize container
n_boundaries = count_required_boundaries(mesh, cell_ids)
boundaries = TreeBoundaryContainer2D{real(elements), eltype(elements)}(n_boundaries,
nvariables(elements),
nnodes(elements))

# Connect elements with boundaries
init_boundaries!(boundaries, elements, mesh)
init_boundaries!(boundaries, elements, mesh, basis)
return boundaries
end

Expand Down Expand Up @@ -390,8 +390,89 @@ function count_required_boundaries(mesh::TreeMesh2D, cell_ids)
return count
end

# For Lobatto points, we can simply use the outer nodes of the elements as boundary nodes.
function set_boundary_node_coordinates!(boundaries, element, count, direction,
elements, mesh::TreeMesh2D,
basis::LobattoLegendreBasis)
el_node_coords = elements.node_coordinates
bnd_node_coords = boundaries.node_coordinates

if direction == 1 # -x direction
@views bnd_node_coords[:, :, count] .= el_node_coords[:, 1, :, element]
elseif direction == 2 # +x direction
@views bnd_node_coords[:, :, count] .= el_node_coords[:, end, :, element]
elseif direction == 3 # -y direction
@views bnd_node_coords[:, :, count] .= el_node_coords[:, :, 1, element]
elseif direction == 4 # +y direction
@views bnd_node_coords[:, :, count] .= el_node_coords[:, :, end, element]
else
error("should not happen")
end

return nothing
end

# For Gauss points, we need to interpolate the boundary node coordinates.
function set_boundary_node_coordinates!(boundaries, element, count, direction,
elements, mesh::TreeMesh2D,
basis::GaussLegendreBasis)
boundary_matrix = basis.boundary_interpolation
el_node_coords = elements.node_coordinates
bnd_node_coords = boundaries.node_coordinates

if direction == 1 # -x direction: interpolate in x for each y node j
for j in eachnode(basis)
for orientation in 1:2 # Need to set both x and y coordinate of boundary node
@views bnd_node_coords[orientation, j, count] = dot(boundary_matrix[:,
1],
el_node_coords[orientation,
:,
j,
element])
end
end
elseif direction == 2 # +x direction: interpolate in x for each y node j
for j in eachnode(basis)
for orientation in 1:2 # Need to set both x and y coordinate of boundary node
@views bnd_node_coords[orientation, j, count] = dot(boundary_matrix[:,
2],
el_node_coords[orientation,
:,
j,
element])
end
end
elseif direction == 3 # -y direction: interpolate in y for each x node i
for i in eachnode(basis)
for orientation in 1:2 # Need to set both x and y coordinate of boundary node
@views bnd_node_coords[orientation, i, count] = dot(boundary_matrix[:,
1],
el_node_coords[orientation,
i,
:,
element])
end
end
elseif direction == 4 # +y direction: interpolate in y for each x node i
for i in eachnode(basis)
for orientation in 1:2 # Need to set both x and y coordinate of boundary node
@views bnd_node_coords[orientation, i, count] = dot(boundary_matrix[:,
2],
el_node_coords[orientation,
i,
:,
element])
end
end
else
error("should not happen")
end

return nothing
end

# Initialize connectivity between elements and boundaries
function init_boundaries!(boundaries, elements, mesh::TreeMesh2D)
function init_boundaries!(boundaries, elements, mesh::TreeMesh2D, basis)
# Exit early if there are no boundaries to initialize
if nboundaries(boundaries) == 0
# In this case n_boundaries_per_direction still needs to be reset!
Expand Down Expand Up @@ -441,24 +522,14 @@ function init_boundaries!(boundaries, elements, mesh::TreeMesh2D)

# Set orientation (x -> 1, y -> 2)
if direction in (1, 2)
boundaries.orientations[count] = 1
boundaries.orientations[count] = 1 # x direction
else
boundaries.orientations[count] = 2
boundaries.orientations[count] = 2 # y direction
end

# Store node coordinates
enc = elements.node_coordinates
if direction == 1 # -x direction
boundaries.node_coordinates[:, :, count] .= enc[:, 1, :, element]
elseif direction == 2 # +x direction
boundaries.node_coordinates[:, :, count] .= enc[:, end, :, element]
elseif direction == 3 # -y direction
boundaries.node_coordinates[:, :, count] .= enc[:, :, 1, element]
elseif direction == 4 # +y direction
boundaries.node_coordinates[:, :, count] .= enc[:, :, end, element]
else
error("should not happen")
end
set_boundary_node_coordinates!(boundaries, element, count, direction,
elements, mesh, basis)
end
end

Expand Down Expand Up @@ -1349,7 +1420,7 @@ function reinitialize_containers!(mesh::Union{TreeMesh{2}, TreeMesh{3}}, equatio
# re-initialize boundaries container
@unpack boundaries = cache
resize!(boundaries, count_required_boundaries(mesh, leaf_cell_ids))
init_boundaries!(boundaries, elements, mesh)
init_boundaries!(boundaries, elements, mesh, dg.basis)

# re-initialize mortars container
@unpack mortars = cache
Expand Down
6 changes: 3 additions & 3 deletions src/solvers/dgsem_tree/containers_3d.jl
Original file line number Diff line number Diff line change
Expand Up @@ -350,15 +350,15 @@ end

# Create boundaries container and initialize boundary data in `elements`.
function init_boundaries(cell_ids, mesh::TreeMesh3D,
elements::TreeElementContainer3D)
elements::TreeElementContainer3D, basis)
# Initialize container
n_boundaries = count_required_boundaries(mesh, cell_ids)
boundaries = TreeBoundaryContainer3D{real(elements), eltype(elements)}(n_boundaries,
nvariables(elements),
nnodes(elements))

# Connect elements with boundaries
init_boundaries!(boundaries, elements, mesh)
init_boundaries!(boundaries, elements, mesh, basis)
return boundaries
end

Expand Down Expand Up @@ -388,7 +388,7 @@ function count_required_boundaries(mesh::TreeMesh3D, cell_ids)
end

# Initialize connectivity between elements and boundaries
function init_boundaries!(boundaries, elements, mesh::TreeMesh3D)
function init_boundaries!(boundaries, elements, mesh::TreeMesh3D, basis)
# Reset boundaries count
count = 0

Expand Down
4 changes: 0 additions & 4 deletions src/solvers/dgsem_tree/dg_1d.jl
Original file line number Diff line number Diff line change
Expand Up @@ -561,12 +561,8 @@ function prolong2boundaries!(cache, u_or_flux_viscous,
end
else # Element in +x direction of boundary => need to evaluate at left boundary node (-1)
for v in eachvariable(equations)
# Interpolate to the boundaries using a local variable for
# the accumulation of values (to reduce global memory operations).
boundary_u_2 = zero(eltype(boundaries.u))
for ii in eachnode(dg)
# Not += to allow `@muladd` to turn these into FMAs
# (see comment at the top of the file)
boundary_u_2 = (boundary_u_2 +
u_or_flux_viscous[v, ii, element] *
boundary_interpolation[ii, 1])
Expand Down
70 changes: 69 additions & 1 deletion src/solvers/dgsem_tree/dg_2d.jl
Original file line number Diff line number Diff line change
Expand Up @@ -20,7 +20,7 @@ function create_cache(mesh::Union{TreeMesh{2}, TreeMesh{3}}, equations,

interfaces = init_interfaces(leaf_cell_ids, mesh, elements)

boundaries = init_boundaries(leaf_cell_ids, mesh, elements)
boundaries = init_boundaries(leaf_cell_ids, mesh, elements, dg.basis)

mortars = init_mortars(leaf_cell_ids, mesh, elements, dg.mortar)

Expand Down Expand Up @@ -699,6 +699,74 @@ function prolong2boundaries!(cache, u,
return nothing
end

function prolong2boundaries!(cache, u,
mesh::TreeMesh{2}, equations,
dg::DG{<:GaussLegendreBasis})
@unpack boundaries = cache
@unpack orientations, neighbor_sides = boundaries
@unpack boundary_interpolation = dg.basis

@threaded for boundary in eachboundary(dg, cache)
element = boundaries.neighbor_ids[boundary]

if orientations[boundary] == 1
# boundary in x-direction
if neighbor_sides[boundary] == 1
# element in -x direction of boundary => interpolate to right boundary node (+1)
for l in eachnode(dg), v in eachvariable(equations)
# Interpolate to the boundaries using a local variable for
# the accumulation of values (to reduce global memory operations).
boundary_u = zero(eltype(boundaries.u))
for ii in eachnode(dg)
# Not += to allow `@muladd` to turn these into FMAs
# (see comment at the top of the file)
boundary_u = (boundary_u +
u[v, ii, l, element] *
boundary_interpolation[ii, 2])
end
boundaries.u[1, v, l, boundary] = boundary_u
end
else # element in +x direction of boundary => interpolate to left boundary node (-1)
for l in eachnode(dg), v in eachvariable(equations)
boundary_u = zero(eltype(boundaries.u))
for ii in eachnode(dg)
boundary_u = (boundary_u +
u[v, ii, l, element] *
boundary_interpolation[ii, 1])
end
boundaries.u[2, v, l, boundary] = boundary_u
end
end
else # if orientations[boundary] == 2
# boundary in y-direction
if neighbor_sides[boundary] == 1
# element in -y direction of boundary => interpolate to right boundary node (+1)
for l in eachnode(dg), v in eachvariable(equations)
boundary_u = zero(eltype(boundaries.u))
for ii in eachnode(dg)
boundary_u = (boundary_u +
u[v, l, ii, element] *
boundary_interpolation[ii, 2])
end
boundaries.u[1, v, l, boundary] = boundary_u
end
else # element in +y direction of boundary => interpolate to left boundary node (-1)
for l in eachnode(dg), v in eachvariable(equations)
boundary_u = zero(eltype(boundaries.u))
for ii in eachnode(dg)
boundary_u = (boundary_u +
u[v, l, ii, element] *
boundary_interpolation[ii, 1])
end
boundaries.u[2, v, l, boundary] = boundary_u
end
end
end
end

return nothing
end

function calc_boundary_flux!(cache, t, boundary_conditions::NamedTuple,
mesh::TreeMesh{2}, equations, surface_integral, dg::DG)
@unpack surface_flux_values = cache.elements
Expand Down
2 changes: 1 addition & 1 deletion src/solvers/dgsem_tree/dg_2d_parallel.jl
Original file line number Diff line number Diff line change
Expand Up @@ -253,7 +253,7 @@ function create_cache(mesh::TreeMeshParallel{2}, equations,

mpi_interfaces = init_mpi_interfaces(leaf_cell_ids, mesh, elements)

boundaries = init_boundaries(leaf_cell_ids, mesh, elements)
boundaries = init_boundaries(leaf_cell_ids, mesh, elements, dg.basis)

mortars = init_mortars(leaf_cell_ids, mesh, elements, dg.mortar)

Expand Down
Loading
Loading