Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Fix bug with L2projection of mixed grid #456

Merged
merged 1 commit into from
Jun 20, 2022
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
4 changes: 2 additions & 2 deletions src/L2_projection.jl
Original file line number Diff line number Diff line change
Expand Up @@ -239,11 +239,11 @@ function _project(vars, proj::L2Projector, fe_values::Values, M::Integer, ::Type
get_data(x::Number, i) = x

## Assemble contributions from each cell
for cellnum in proj.set
for (ic,cellnum) in enumerate(proj.set)
celldofs!(cell_dofs, proj.dh, cellnum)
fill!(fe, 0)
Xe = getcoordinates(proj.dh.grid, cellnum)
cell_vars = vars[cellnum]
cell_vars = vars[ic]
reinit!(fe_values, Xe)

for q_point = 1:nqp
Expand Down
35 changes: 35 additions & 0 deletions test/test_l2_projection.jl
Original file line number Diff line number Diff line change
Expand Up @@ -129,6 +129,8 @@ function test_projection_mixedgrid()
push!(cells, Triangle((2,3,6)))
push!(cells, Triangle((2,6,5)))

quadset = 1:1
triaset = 2:3
mesh = Grid(cells, nodes)

order = 2
Expand Down Expand Up @@ -175,6 +177,39 @@ function test_projection_mixedgrid()
@test isnan(point_vars_4[node][d1, d2])
end
end


#
#Do the same thing but for the triangle set
#
order = 2
ip = Lagrange{dim, RefTetrahedron, order}()
ip_geom = Lagrange{dim, RefTetrahedron, 1}()
qr = QuadratureRule{dim, RefTetrahedron}(4)
cv = CellScalarValues(qr, ip, ip_geom)
nqp = getnquadpoints(cv)

qp_values_tria = [zeros(SymmetricTensor{2,2}, nqp) for _ in triaset]
qp_values_matrix_tria = [zero(SymmetricTensor{2,2}) for _ in 1:nqp, _ in triaset]
for (ic, cellid) in enumerate(triaset)
xe = getcoordinates(mesh, cellid)
ae = compute_vertex_values(mesh, f)
# analytical values
qp_values = [f(spatial_coordinate(cv, qp, xe)) for qp in 1:getnquadpoints(cv)]
qp_values_tria[ic] = qp_values
qp_values_matrix_tria[:, ic] .= qp_values
end

#tria
proj = L2Projector(ip, mesh; geom_ip=ip_geom, set=triaset)
point_vars = project(proj, qp_values_tria, qr)
point_vars_2 = project(proj, qp_values_matrix_tria, qr)
for cellid in triaset
for node in mesh.cells[cellid].nodes
@test ae[node] ≈ point_vars[node] ≈ point_vars_2[node]
end
end

end

function test_node_reordering()
Expand Down