Skip to content

Commit

Permalink
fix P[Curl]Grad coverage and API for 1D + isHierarchical cov
Browse files Browse the repository at this point in the history
  • Loading branch information
Antoinemarteau committed Jan 6, 2025
1 parent 6ecc309 commit 9555fe9
Show file tree
Hide file tree
Showing 8 changed files with 44 additions and 7 deletions.
5 changes: 5 additions & 0 deletions src/Polynomials/MonomialBases.jl
Original file line number Diff line number Diff line change
Expand Up @@ -26,6 +26,11 @@ MonomialBasis(args...) = UniformPolyBasis(Monomial, args...)
function PGradBasis(::Type{Monomial},::Val{D},::Type{T},order::Int) where {D,T}
NedelecPolyBasisOnSimplex{D}(Monomial,T,order)
end
function PGradBasis(::Type{Monomial},::Val{1},::Type{T},order::Int) where T
@check T<:Real "T needs to be <:Real since represents the type of the components of the vector value"
V = VectorValue{1,T}
UniformPolyBasis(Monomial, Val(1), V, order+1)
end

# 1D evaluation implementation

Expand Down
4 changes: 2 additions & 2 deletions src/Polynomials/NedelecPolyBases.jl
Original file line number Diff line number Diff line change
Expand Up @@ -2,6 +2,7 @@ struct NedelecPolyBasisOnSimplex{D,V,K,PT} <: PolynomialBasis{D,V,K,PT}
order::Int
function NedelecPolyBasisOnSimplex{D}(::Type{PT},::Type{T},order::Integer) where {D,PT,T}
@check T<:Real "T needs to be <:Real since represents the type of the components of the vector value"
@notimplementedif !(D in (2,3))
K = Int(order)+1
V = VectorValue{D,T}
new{D,V,K,PT}(Int(order))
Expand Down Expand Up @@ -242,7 +243,7 @@ Return a basis of
ℕ𝔻ᴰₙ(△) = (ℙₙ)ᴰ ⊕ x × (ℙₙ \\ ℙₙ₋₁)ᴰ
with n=`order`, the polynomial space for Nedelec elements on `D`-dimensional
simplices with scalar type `T`.
simplices with scalar type `T`. `D` must be 1, 2 or 3.
The `order`=n argument has the following meaning: the curl of the functions in
this basis is in ℙₙ.
Expand All @@ -262,4 +263,3 @@ function PGradBasis(::Type{PT},::Val{D},::Type{T},order::Int) where {PT,D,T}
@notimplemented "only implemented for monomials"
end


6 changes: 6 additions & 0 deletions test/PolynomialsTests/BernsteinBasesTests.jl
Original file line number Diff line number Diff line change
Expand Up @@ -5,6 +5,8 @@ using Gridap.TensorValues
using Gridap.Fields
using Gridap.Polynomials

@test isHierarchical(Bernstein) == false

np = 3
x = [Point(0.),Point(1.),Point(.4)]
xi = x[1]
Expand All @@ -15,6 +17,7 @@ V = Float64
G = gradient_type(V,xi)
H = gradient_type(G,xi)


# order 0 degenerated case

order = 0
Expand All @@ -32,6 +35,7 @@ Hbx = repeat(permutedims(h),np)
test_field_array(b,x,bx,grad=∇bx,gradgrad=Hbx)
test_field_array(b,x[1],bx[1,:],grad=∇bx[1,:],gradgrad=Hbx[1,:])


# Order 1

order = 1
Expand All @@ -52,6 +56,7 @@ Hbx = H[ 0. 0.
test_field_array(b,x,bx,grad=∇bx,gradgrad=Hbx)
test_field_array(b,x[1],bx[1,:],grad=∇bx[1,:],gradgrad=Hbx[1,:])


# Order 2

order = 2
Expand All @@ -73,6 +78,7 @@ Hbx = H[ 2. -4. 2.
test_field_array(b,x,bx,, grad=∇bx,gradgrad=Hbx)
test_field_array(b,x[1],bx[1,:],grad=∇bx[1,:],gradgrad=Hbx[1,:])


# Order 3

order = 3
Expand Down
2 changes: 2 additions & 0 deletions test/PolynomialsTests/LegendreBasesTests.jl
Original file line number Diff line number Diff line change
Expand Up @@ -6,6 +6,8 @@ using Gridap.Fields
using Gridap.Fields: Broadcasting
using Gridap.Polynomials

@test isHierarchical(Legendre) == true

# Real-valued Q space with isotropic order

x1 = Point(0.0)
Expand Down
2 changes: 2 additions & 0 deletions test/PolynomialsTests/MonomialBasesTests.jl
Original file line number Diff line number Diff line change
Expand Up @@ -7,6 +7,8 @@ using Gridap.Polynomials

using Gridap.Polynomials: _q_filter, _qs_filter, _p_filter, _ps_filter

@test isHierarchical(Monomial) == true

xi = Point(2,3)
np = 5
x = fill(xi,np)
Expand Down
20 changes: 18 additions & 2 deletions test/PolynomialsTests/PCurlGradBasesTests.jl
Original file line number Diff line number Diff line change
Expand Up @@ -10,12 +10,14 @@ xi = Point(4,2)
np = 1
x = fill(xi,np)

PT = Monomial

order = 2
D = 2
T = Float64
V = VectorValue{D,T}
G = gradient_type(V,xi)
b = PCurlGradBasis(Monomial, Val(D),T,order)
b = PCurlGradBasis(PT, Val(D),T,order)

v = V[
(1.0, 0.0), (4.0, 0.0), (16.0, 0.0), (2.0, 0.0), (8.0, 0.0), (4.0, 0.0), # pterm ex
Expand Down Expand Up @@ -57,9 +59,23 @@ D = 3
T = Float64
V = VectorValue{D,T}
G = gradient_type(V,xi)
b = PCurlGradBasis(Monomial, Val(D),T,order)
b = PCurlGradBasis(PT, Val(D),T,order)

@test length(b) == 15
@test get_order(b) == 2


# 1D

order = 0
D = 1
T = Float64
V = VectorValue{D,T}
b = PCurlGradBasis(PT,Val(D),T,order)

@test b isa UniformPolyBasis{D,V,order+1,PT}

@test_throws AssertionError PCurlGradBasis(PT,Val(D),V,order)


end # module
11 changes: 9 additions & 2 deletions test/PolynomialsTests/PGradBasesTests.jl
Original file line number Diff line number Diff line change
Expand Up @@ -11,6 +11,11 @@ np = 3
x = fill(xi,np)
T = Float64

D = 2
for PT in (Legendre, Chebyshev, ModalC0, Bernstein)
@test_throws ErrorException PGradBasis(PT,Val(D),T,0)
end

PT = Monomial

order = 0
Expand Down Expand Up @@ -46,20 +51,22 @@ bx = repeat(permutedims(v),np)
test_field_array(b,x,bx,grad=∇bx)
test_field_array(b,x[1],bx[1,:],grad=∇bx[1,:])


# 1D

order = 0
D = 1
T = Float64
V = VectorValue{D,T}
b = QGradBasis(PT,Val(D),T,order)
b = PGradBasis(PT,Val(D),T,order)

@test b isa UniformPolyBasis{D,V,order+1,PT}

@test_throws AssertionError QGradBasis(PT,Val(D),V,order)
@test_throws AssertionError PGradBasis(PT,Val(D),V,order)

# D > 3 not implemented

@test_throws ErrorException NedelecPolyBasisOnSimplex{4}(PT, T, order)


end # module
1 change: 0 additions & 1 deletion test/PolynomialsTests/QCurlGradBasesTests.jl
Original file line number Diff line number Diff line change
Expand Up @@ -92,7 +92,6 @@ order = 0
D = 1
T = Float64
V = VectorValue{D,T}
G = gradient_type(V,xi)
b = QCurlGradBasis(PT,Val(D),T,order)

@test b isa UniformPolyBasis{D,V,order+1,PT}
Expand Down

0 comments on commit 9555fe9

Please sign in to comment.