comparison +sbp/+implementations/d1_noneq_minimal_10.m @ 423:a2cb0d4f4a02 feature/grids

Merge in default.
author Jonatan Werpers <jonatan@werpers.com>
date Tue, 07 Feb 2017 15:47:51 +0100
parents f7ac3cd6eeaa
children 4cb627c7fb90
comparison
equal deleted inserted replaced
218:da058ce66876 423:a2cb0d4f4a02
1 function [D1,H,x,h] = d1_noneq_minimal_10(N,L)
2
3 % L: Domain length
4 % N: Number of grid points
5 if(nargin < 2)
6 L = 1;
7 end
8
9 if(N<16)
10 error('Operator requires at least 16 grid points');
11 end
12
13 % BP: Number of boundary points
14 % m: Number of nonequidistant spacings
15 % order: Accuracy of interior stencil
16 BP = 8;
17 m = 3;
18 order = 10;
19
20 %%%% Non-equidistant grid points %%%%%
21 x0 = 0.0000000000000e+00;
22 x1 = 5.8556160757529e-01;
23 x2 = 1.7473267488572e+00;
24 x3 = 3.0000000000000e+00;
25 x4 = 4.0000000000000e+00;
26 x5 = 5.0000000000000e+00;
27 x6 = 6.0000000000000e+00;
28 x7 = 7.0000000000000e+00;
29 x8 = 8.0000000000000e+00;
30
31 xb = sparse(m+1,1);
32 for i = 0:m
33 xb(i+1) = eval(['x' num2str(i)]);
34 end
35 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
36
37 %%%% Compute h %%%%%%%%%%
38 h = L/(2*xb(end) + N-1-2*m);
39 %%%%%%%%%%%%%%%%%%%%%%%%%
40
41 %%%% Define grid %%%%%%%%
42 x = h*[xb; linspace(xb(end)+1,L/h-xb(end)-1,N-2*(m+1))'; L/h-flip(xb) ];
43 %%%%%%%%%%%%%%%%%%%%%%%%%
44
45 %%%% Norm matrix %%%%%%%%
46 P = sparse(BP,1);
47 %#ok<*NASGU>
48 P0 = 1.6717213975289e-01;
49 P1 = 9.3675739171278e-01;
50 P2 = 1.3035532379753e+00;
51 P3 = 1.1188461804303e+00;
52 P4 = 9.6664345922660e-01;
53 P5 = 1.0083235564392e+00;
54 P6 = 9.9858767377362e-01;
55 P7 = 1.0001163606893e+00;
56
57 for i = 0:BP-1
58 P(i+1) = eval(['P' num2str(i)]);
59 end
60
61 H = ones(N,1);
62 H(1:BP) = P;
63 H(end-BP+1:end) = flip(P);
64 H = spdiags(h*H,0,N,N);
65 %%%%%%%%%%%%%%%%%%%%%%%%%
66
67 %%%% Q matrix %%%%%%%%%%%
68
69 % interior stencil
70 switch order
71 case 2
72 d = [-1/2,0,1/2];
73 case 4
74 d = [1/12,-2/3,0,2/3,-1/12];
75 case 6
76 d = [-1/60,3/20,-3/4,0,3/4,-3/20,1/60];
77 case 8
78 d = [1/280,-4/105,1/5,-4/5,0,4/5,-1/5,4/105,-1/280];
79 case 10
80 d = [-1/1260,5/504,-5/84,5/21,-5/6,0,5/6,-5/21,5/84,-5/504,1/1260];
81 case 12
82 d = [1/5544,-1/385,1/56,-5/63,15/56,-6/7,0,6/7,-15/56,5/63,-1/56,1/385,-1/5544];
83 end
84 d = repmat(d,N,1);
85 Q = spdiags(d,-order/2:order/2,N,N);
86
87 % Boundaries
88 Q0_0 = -5.0000000000000e-01;
89 Q0_1 = 6.7349296966214e-01;
90 Q0_2 = -2.5186401896559e-01;
91 Q0_3 = 8.3431385420901e-02;
92 Q0_4 = 2.5480326895984e-02;
93 Q0_5 = -4.5992420658252e-02;
94 Q0_6 = 1.7526412909003e-02;
95 Q0_7 = -2.0746552641799e-03;
96 Q0_8 = 0.0000000000000e+00;
97 Q0_9 = 0.0000000000000e+00;
98 Q0_10 = 0.0000000000000e+00;
99 Q0_11 = 0.0000000000000e+00;
100 Q0_12 = 0.0000000000000e+00;
101 Q1_0 = -6.7349296966214e-01;
102 Q1_1 = 0.0000000000000e+00;
103 Q1_2 = 9.1982892384044e-01;
104 Q1_3 = -2.7262271754043e-01;
105 Q1_4 = -5.0992113348238e-02;
106 Q1_5 = 1.1814647281129e-01;
107 Q1_6 = -4.6693123378079e-02;
108 Q1_7 = 5.8255272771571e-03;
109 Q1_8 = 0.0000000000000e+00;
110 Q1_9 = 0.0000000000000e+00;
111 Q1_10 = 0.0000000000000e+00;
112 Q1_11 = 0.0000000000000e+00;
113 Q1_12 = 0.0000000000000e+00;
114 Q2_0 = 2.5186401896559e-01;
115 Q2_1 = -9.1982892384044e-01;
116 Q2_2 = 0.0000000000000e+00;
117 Q2_3 = 7.8566746772741e-01;
118 Q2_4 = -2.4097806629929e-02;
119 Q2_5 = -1.5312168858669e-01;
120 Q2_6 = 6.9451518963875e-02;
121 Q2_7 = -9.9345865998262e-03;
122 Q2_8 = 0.0000000000000e+00;
123 Q2_9 = 0.0000000000000e+00;
124 Q2_10 = 0.0000000000000e+00;
125 Q2_11 = 0.0000000000000e+00;
126 Q2_12 = 0.0000000000000e+00;
127 Q3_0 = -8.3431385420901e-02;
128 Q3_1 = 2.7262271754043e-01;
129 Q3_2 = -7.8566746772741e-01;
130 Q3_3 = 0.0000000000000e+00;
131 Q3_4 = 6.2047871210535e-01;
132 Q3_5 = 1.4776775176509e-02;
133 Q3_6 = -4.6889652372990e-02;
134 Q3_7 = 7.3166499053672e-03;
135 Q3_8 = 7.9365079365079e-04;
136 Q3_9 = 0.0000000000000e+00;
137 Q3_10 = 0.0000000000000e+00;
138 Q3_11 = 0.0000000000000e+00;
139 Q3_12 = 0.0000000000000e+00;
140 Q4_0 = -2.5480326895984e-02;
141 Q4_1 = 5.0992113348238e-02;
142 Q4_2 = 2.4097806629929e-02;
143 Q4_3 = -6.2047871210535e-01;
144 Q4_4 = 0.0000000000000e+00;
145 Q4_5 = 6.9425006383507e-01;
146 Q4_6 = -1.5686345740485e-01;
147 Q4_7 = 4.2609496719925e-02;
148 Q4_8 = -9.9206349206349e-03;
149 Q4_9 = 7.9365079365079e-04;
150 Q4_10 = 0.0000000000000e+00;
151 Q4_11 = 0.0000000000000e+00;
152 Q4_12 = 0.0000000000000e+00;
153 Q5_0 = 4.5992420658252e-02;
154 Q5_1 = -1.1814647281129e-01;
155 Q5_2 = 1.5312168858669e-01;
156 Q5_3 = -1.4776775176509e-02;
157 Q5_4 = -6.9425006383507e-01;
158 Q5_5 = 0.0000000000000e+00;
159 Q5_6 = 8.0719535654891e-01;
160 Q5_7 = -2.2953297936781e-01;
161 Q5_8 = 5.9523809523809e-02;
162 Q5_9 = -9.9206349206349e-03;
163 Q5_10 = 7.9365079365079e-04;
164 Q5_11 = 0.0000000000000e+00;
165 Q5_12 = 0.0000000000000e+00;
166 Q6_0 = -1.7526412909003e-02;
167 Q6_1 = 4.6693123378079e-02;
168 Q6_2 = -6.9451518963875e-02;
169 Q6_3 = 4.6889652372990e-02;
170 Q6_4 = 1.5686345740485e-01;
171 Q6_5 = -8.0719535654891e-01;
172 Q6_6 = 0.0000000000000e+00;
173 Q6_7 = 8.3142546796428e-01;
174 Q6_8 = -2.3809523809524e-01;
175 Q6_9 = 5.9523809523809e-02;
176 Q6_10 = -9.9206349206349e-03;
177 Q6_11 = 7.9365079365079e-04;
178 Q6_12 = 0.0000000000000e+00;
179 Q7_0 = 2.0746552641799e-03;
180 Q7_1 = -5.8255272771571e-03;
181 Q7_2 = 9.9345865998262e-03;
182 Q7_3 = -7.3166499053672e-03;
183 Q7_4 = -4.2609496719925e-02;
184 Q7_5 = 2.2953297936781e-01;
185 Q7_6 = -8.3142546796428e-01;
186 Q7_7 = 0.0000000000000e+00;
187 Q7_8 = 8.3333333333333e-01;
188 Q7_9 = -2.3809523809524e-01;
189 Q7_10 = 5.9523809523809e-02;
190 Q7_11 = -9.9206349206349e-03;
191 Q7_12 = 7.9365079365079e-04;
192 for i = 1:BP
193 for j = 1:BP
194 Q(i,j) = eval(['Q' num2str(i-1) '_' num2str(j-1)]);
195 Q(N+1-i,N+1-j) = -eval(['Q' num2str(i-1) '_' num2str(j-1)]);
196 end
197 end
198 %%%%%%%%%%%%%%%%%%%%%%%%%%%
199
200 %%%% Difference operator %%
201 D1 = H\Q;
202 %%%%%%%%%%%%%%%%%%%%%%%%%%%