comparison +scheme/LaplaceCurvilinear.m @ 704:111fcbcff2e9 feature/optim

merg with featuew grids
author Ylva Rydin <ylva.rydin@telia.com>
date Fri, 03 Nov 2017 10:53:15 +0100
parents 33b962620e24
children 07f8311374c6
comparison
equal deleted inserted replaced
703:027f606fa691 704:111fcbcff2e9
1 classdef LaplaceCurvilinear < scheme.Scheme
2 properties
3 m % Number of points in each direction, possibly a vector
4 h % Grid spacing
5
6 grid
7
8 order % Order accuracy for the approximation
9
10 a,b % Parameters of the operator
11
12
13 % Inner products and operators for physical coordinates
14 D % Laplace operator
15 H, Hi % Inner product
16 e_w, e_e, e_s, e_n
17 d_w, d_e, d_s, d_n % Normal derivatives at the boundary
18 H_w, H_e, H_s, H_n % Boundary inner products
19 Dx, Dy % Physical derivatives
20 M % Gradient inner product
21
22 % Metric coefficients
23 J, Ji
24 a11, a12, a22
25 x_u
26 x_v
27 y_u
28 y_v
29
30 % Inner product and operators for logical coordinates
31 H_u, H_v % Norms in the x and y directions
32 Hi_u, Hi_v
33 Hu,Hv % Kroneckerd norms. 1'*Hx*v corresponds to integration in the x dir.
34 Hiu, Hiv
35 du_w, dv_w
36 du_e, dv_e
37 du_s, dv_s
38 du_n, dv_n
39 gamm_u, gamm_v
40 lambda
41 end
42
43 methods
44 % Implements a*div(b*grad(u)) as a SBP scheme
45 % TODO: Implement proper H, it should be the real physical quadrature, the logic quadrature may be but in a separate variable (H_logic?)
46
47 function obj = LaplaceCurvilinear(g ,order, a, b, opSet)
48 default_arg('opSet',@sbp.D2Variable);
49 default_arg('a', 1);
50 default_arg('b', 1);
51
52 if b ~=1
53 error('Not implemented yet')
54 end
55
56 assert(isa(g, 'grid.Curvilinear'))
57
58 m = g.size();
59 m_u = m(1);
60 m_v = m(2);
61 m_tot = g.N();
62
63 h = g.scaling();
64 h_u = h(1);
65 h_v = h(2);
66
67
68 % 1D operators
69 ops_u = opSet(m_u, {0, 1}, order);
70 ops_v = opSet(m_v, {0, 1}, order);
71
72 I_u = speye(m_u);
73 I_v = speye(m_v);
74
75 D1_u = ops_u.D1;
76 D2_u = ops_u.D2;
77 H_u = ops_u.H;
78 Hi_u = ops_u.HI;
79 e_l_u = ops_u.e_l;
80 e_r_u = ops_u.e_r;
81 d1_l_u = ops_u.d1_l;
82 d1_r_u = ops_u.d1_r;
83
84 D1_v = ops_v.D1;
85 D2_v = ops_v.D2;
86 H_v = ops_v.H;
87 Hi_v = ops_v.HI;
88 e_l_v = ops_v.e_l;
89 e_r_v = ops_v.e_r;
90 d1_l_v = ops_v.d1_l;
91 d1_r_v = ops_v.d1_r;
92
93
94 % Logical operators
95 Du = kr(D1_u,I_v);
96 Dv = kr(I_u,D1_v);
97 obj.Hu = kr(H_u,I_v);
98 obj.Hv = kr(I_u,H_v);
99 obj.Hiu = kr(Hi_u,I_v);
100 obj.Hiv = kr(I_u,Hi_v);
101
102 e_w = kr(e_l_u,I_v);
103 e_e = kr(e_r_u,I_v);
104 e_s = kr(I_u,e_l_v);
105 e_n = kr(I_u,e_r_v);
106 obj.du_w = kr(d1_l_u,I_v);
107 obj.dv_w = (e_w'*Dv)';
108 obj.du_e = kr(d1_r_u,I_v);
109 obj.dv_e = (e_e'*Dv)';
110 obj.du_s = (e_s'*Du)';
111 obj.dv_s = kr(I_u,d1_l_v);
112 obj.du_n = (e_n'*Du)';
113 obj.dv_n = kr(I_u,d1_r_v);
114
115
116 % Metric coefficients
117 coords = g.points();
118 x = coords(:,1);
119 y = coords(:,2);
120
121 x_u = Du*x;
122 x_v = Dv*x;
123 y_u = Du*y;
124 y_v = Dv*y;
125
126 J = x_u.*y_v - x_v.*y_u;
127 a11 = 1./J .* (x_v.^2 + y_v.^2);
128 a12 = -1./J .* (x_u.*x_v + y_u.*y_v);
129 a22 = 1./J .* (x_u.^2 + y_u.^2);
130 lambda = 1/2 * (a11 + a22 - sqrt((a11-a22).^2 + 4*a12.^2));
131
132 obj.x_u = x_u;
133 obj.x_v = x_v;
134 obj.y_u = y_u;
135 obj.y_v = y_v;
136
137
138 % Assemble full operators
139 L_12 = spdiag(a12);
140 Duv = Du*L_12*Dv;
141 Dvu = Dv*L_12*Du;
142
143 Duu = sparse(m_tot);
144 Dvv = sparse(m_tot);
145 ind = grid.funcToMatrix(g, 1:m_tot);
146
147 for i = 1:m_v
148 D = D2_u(a11(ind(:,i)));
149 p = ind(:,i);
150 Duu(p,p) = D;
151 end
152
153 for i = 1:m_u
154 D = D2_v(a22(ind(i,:)));
155 p = ind(i,:);
156 Dvv(p,p) = D;
157 end
158
159
160 % Physical operators
161 obj.J = spdiag(J);
162 obj.Ji = spdiag(1./J);
163
164 obj.D = obj.Ji*a*(Duu + Duv + Dvu + Dvv);
165 obj.H = obj.J*kr(H_u,H_v);
166 obj.Hi = obj.Ji*kr(Hi_u,Hi_v);
167
168 obj.e_w = e_w;
169 obj.e_e = e_e;
170 obj.e_s = e_s;
171 obj.e_n = e_n;
172
173 %% normal derivatives
174 I_w = ind(1,:);
175 I_e = ind(end,:);
176 I_s = ind(:,1);
177 I_n = ind(:,end);
178
179 a11_w = spdiag(a11(I_w));
180 a12_w = spdiag(a12(I_w));
181 a11_e = spdiag(a11(I_e));
182 a12_e = spdiag(a12(I_e));
183 a22_s = spdiag(a22(I_s));
184 a12_s = spdiag(a12(I_s));
185 a22_n = spdiag(a22(I_n));
186 a12_n = spdiag(a12(I_n));
187
188 s_w = sqrt((e_w'*x_v).^2 + (e_w'*y_v).^2);
189 s_e = sqrt((e_e'*x_v).^2 + (e_e'*y_v).^2);
190 s_s = sqrt((e_s'*x_u).^2 + (e_s'*y_u).^2);
191 s_n = sqrt((e_n'*x_u).^2 + (e_n'*y_u).^2);
192
193 obj.d_w = -1*(spdiag(1./s_w)*(a11_w*obj.du_w' + a12_w*obj.dv_w'))';
194 obj.d_e = (spdiag(1./s_e)*(a11_e*obj.du_e' + a12_e*obj.dv_e'))';
195 obj.d_s = -1*(spdiag(1./s_s)*(a22_s*obj.dv_s' + a12_s*obj.du_s'))';
196 obj.d_n = (spdiag(1./s_n)*(a22_n*obj.dv_n' + a12_n*obj.du_n'))';
197
198 obj.Dx = spdiag( y_v./J)*Du + spdiag(-y_u./J)*Dv;
199 obj.Dy = spdiag(-x_v./J)*Du + spdiag( x_u./J)*Dv;
200
201 %% Boundary inner products
202 obj.H_w = H_v*spdiag(s_w);
203 obj.H_e = H_v*spdiag(s_e);
204 obj.H_s = H_u*spdiag(s_s);
205 obj.H_n = H_u*spdiag(s_n);
206
207 % Misc.
208 obj.m = m;
209 obj.h = [h_u h_v];
210 obj.order = order;
211 obj.grid = g;
212
213 obj.a = a;
214 obj.b = b;
215 obj.a11 = a11;
216 obj.a12 = a12;
217 obj.a22 = a22;
218 obj.lambda = lambda;
219
220 obj.gamm_u = h_u*ops_u.borrowing.M.d1;
221 obj.gamm_v = h_v*ops_v.borrowing.M.d1;
222 end
223
224
225 % Closure functions return the opertors applied to the own doamin to close the boundary
226 % Penalty functions return the opertors to force the solution. In the case of an interface it returns the operator applied to the other doamin.
227 % boundary is a string specifying the boundary e.g. 'l','r' or 'e','w','n','s'.
228 % type is a string specifying the type of boundary condition if there are several.
229 % data is a function returning the data that should be applied at the boundary.
230 % neighbour_scheme is an instance of Scheme that should be interfaced to.
231 % neighbour_boundary is a string specifying which boundary to interface to.
232 function [closure, penalty] = boundary_condition(obj, boundary, type, parameter)
233 default_arg('type','neumann');
234 default_arg('parameter', []);
235
236 [e, d, gamm, H_b, ~] = obj.get_boundary_ops(boundary);
237 switch type
238 % Dirichlet boundary condition
239 case {'D','d','dirichlet'}
240 tuning = 1.2;
241 % tuning = 20.2;
242
243 b1 = gamm*obj.lambda./obj.a11.^2;
244 b2 = gamm*obj.lambda./obj.a22.^2;
245
246 tau1 = tuning * spdiag(-1./b1 - 1./b2);
247 tau2 = 1;
248
249 tau = (tau1*e + tau2*d)*H_b;
250
251 closure = obj.a*obj.Hi*tau*e';
252 penalty = -obj.a*obj.Hi*tau;
253
254
255 % Neumann boundary condition
256 case {'N','n','neumann'}
257 tau1 = -1;
258 tau2 = 0;
259 tau = (tau1*e + tau2*d)*H_b;
260
261 closure = obj.a*obj.Hi*tau*d';
262 penalty = -obj.a*obj.Hi*tau;
263
264
265 % Unknown, boundary condition
266 otherwise
267 error('No such boundary condition: type = %s',type);
268 end
269 end
270
271 function [closure, penalty] = interface(obj,boundary,neighbour_scheme,neighbour_boundary)
272 % u denotes the solution in the own domain
273 % v denotes the solution in the neighbour domain
274 tuning = 1.2;
275 % tuning = 20.2;
276 [e_u, d_u, gamm_u, H_b_u, I_u] = obj.get_boundary_ops(boundary);
277 [e_v, d_v, gamm_v, H_b_v, I_v] = neighbour_scheme.get_boundary_ops(neighbour_boundary);
278
279 u = obj;
280 v = neighbour_scheme;
281
282 b1_u = gamm_u*u.lambda(I_u)./u.a11(I_u).^2;
283 b2_u = gamm_u*u.lambda(I_u)./u.a22(I_u).^2;
284 b1_v = gamm_v*v.lambda(I_v)./v.a11(I_v).^2;
285 b2_v = gamm_v*v.lambda(I_v)./v.a22(I_v).^2;
286
287 tau1 = -1./(4*b1_u) -1./(4*b1_v) -1./(4*b2_u) -1./(4*b2_v);
288 tau1 = tuning * spdiag(tau1);
289 tau2 = 1/2;
290
291 sig1 = -1/2;
292 sig2 = 0;
293
294 tau = (e_u*tau1 + tau2*d_u)*H_b_u;
295 sig = (sig1*e_u + sig2*d_u)*H_b_u;
296
297 closure = obj.a*obj.Hi*( tau*e_u' + sig*d_u');
298 penalty = obj.a*obj.Hi*(-tau*e_v' + sig*d_v');
299 end
300
301 % Ruturns the boundary ops and sign for the boundary specified by the string boundary.
302 % The right boundary is considered the positive boundary
303 %
304 % I -- the indecies of the boundary points in the grid matrix
305 function [e, d, gamm, H_b, I] = get_boundary_ops(obj, boundary)
306
307 % gridMatrix = zeros(obj.m(2),obj.m(1));
308 % gridMatrix(:) = 1:numel(gridMatrix);
309
310 ind = grid.funcToMatrix(obj.grid, 1:prod(obj.m));
311
312 switch boundary
313 case 'w'
314 e = obj.e_w;
315 d = obj.d_w;
316 H_b = obj.H_w;
317 I = ind(1,:);
318 case 'e'
319 e = obj.e_e;
320 d = obj.d_e;
321 H_b = obj.H_e;
322 I = ind(end,:);
323 case 's'
324 e = obj.e_s;
325 d = obj.d_s;
326 H_b = obj.H_s;
327 I = ind(:,1)';
328 case 'n'
329 e = obj.e_n;
330 d = obj.d_n;
331 H_b = obj.H_n;
332 I = ind(:,end)';
333 otherwise
334 error('No such boundary: boundary = %s',boundary);
335 end
336
337 switch boundary
338 case {'w','e'}
339 gamm = obj.gamm_u;
340 case {'s','n'}
341 gamm = obj.gamm_v;
342 end
343 end
344
345 function N = size(obj)
346 N = prod(obj.m);
347 end
348 end
349 end