Skip to content

Commit 42de0a4

Browse files
authored
Merge pull request #64 from WIAS-PDELib/chmerdon/curl4hdiv
Curl2D and Curl3D evaluations for Hdiv
2 parents 8b9bde0 + b94b896 commit 42de0a4

4 files changed

Lines changed: 152 additions & 3 deletions

File tree

CHANGELOG.md

Lines changed: 4 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -1,5 +1,9 @@
11
# CHANGES
22

3+
## v1.6.0
4+
- added Curl2D and Curl3D evaluation for Hdiv
5+
6+
37
## v1.5.0 November 21, 2025
48
- added type parameter to `PointEvaluator` so `lazy_interpolate!` is now accessible for automatic differentiation
59
- added new `compute_lazy_interpolation_jacobian` function for computing Jacobians of interpolations between finite element spaces

Project.toml

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -1,6 +1,6 @@
11
name = "ExtendableFEMBase"
22
uuid = "12fb9182-3d4c-4424-8fd1-727a0899810c"
3-
version = "1.5.1"
3+
version = "1.6.0"
44
authors = ["Christian Merdon <merdon@wias-berlin.de>", "Patrick Jaap <patrick.jaap@wias-berlin.de>"]
55

66
[deps]

src/feevaluator_hdiv.jl

Lines changed: 46 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -106,3 +106,49 @@ function update_basis!(FEBE::SingleFEEvaluator{<:Real, <:Real, <:Integer, <:Grad
106106
end
107107
return nothing
108108
end
109+
110+
111+
# CURL2D HDIV
112+
function update_basis!(FEBE::SingleFEEvaluator{<:Real, <:Real, <:Integer, <:Curl2D, <:AbstractHdivFiniteElement})
113+
L2GAinv = _update_trafo!(FEBE)
114+
L2GM = _update_piola!(FEBE)
115+
subset = _update_subset!(FEBE)
116+
coefficients = _update_coefficients!(FEBE)
117+
cvals = FEBE.cvals
118+
offsets2 = FEBE.offsets2
119+
refbasisderivvals = FEBE.refbasisderivvals
120+
fill!(cvals, 0)
121+
det = FEBE.L2G.det # 1 alloc
122+
for i in 1:size(cvals, 3), dof_i in 1:size(cvals, 2)
123+
for j in 1:size(L2GM, 2), m in 1:size(L2GAinv, 2)
124+
cvals[1, dof_i, i] -= L2GAinv[2, m] * L2GM[1, j] * refbasisderivvals[subset[dof_i] + offsets2[j], m, i] * coefficients[1, dof_i] / det
125+
cvals[1, dof_i, i] += L2GAinv[1, m] * L2GM[2, j] * refbasisderivvals[subset[dof_i] + offsets2[j], m, i] * coefficients[2, dof_i] / det
126+
end
127+
end
128+
return nothing
129+
end
130+
131+
132+
# CURL3D HDIV
133+
function update_basis!(FEBE::SingleFEEvaluator{<:Real, <:Real, <:Integer, <:Curl3D, <:AbstractHdivFiniteElement})
134+
L2GAinv = _update_trafo!(FEBE)
135+
L2GM = _update_piola!(FEBE)
136+
subset = _update_subset!(FEBE)
137+
coefficients = _update_coefficients!(FEBE)
138+
cvals = FEBE.cvals
139+
offsets2 = FEBE.offsets2
140+
refbasisderivvals = FEBE.refbasisderivvals
141+
fill!(cvals, 0)
142+
det = FEBE.L2G.det # 1 alloc
143+
for i in 1:size(cvals, 3), dof_i in 1:size(cvals, 2)
144+
for j in 1:size(L2GM, 2), m in 1:size(L2GAinv, 2)
145+
cvals[1, dof_i, i] -= L2GAinv[3, m] * L2GM[2, j] * refbasisderivvals[subset[dof_i] + offsets2[j], m, i] * coefficients[2, dof_i] / det
146+
cvals[1, dof_i, i] += L2GAinv[2, m] * L2GM[3, j] * refbasisderivvals[subset[dof_i] + offsets2[j], m, i] * coefficients[3, dof_i] / det
147+
cvals[2, dof_i, i] -= L2GAinv[1, m] * L2GM[3, j] * refbasisderivvals[subset[dof_i] + offsets2[j], m, i] * coefficients[3, dof_i] / det
148+
cvals[2, dof_i, i] += L2GAinv[3, m] * L2GM[1, j] * refbasisderivvals[subset[dof_i] + offsets2[j], m, i] * coefficients[1, dof_i] / det
149+
cvals[3, dof_i, i] -= L2GAinv[2, m] * L2GM[1, j] * refbasisderivvals[subset[dof_i] + offsets2[j], m, i] * coefficients[1, dof_i] / det
150+
cvals[3, dof_i, i] += L2GAinv[1, m] * L2GM[2, j] * refbasisderivvals[subset[dof_i] + offsets2[j], m, i] * coefficients[2, dof_i] / det
151+
end
152+
end
153+
return nothing
154+
end

test/test_operators.jl

Lines changed: 101 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -6,8 +6,12 @@ function run_operator_tests()
66
println("============================")
77
error = test_derivatives2D()
88
@test error < 1.0e-14
9+
error = test_derivatives2D_hdiv()
10+
@test error < 1.0e-14
911
error = test_derivatives3D()
1012
@test error < 1.0e-14
13+
error = test_derivatives3D_hdiv()
14+
@test error < 1.0e-14
1115
test_reconstructions()
1216
end
1317
end
@@ -79,8 +83,7 @@ function test_derivatives2D()
7983
xgrid = grid_triangle([-1.0 0.0; 1.0 0.0; 0.0 1.0]')
8084

8185
## define P2-Courant finite element space
82-
FEType = H1P2{2, 2}
83-
FES = FESpace{FEType}(xgrid)
86+
FES = FESpace{H1P2{2, 2}}(xgrid)
8487
show(devnull, FES)
8588

8689
## get midpoint quadrature rule for constants
@@ -140,6 +143,52 @@ function test_derivatives2D()
140143
end
141144

142145

146+
function test_derivatives2D_hdiv()
147+
## define test function and expected operator evals
148+
function testf(result, qpinfo)
149+
x = qpinfo.x
150+
result[1] = x[1]^2
151+
return result[2] = 3 * x[2]^2 + x[1] * x[2]
152+
end
153+
154+
## expected values of operators in cell midpoint
155+
expected_curl2 = [1 / 3]
156+
157+
## define grid = a single non-refenrece triangle
158+
xgrid = grid_triangle([-1.0 0.0; 1.0 0.0; 0.0 1.0]')
159+
160+
## define P2-Courant finite element space
161+
FES = FESpace{HDIVBDM2{2}}(xgrid)
162+
show(devnull, FES)
163+
164+
## get midpoint quadrature rule for constants
165+
qf = QuadratureRule{Float64, Triangle2D}(0)
166+
167+
## define FE basis Evaluator for Hessian
168+
FEBE_curl2 = FEEvaluator(FES, Curl2D, qf)
169+
170+
## update on cell 1
171+
update_basis!(FEBE_curl2, 1)
172+
173+
## interpolate quadratic testfunction
174+
Iu = FEVector(FES)
175+
interpolate!(Iu[1], testf)
176+
177+
## check if operator evals have the correct length
178+
@assert size(FEBE_curl2.cvals, 1) == length(expected_curl2)
179+
180+
# compute derivatives
181+
curl2 = zeros(Float64, 1)
182+
eval_febe!(curl2, FEBE_curl2, Iu.entries[FES[CellDofs][:, 1]], 1)
183+
184+
## compute errors to expected values
185+
error_curl2 = sqrt(sum((curl2 - expected_curl2) .^ 2))
186+
println("EG = Triangle2D | operator = Curl2 | error = $error_curl2")
187+
188+
return maximum([error_curl2])
189+
end
190+
191+
143192
function test_derivatives3D()
144193
## define test function and expected operator evals
145194
function testf(result, qpinfo)
@@ -220,3 +269,53 @@ function test_derivatives3D()
220269

221270
return maximum([error_curl3, error_L, error_H, error_symH, error_symH2])
222271
end
272+
273+
274+
function test_derivatives3D_hdiv()
275+
## define test function and expected operator evals
276+
function testf(result, qpinfo)
277+
x = qpinfo.x
278+
result[1] = x[2]
279+
result[2] = 3 * x[3]
280+
return result[3] = x[2]
281+
end
282+
283+
## expected values of operators in cell midpoint
284+
expected_curl3 = [1 - 3, 0, -1]
285+
286+
## define grid = a single non-refenrece triangle
287+
xgrid = reference_domain(Tetrahedron3D)
288+
xgrid[Coordinates][:, 2] = [2, 0, 0]
289+
290+
## define P2-Courant finite element space
291+
FES = FESpace{HDIVRT1{3}}(xgrid)
292+
show(devnull, FES)
293+
294+
## get midpoint quadrature rule for constants
295+
qf = QuadratureRule{Float64, Tetrahedron3D}(0)
296+
297+
## define FE basis Evaluator for Hessian
298+
FEBE_curl3 = FEEvaluator(FES, Curl3D, qf)
299+
300+
## update on cell 1
301+
update_basis!(FEBE_curl3, 1)
302+
303+
## interpolate quadratic testfunction
304+
Iu = FEVector(FES)
305+
interpolate!(Iu[1], testf)
306+
307+
## check if operator evals have the correct length
308+
@assert size(FEBE_curl3.cvals, 1) == length(expected_curl3)
309+
310+
## eval 2nd order derivatives at only quadrature point 1
311+
## since function is quadratic this should be constant
312+
curl3 = zeros(Float64, 3)
313+
eval_febe!(curl3, FEBE_curl3, Iu.entries[FES[CellDofs][:, 1]], 1)
314+
@info curl3
315+
316+
## compute errors to expected values
317+
error_curl3 = sqrt(sum((curl3 - expected_curl3) .^ 2))
318+
println("EG = Tetrahedron3D | operator = Curl3 | error = $error_curl3")
319+
320+
return maximum([error_curl3])
321+
end

0 commit comments

Comments
 (0)