0%

Chebyshev 谱方法及其算例

谱方法的正交多项式为 Chebyshev 多项式,为 Chebyshev 谱方法。这里介绍一维问题的 Chebyshev 谱方法,然后 MATLAB 编程算出数值解以及误差。

两点边值问题

考虑两点边值问题:

Chebyshev-Galerkin 方法

Chebyshev 多项式的权函数 $\omega=\left(1-x^{2}\right)^{-1 / 2}$.

令 $\phi_{k}(x)=T_{k}(x)+a_{k} T_{k+1}(x)+b_{k} T_{k+2}(x)$ 满足边界条件,[Shen2011, P146] 可得 $a_{k}=0$, $b_{k}=-1$,所以

设 $\{\phi_{k}\}_{k=0}^{N-2}$ 是 $X_{N}$ 的一组基函数,我们将逼近解展开为

谱格式:

将 $u_{N}(x)$ 的展开式代入谱格式方程。令

取 $v_N=\phi_{k}$,可得到以下方程组

质量矩阵 $M=(m_{jk})$ 是一个对称正定的五对角矩阵 [P149-4.29],它的非零元素为

其中 $c_0=2$,当 $k \geq 1$ 时,$c_k=1$.

刚度矩阵 $S=(s_{jk})$ 是一个上三角矩阵 [P149-4.30],其元素如下:

数值例子

对于前面提到的两点边值问题,取 $\alpha=1$,如果取精确解为 $\sin (k \pi x)$,则右端项为 $f=k^2 \pi^2 \sin(k\pi x)+\sin(k\pi x)$.

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
% ChebSM1.m
% Chebyshev spectral-Galerkin method for the model equation
% -u_xx+u=f in (-1,1)
% boundary condition: u(-1)=u(1)=0;
% exact solution: u=sin(kw*pi*x);
% RHS: f=kw*kw*pi^2*sin(kw*pi*x)+sin(kw*pi*x);
% Rmk: Use routines japoly(); jags(); japolym();
clear all; close all;
kw=10;
Nvec=32:2:76;
% Initialization for error
L2_Err=[]; Max_Err=[];
% Loop for various modes N to calculate numerical errors
for N=Nvec
[xv,wv]=jags(N+1,-1/2,-1/2); % Chebyshev-Gauss points and weights
Cm=japolym(N,-1/2,-1/2,xv)./japolym(N,-1/2,-1/2,1); % matrix of Chebyshev polynomals
u=sin(kw*pi*xv); % exact solution
f=kw*kw*pi^2*sin(kw*pi*xv)+sin(kw*pi*xv); % RHS
% Calculting coefficient matrix
S=zeros(N-1); % stiff matrix
for k=1:N-1
for j=1:N-1
if k==j
S(k,j)=2*pi*k*(k+1);
elseif (k<j && mod(j-k,2)==0)
S(k,j)=4*pi*k;
else
S(k,j)=0;
end
end
end

M=diag([3/2*pi, pi*ones(1,N-2)])+diag(-pi/2*ones(1,N-3),2)...
+diag(-pi/2*ones(1,N-3),-2); % mass matrix
A=S+M;

% Solving the linear system
Pm=Cm(1:end-2,:)-Cm(3:end,:); % matrix of Phi(x)
b=Pm*diag(wv)*f; % Solving RHS
uh=A\b; % expansion coefficients of u_N(x)
un=Pm'*uh; % compositing the numerical solution

L2_err=sqrt(((un-u).^2)'*wv); % L^2 error
Max_err=norm(abs(un-u),inf); % maximum pointwise error
L2_Err=[L2_Err;L2_err];
Max_Err=[Max_Err;Max_err];
end
% plot L^2 and maximum pointwise error
plot(Nvec,log10(L2_Err),'bo-','MarkerFaceColor','w','LineWidth',1)
hold on
plot(Nvec,log10(Max_Err),'rd-','MarkerFaceColor','w','LineWidth',1)
grid on
legend('L^2 error','L^{\infty} error','location','NorthEast')
% title('Convergence of Chebyshev-Galerkin method','fontsize',12)
set(gca,'fontsize',12)
xlabel('N','fontsize', 14), ylabel('log_{10}Error','fontsize',14)

% sets axis tick and axis limits
xticks(30:10:80)
yticks(-14:2:0)
xlim([30 80])
ylim([-14 0])

% print -dpng -r600 ChebSM1.png
% print -depsc2 ChebSM1.eps

ChebSM1.png

参考书籍

  • Spectral Methods: Algorithms, Analysis and Applications