Skip to content

Commit 16fe1a5

Browse files
authored
add docstrings for Jacobi and Legendre (#215)
1 parent 9d02cab commit 16fe1a5

File tree

4 files changed

+216
-26
lines changed

4 files changed

+216
-26
lines changed

docs/make.jl

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

3+
DocMeta.setdocmeta!(ClassicalOrthogonalPolynomials,
4+
:DocTestSetup,
5+
:(using ClassicalOrthogonalPolynomials))
6+
37
makedocs(
48
modules = [ClassicalOrthogonalPolynomials],
59
sitename="ClassicalOrthogonalPolynomials.jl",

docs/src/index.md

Lines changed: 16 additions & 24 deletions
Original file line numberDiff line numberDiff line change
@@ -126,17 +126,15 @@ U\T
126126

127127
```@docs
128128
ClassicalOrthogonalPolynomials.Chebyshev
129-
```
130-
```@docs
131129
ClassicalOrthogonalPolynomials.chebyshevt
132-
```
133-
```@docs
134130
ClassicalOrthogonalPolynomials.chebyshevu
135131
```
136132
```@docs
133+
ClassicalOrthogonalPolynomials.Legendre
137134
ClassicalOrthogonalPolynomials.legendrep
138135
```
139136
```@docs
137+
ClassicalOrthogonalPolynomials.Jacobi
140138
ClassicalOrthogonalPolynomials.jacobip
141139
```
142140
```@docs
@@ -147,8 +145,6 @@ ClassicalOrthogonalPolynomials.hermiteh
147145
```
148146

149147

150-
151-
152148
### Weights
153149

154150
```@docs
@@ -161,7 +157,9 @@ ClassicalOrthogonalPolynomials.HermiteWeight
161157
ClassicalOrthogonalPolynomials.Weighted
162158
```
163159
```@docs
160+
ClassicalOrthogonalPolynomials.LegendreWeight
164161
ClassicalOrthogonalPolynomials.ChebyshevWeight
162+
ClassicalOrthogonalPolynomials.JacobiWeight
165163
```
166164
```@docs
167165
ClassicalOrthogonalPolynomials.LaguerreWeight
@@ -170,6 +168,13 @@ ClassicalOrthogonalPolynomials.LaguerreWeight
170168
ClassicalOrthogonalPolynomials.HalfWeighted
171169
```
172170

171+
### Affine-mapped
172+
```@docs
173+
ClassicalOrthogonalPolynomials.legendre
174+
ClassicalOrthogonalPolynomials.jacobi
175+
ClassicalOrthogonalPolynomials.legendreweight
176+
ClassicalOrthogonalPolynomials.jacobiweight
177+
```
173178

174179
### Recurrences
175180

@@ -190,32 +195,22 @@ ClassicalOrthogonalPolynomials.recurrencecoefficients
190195
```
191196

192197

193-
194198
### Internal
195199

196200
```@docs
197-
ClassicalOrthogonalPolynomials.ShuffledIR2HC
198-
```
199-
```@docs
200-
ClassicalOrthogonalPolynomials.ShuffledR2HC
201-
```
202-
```@docs
201+
ClassicalOrthogonalPolynomials.ShuffledFFT
203202
ClassicalOrthogonalPolynomials.ShuffledIFFT
203+
ClassicalOrthogonalPolynomials.ShuffledR2HC
204+
ClassicalOrthogonalPolynomials.ShuffledIR2HC
204205
```
205206
```@docs
206207
ClassicalOrthogonalPolynomials.qr_jacobimatrix
207-
```
208-
```@docs
209-
ClassicalOrthogonalPolynomials.MappedOPLayout
210-
```
211-
```@docs
212208
ClassicalOrthogonalPolynomials.cholesky_jacobimatrix
213209
```
214210
```@docs
215211
ClassicalOrthogonalPolynomials.AbstractNormalizedOPLayout
216-
```
217-
```@docs
218-
ClassicalOrthogonalPolynomials.ShuffledFFT
212+
ClassicalOrthogonalPolynomials.MappedOPLayout
213+
ClassicalOrthogonalPolynomials.WeightedOPLayout
219214
```
220215
```@docs
221216
ClassicalOrthogonalPolynomials.legendre_grammatrix
@@ -224,9 +219,6 @@ ClassicalOrthogonalPolynomials.legendre_grammatrix
224219
ClassicalOrthogonalPolynomials.weightedgrammatrix
225220
```
226221
```@docs
227-
ClassicalOrthogonalPolynomials.WeightedOPLayout
228-
```
229-
```@docs
230222
ClassicalOrthogonalPolynomials.interlace!
231223
```
232224
```@docs

src/classical/jacobi.jl

Lines changed: 98 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -21,13 +21,48 @@ broadcasted(::LazyQuasiArrayStyle{1}, ::typeof(sqrt), w::AbstractJacobiWeight) =
2121
broadcasted(::LazyQuasiArrayStyle{1}, ::typeof(Base.literal_pow), ::Base.RefValue{typeof(^)}, w::AbstractJacobiWeight, ::Base.RefValue{Val{k}}) where k =
2222
JacobiWeight(k * w.a, k * w.b)
2323

24+
"""
25+
JacobiWeight{T}(a,b)
26+
JacobiWeight(a,b)
27+
28+
The quasi-vector representing the Jacobi weight function ``(1-x)^a (1+x)^b`` on ``[-1,1]``. See also [`jacobiweight`](@ref) and [`Jacobi`](@ref).
29+
# Examples
30+
```jldoctest
31+
julia> J=JacobiWeight(1.0,1.0)
32+
(1-x)^1.0 * (1+x)^1.0 on -1..1
33+
34+
julia> J[0.5]
35+
0.75
36+
37+
julia> axes(J)
38+
(Inclusion(-1.0 .. 1.0 (Chebyshev)),)
39+
```
40+
"""
2441
struct JacobiWeight{T} <: AbstractJacobiWeight{T}
2542
a::T
2643
b::T
2744
JacobiWeight{T}(a, b) where T = new{T}(convert(T,a), convert(T,b))
2845
end
2946

3047
JacobiWeight(a::V, b::T) where {T,V} = JacobiWeight{promote_type(T,V)}(a,b)
48+
49+
"""
50+
jacobiweight(a,b, d::AbstractInterval)
51+
52+
The [`JacobiWeight`](@ref) affine-mapped to interval `d`.
53+
54+
# Examples
55+
```jldoctest
56+
julia> J = jacobiweight(1, 1, 0..1)
57+
(1-x)^1 * (1+x)^1 on -1..1 affine mapped to 0 .. 1
58+
59+
julia> axes(J)
60+
(Inclusion(0 .. 1),)
61+
62+
julia> J[0.5]
63+
1.0
64+
```
65+
"""
3166
jacobiweight(a,b, d::AbstractInterval{T}) where T = JacobiWeight(a,b)[affine(d,ChebyshevInterval{T}())]
3267

3368
AbstractQuasiArray{T}(w::JacobiWeight) where T = JacobiWeight{T}(w.a, w.b)
@@ -92,7 +127,41 @@ include("legendre.jl")
92127
singularitiesbroadcast(::typeof(*), ::LegendreWeight, b::AbstractJacobiWeight) = b
93128
singularitiesbroadcast(::typeof(*), a::AbstractJacobiWeight, ::LegendreWeight) = a
94129

95-
130+
"""
131+
Jacobi{T}(a,b)
132+
Jacobi(a,b)
133+
134+
The quasi-matrix representing Jacobi polynomials, where the first axes represents the interval and the second axes represents the polynomial index (starting from 1). See also [`jacobi`](@ref), [`jacobip`](@ref) and [`JacobiWeight`](@ref).
135+
136+
The eltype, when not specified, will be converted to a floating point data type.
137+
# Examples
138+
```jldoctest
139+
julia> J=Jacobi(0,0) # The eltype will be converted to float
140+
Jacobi(0.0, 0.0)
141+
142+
julia> axes(J)
143+
(Inclusion(-1.0 .. 1.0 (Chebyshev)), OneToInf())
144+
145+
julia> J[0,:] # Values of polynomials at x=0
146+
ℵ₀-element view(::Jacobi{Float64}, 0.0, :) with eltype Float64 with indices OneToInf():
147+
1.0
148+
0.0
149+
-0.5
150+
-0.0
151+
0.375
152+
0.0
153+
-0.3125
154+
-0.0
155+
0.2734375
156+
0.0
157+
158+
159+
julia> J0=J[:,1]; # J0 is the first Jacobi polynomial which is constant.
160+
161+
julia> J0[0],J0[0.5]
162+
(1.0, 1.0)
163+
```
164+
"""
96165
struct Jacobi{T} <: AbstractJacobi{T}
97166
a::T
98167
b::T
@@ -104,7 +173,34 @@ Jacobi(a::V, b::T) where {T,V} = Jacobi{float(promote_type(T,V))}(a, b)
104173
AbstractQuasiArray{T}(w::Jacobi) where T = Jacobi{T}(w.a, w.b)
105174
AbstractQuasiMatrix{T}(w::Jacobi) where T = Jacobi{T}(w.a, w.b)
106175

107-
176+
"""
177+
jacobi(a,b, d::AbstractInterval)
178+
179+
The [`Jacobi`](@ref) polynomials affine-mapped to interval `d`.
180+
181+
# Examples
182+
```jldoctest
183+
julia> J = jacobi(1, 1, 0..1)
184+
Jacobi(1.0, 1.0) affine mapped to 0 .. 1
185+
186+
julia> axes(J)
187+
(Inclusion(0 .. 1), OneToInf())
188+
189+
julia> J[0,:]
190+
ℵ₀-element view(::Jacobi{Float64}, -1.0, :) with eltype Float64 with indices OneToInf():
191+
1.0
192+
-2.0
193+
3.0
194+
-4.0
195+
5.0
196+
-6.0
197+
7.0
198+
-8.0
199+
9.0
200+
-10.0
201+
202+
```
203+
"""
108204
jacobi(a,b) = Jacobi(a,b)
109205
jacobi(a,b, d::AbstractInterval{T}) where T = Jacobi{float(promote_type(eltype(a),eltype(b),T))}(a,b)[affine(d,ChebyshevInterval{T}()), :]
110206

src/classical/legendre.jl

Lines changed: 98 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -1,5 +1,40 @@
1+
"""
2+
LegendreWeight{T}()
3+
LegendreWeight()
4+
5+
The quasi-vector representing the Legendre weight function (const 1) on [-1,1]. See also [`legendreweight`](@ref) and [`Legendre`](@ref).
6+
# Examples
7+
```jldoctest
8+
julia> w = LegendreWeight()
9+
LegendreWeight{Float64}()
10+
11+
julia> w[0.5]
12+
1.0
13+
14+
julia> axes(w)
15+
(Inclusion(-1.0 .. 1.0 (Chebyshev)),)
16+
```
17+
"""
118
struct LegendreWeight{T} <: AbstractJacobiWeight{T} end
219
LegendreWeight() = LegendreWeight{Float64}()
20+
21+
"""
22+
legendreweight(d::AbstractInterval)
23+
24+
The [`LegendreWeight`](@ref) affine-mapped to interval `d`.
25+
26+
# Examples
27+
```jldoctest
28+
julia> w = legendreweight(0..1)
29+
LegendreWeight{Float64}() affine mapped to 0 .. 1
30+
31+
julia> axes(w)
32+
(Inclusion(0 .. 1),)
33+
34+
julia> w[0]
35+
1.0
36+
```
37+
"""
338
legendreweight(d::AbstractInterval{T}) where T = LegendreWeight{float(T)}()[affine(d,ChebyshevInterval{T}())]
439

540
AbstractQuasiArray{T}(::LegendreWeight) where T = LegendreWeight{T}()
@@ -46,6 +81,41 @@ singularities(::AbstractFillLayout, P) = LegendreWeight{eltype(P)}()
4681
basis_singularities(::LegendreWeight{T}) where T = Legendre{T}()
4782
basis_singularities(v::SubQuasiArray) = view(basis_singularities(parent(v)), parentindices(v)[1], :)
4883

84+
"""
85+
Legendre{T=Float64}(a,b)
86+
87+
The quasi-matrix representing Legendre polynomials, where the first axes represents the interval and the second axes represents the polynomial index (starting from 1). See also [`legendre`](@ref), [`legendrep`](@ref), [`LegendreWeight`](@ref) and [`Jacobi`](@ref).
88+
# Examples
89+
```jldoctest
90+
julia> P = Legendre()
91+
Legendre()
92+
93+
julia> typeof(P) # default eltype
94+
Legendre{Float64}
95+
96+
julia> axes(P)
97+
(Inclusion(-1.0 .. 1.0 (Chebyshev)), OneToInf())
98+
99+
julia> P[0,:] # Values of polynomials at x=0
100+
ℵ₀-element view(::Legendre{Float64}, 0.0, :) with eltype Float64 with indices OneToInf():
101+
1.0
102+
0.0
103+
-0.5
104+
-0.0
105+
0.375
106+
0.0
107+
-0.3125
108+
-0.0
109+
0.2734375
110+
0.0
111+
112+
113+
julia> P₀=P[:,1]; # P₀ is the first Legendre polynomial which is constant.
114+
115+
julia> P₀[0],P₀[0.5]
116+
(1.0, 1.0)
117+
```
118+
"""
49119
struct Legendre{T} <: AbstractJacobi{T} end
50120
Legendre() = Legendre{Float64}()
51121

@@ -57,6 +127,34 @@ weighted(P::Normalized{<:Any,<:Legendre}) = P
57127
weighted(P::SubQuasiArray{<:Any,2,<:Legendre}) = P
58128
weighted(P::SubQuasiArray{<:Any,2,<:Normalized{<:Any,<:Legendre}}) = P
59129

130+
"""
131+
legendre(d::AbstractInterval)
132+
133+
The [`Legendre`](@ref) polynomials affine-mapped to interval `d`.
134+
135+
# Examples
136+
```jldoctest
137+
julia> P = legendre(0..1)
138+
Legendre() affine mapped to 0 .. 1
139+
140+
julia> axes(P)
141+
(Inclusion(0 .. 1), OneToInf())
142+
143+
julia> P[0.5,:]
144+
ℵ₀-element view(::Legendre{Float64}, 0.0, :) with eltype Float64 with indices OneToInf():
145+
1.0
146+
0.0
147+
-0.5
148+
-0.0
149+
0.375
150+
0.0
151+
-0.3125
152+
-0.0
153+
0.2734375
154+
0.0
155+
156+
```
157+
"""
60158
legendre() = Legendre()
61159
legendre(d::AbstractInterval{T}) where T = Legendre{float(T)}()[affine(d,ChebyshevInterval{T}()), :]
62160
legendre(d::ChebyshevInterval{T}) where T = Legendre{float(T)}()

0 commit comments

Comments
 (0)