-
Notifications
You must be signed in to change notification settings - Fork 2
/
Copy pathchebpts.m
230 lines (215 loc) · 7.41 KB
/
chebpts.m
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
function [x w v] = chebpts(n,d,kind)
%CHEBPTS Chebyshev points in [-1,1].
% CHEBPTS(N) returns N Chebyshev points of the 2nd-kind in [-1,1].
%
% CHEBPTS(N,D), where D is vector of length 2 and N is a scalar integer,
% scales the nodes and weights for the interval [D(1) D(2)]. If the
% interval is infinite, the map is chosen to be the default 'unbounded
% map' with mappref('parinf') = [1 0] and mappref('adaptinf') = 0. If
% length(D) > 2 and N a vector of length(D)-1, then CHEBPTS returns a
% column vector of the stacked N(k) Chebyshev points on the subintervals
% D(k:k+1). If length(N) is 1, then D is treated as [D(1) D(end)].
%
% [X W] = CHEBPTS(N,D) returns also a row vector of the (scaled) weights
% for Clenshaw-Curtis quadrature (computed using [1]). (For nodes and
% weights of Gauss-Chebyshev quadrature, use [X W] = JACPTS(N,-.5,-.5,D))
%
% [X W V] = CHEBPTS(N,D) returns, in addition to X and W, the barycentric
% weights V corresponding to the Chebyshev points X.
%
% [X W V] = CHEBPTS(F) returns the Chebyshev nodes and weights
% corresponding to the domain and length of the chebfun F.
%
% [X W V] = CHEBPTS(N,KIND) or CHEBPTS(N,D,KIND) returns Chebyshev points
% and weights of the 1st-kind if KIND = 1 and 2nd-kind if KIND = 2
% (default). (Note that if KIND is not supplied, chebpts will always
% return 2nd-kind points, regardless of the value of 'chebkind' in
% chebfunpref.).
%
% See also legpts, jacpts, lagpts, and hermpts.
% Copyright 2011 by The University of Oxford and The Chebfun Developers.
% See http://www.maths.ox.ac.uk/chebfun/ for Chebfun information.
% [1] Jörg Waldvogel, "Fast construction of the Fejér and Clenshaw-Curtis
% quadrature rules", BIT Numerical Mathematics 46 (2006), pp 195-202.
% Intialise
x = []; w = []; v = []; scale = false;
% Shortcut for simple input (2nd-kind points on [-1 1]).
if nargin == 1 && isnumeric(n) && numel(n) == 1 && n > 0
if n == 1, x = 0; w = 2; v = 1; return, end % Special case
m = n-1;
x = sin(pi*(-m:2:m)/(2*m)).'; % Chebyshev points
if nargout > 1 % Quadrature weights
w = weights2(n);
end
if nargout > 2 % Barycentric weights
v = [.5 ; ones(n-1,1)]; v(2:2:end) = -1; v(end) = .5*v(end);
end
return
end
% Parse the inputs
if isa(n,'chebfun')
if numel(n) > 1
error('CHEBFUN:chebpts:quasi',...
'chebpts does not support quasi-matrices.');
end
nn = zeros(n.nfuns,1);
for k = 1:n.nfuns
nn(k) = n.funs(k).n;
end
if nargin == 1, d = 2; end % Default to 2nd-kind
if nargout == 1
x = chebpts(nn,n.ends,d);
elseif nargout == 2
[x w] = chebpts(nn,n.ends,d);
else
[x w v] = chebpts(nn,n.ends,d);
end
return
elseif nargin == 1
d = [-1 1];
kind = 2;
elseif nargin == 2
if isa(d,'domain')
scale = true;
kind = 2;
elseif length(d) == 1
kind = d;
d = [-1 1];
else
scale = true;
kind = 2;
end
elseif nargin == 3
scale = true;
end
if isa(d,'domain')
d = d.endsandbreaks;
end
if isempty(d) || ~any(n)
return % Return empty vector if domain is empty or n == 0
end
if numel(n) == 1
d = d([1 end]);
end
% Deal with the piecewise case (where d has breakpoints and n is a vector).
numints = numel(d)-1;
if numints > 1
if numel(n) ~= numints
error('CHEBFUN:chebpts:numints', ...
'Vector N does not match domain D.');
end
csn = cumsum([0 ; n(:)]);
x = zeros(csn(end),1);
if nargout == 1
for k = 1:numints
idxk = csn(k)+1:csn(k+1);
x(idxk) = chebpts(n(k),d(k:k+1),kind);
end
elseif nargout == 2
w = zeros(1,csn(end));
for k = 1:numints
idxk = csn(k)+1:csn(k+1);
[x(idxk) w(idxk)] = chebpts(n(k),d(k:k+1),kind);
end
else
w = zeros(1,csn(end)); v = zeros(csn(end),1);
for k = 1:numints
idxk = csn(k)+1:csn(k+1);
[x(idxk) w(idxk) v(idxk)] = chebpts(n(k),d(k:k+1),kind);
end
end
return
end
if numel(n) > 1,
error('CHEBFUN:chebpts:vecn','Vector N does not match domain D.');
end
% Avoid unnecessary scaling
if (d(1)==-1 && d(2)==1), scale = false; end
% Allow strings to determine which kind of points
if ischar(kind)
if strcmpi(kind,'1st'), kind = 1;
elseif strcmpi(kind,'2nd'), kind = 2; end
end
if n < 0,
error('CHEBFUN:chebpts:posinpt',...
'Input should be a nonnegative number');
elseif n == 1,
x = 0; w = 2; v = 1;
else
m = n-1;
if kind == 1
x = sin(pi*(-m:2:m)/(2*m+2)).'; % 1st-kind Chebyshev points
if nargout > 1 % Quadrature weights
w = weights1(n);
end
if nargout > 2 % Barycentric weights
v = sin((2*(0:n-1)+1)*pi/(2*n)).'; v(2:2:end) = -v(2:2:end);
if ~mod(n,2), v = v./max(abs(v)); end
end
else
x = sin(pi*(-m:2:m)/(2*m)).'; % 2nd-kind Chebyshev points
if nargout > 1 % Quadrature weights
w = weights2(n);
end
if nargout > 2 % Barycentric weights
v = [.5 ; ones(n-1,1)]; v(2:2:end) = -1; v(end) = .5*v(end);
end
end
end
% Rescale if d is provided:
if scale
if ~any(isinf(d)) % Finite interval
dab05 = .5*diff(d);
x = x*dab05 + (d(1) + dab05);
if ( kind == 2 )
x([1,end]) = d([1,end]);
end
w = dab05*w;
else % Infinite interval
m = maps(fun,{'unbounded'},d); % Use default map
if nargout > 1 % Quadrature weights
w = w.*m.der(x.');
if isinf(d(1)), w(1) = 0; end
if isinf(d(end)), w(end) = 0; end
end
x = m.for(x);
if kind == 2 % Force endpoints for 2nd-kind points
x([1 end]) = d([1 end]);
end
end
end
function w = weights1(n) % 1st-kind Chebyshev weights
% Jörg Waldvogel, "Fast construction of the Fejér and Clenshaw-Curtis
% quadrature rules", BIT Numerical Mathematics 43 (1), p. 001-018 (2004).
% http://www2.maths.ox.ac.uk/chebfun/and_beyond/programme/slides/wald.pdf
if n == 1
w = 2;
else
% new
L = 0:n-1; r = 2./(1-4*min(L,n-L).^2); s1 = sign(n/2-L); % Aux vecs
w = real(ifft(s1.*r.*exp(1i*pi/n*L))); % Fejer weights
% % old
% l = floor(n/2)+1;
% K = 0:n-l;
% v = [2*exp(1i*pi*K/n)./(1-4*K.^2) zeros(1,l)];
% w = real(ifft(v(1:n) + conj(v(n+1:-1:2))));
end
function w = weights2(n) % 2nd-kind Chebyshev wieghts
% Jörg Waldvogel, "Fast construction of the Fejér and Clenshaw-Curtis
% quadrature rules", BIT Numerical Mathematics 43 (1), p. 001-018 (2004).
% http://www2.maths.ox.ac.uk/chebfun/and_beyond/programme/slides/wald.pdf
if n == 1
w = 2;
else
% new
n = n-1;
u0 = 1/(n^2-1+mod(n,2)); % Boundary weights
L = 0:n-1; r = 2./(1-4*min(L,n-L).^2); % Auxiliary vectors
w = [ifft(r-u0) u0]; % C-C weights
% % old
% m = n-1;
% c = zeros(1,n);
% c(1:2:n) = 2./[1 1-(2:2:m).^2 ];
% f = real(ifft([c(1:n) c(m:-1:2)]));
% w = [f(1) 2*f(2:m) f(n)];
end