OptFEM2D Toolbox for Matlab  V1.2b1
Matlab/Octave Optimized P1-Lagrange Finite Element Method in 2D
 All Files Functions Pages
validStiffElasP1.m
Go to the documentation of this file.
1 function validStiffElasP1(Num)
2 % function validStiffElasP1(Num)
3 % Validation function for the Assembly of the Stiffness Elasticity
4 % Matrix for `P_1`-Lagrange finite elements.
5 %
6 % The Stiffness Elasticity Matrix is given by
7 % ``\StiffElas_{m,l}=\int_{\DOMH} \Odv^t(\BasisFuncTwoD_m) \Ocv(\BasisFuncTwoD_l)dT, \ \forall (m,l)\in\ENS{1}{2\,\nq}^2,``
8 % where `\BasisFuncTwoD_m` are vectors field `P_1`-Lagrange basis functions.
9 % Here `\Ocv=(\Occ_{xx},\Occ_{yy},\Occ_{xy})^t` and `\Odv=(\Odc_{xx},\Odc_{yy},2\Odc_{xy})^t`
10 % are the elastic stress and strain tensors respectively.
11 %
12 % This Matrix is computed by functions StiffElasAssemblingP1{Version} where {Version} is one of
13 % 'base', 'OptV0', 'OptV1' and 'OptV2'.
14 % - Test 1: computation of the Elastic Stiffness Matrix using previous functions giving errors and cputimes
15 % - Test 2: compute ``\int_\DOM \Odv^t(\vecb{u}) \Ocv(\vecb{v})d\q \approx \DOT{\StiffElas \vecb{U}}{\vecb{V}}``
16 % with `\foncdefsmall{\vecb{u}=(u_1,u_2)}{\DOM}{\R^2}` and `\foncdefsmall{\vecb{v}=(v_1,v_2)}{\DOM}{\R^2}`.
17 % We have `\vecb{U}_{2i-1}=u_1(\q^i),` `\vecb{U}_{2i}=u_2(\q^i),` `\vecb{V}_{2i-1}=v_1(\q^i),` `\vecb{V}_{2i}=v_2(\q^i)`
18 % if Num in {0,2} and `\vecb{U}_{i}=u_1(\q^i),` `\vecb{U}_{i+\nq}=u_2(\q^i),` `\vecb{V}_{i}=v_1(\q^i),` `\vecb{V}_{i+\nq}=v_2(\q^i)`
19 % if Num in {1,3}. Use functions `\vecb{u}` and `\vecb{u}` defined in #valid_StiffElas.
20 % - Test 3: retrieve order 2 of `P_1`-Lagrange integration
21 % ``|\int_\DOM \Odv^t(\vecb{u}) \Ocv(\vecb{v})d\q - \DOT{\StiffElas \vecb{U}}{\vecb{V}}| \leq C h^2``
22 %
23 % Num :
24 % - 0 global alternate numbering with local alternate numbering (classical method),
25 % - 1 global block numbering with local alternate numbering,
26 % - 2 global alternate numbering with local block numbering,
27 % - 3 global block numbering with local block numbering.
28 % See also:
29 % #StiffElasAssemblingP1base, #StiffElasAssemblingP1OptV0,
30 % #StiffElasAssemblingP1OptV1, #StiffElasAssemblingP1OptV2, #SquareMesh
31 %
32 % @author François Cuvelier @date 2012-11-26
33 %
34 % OptFEM2DP1 [V1.2b1] - Copyright (C) 2013 CJS (LAGA)
35 %
36 % This file is part of OptFEM2DP1.
37 % OptFEM2DP1 is free software: you can redistribute it and/or modify
38 % it under the terms of the GNU General Public License as published by
39 % the Free Software Foundation, either version 3 of the License, or
40 % (at your option) any later version.
41 %
42 % OptFEM2DP1 is distributed in the hope that it will be useful,
43 % but WITHOUT ANY WARRANTY; without even the implied warranty of
44 % MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
45 % GNU General Public License for more details.
46 %
47 % You should have received a copy of the GNU General Public License
48 % along with this program. If not, see <http://www.gnu.org/licenses/>.
49 
50  disp('********************************************')
51  disp('* Stiff Elas Assembling P1 validations *')
52  disp('********************************************')
53  switch Num
54  case 0
55  s=sprintf('Global alternate numbering / local alternate numbering');
56  case 1
57  s=sprintf('Global block numbering / local alternate numbering');
58  case 2
59  s=sprintf('Global alternate numbering / local block numbering');
60  case 3
61  s=sprintf('Global block numbering / local block numbering');
62  otherwise
63  error('invalid Num value')
64  end
65  fprintf(' Numbering Choice : %s\n',s);
66 
67 
68 
69  Th=SquareMesh(20);
70 
71 % TEST 1
72  disp('-----------------------------------------')
73  disp(' Test 1: Matrices errors and CPU times ')
74  disp('-----------------------------------------')
75  tic();
76  Mbase=StiffElasAssemblingP1base(Th.nq,Th.nme,Th.q,Th.me,Th.areas,1,1,Num);
77  T(1)=toc();
78  tic();
79  MOptV0=StiffElasAssemblingP1OptV0(Th.nq,Th.nme,Th.q,Th.me,Th.areas,1,1,Num);
80  T(2)=toc();
81  Test1.error(1)=norm(Mbase-MOptV0,Inf);
82  Test1.name{1}='StiffElasAssemblingP1OptV0';
83  fprintf(' Error P1base vs OptV0 : %e\n',Test1.error(1))
84  tic();
85  MOptV1=StiffElasAssemblingP1OptV1(Th.nq,Th.nme,Th.q,Th.me,Th.areas,1,1,Num);
86  T(3)=toc();
87  Test1.error(2)=norm(Mbase-MOptV1,Inf);
88  Test1.name{2}='StiffElasAssemblingP1OptV1';
89  fprintf(' Error P1base vs OptV1 : %e\n',Test1.error(2))
90  tic();
91  MOptV2=StiffElasAssemblingP1OptV2(Th.nq,Th.nme,Th.q,Th.me,Th.areas,1,1,Num);
92  T(4)=toc();
93  Test1.error(3)=norm(Mbase-MOptV2,Inf);
94  Test1.name{3}='StiffElasAssemblingP1OptV2';
95  fprintf(' Error P1base vs OptV2 : %e\n',Test1.error(3))
96 
97  fprintf(' CPU times base (ref) : %3.4f (s)\n',T(1))
98  fprintf(' CPU times OptV0 : %3.4f (s) - Speed Up X%3.3f\n',T(2),T(1)/T(2))
99  fprintf(' CPU times OptV1 : %3.4f (s) - Speed Up X%3.3f\n',T(3),T(1)/T(3))
100  fprintf(' CPU times OptV2 : %3.4f (s) - Speed Up X%3.3f\n',T(4),T(1)/T(4))
101  checkTest1(Test1)
102 
103 
104 
105 
106 % TEST 2
107  disp('-----------------------------------------------------')
108  disp(' Test 2: Validations by integration on [0,1]x[0,1] ')
109  disp('-----------------------------------------------------')
110  Test=valid_StiffElas();
111  qf=Th.q;
112  for kk=1:length(Test)
113  KK=StiffElasAssemblingP1OptV2(Th.nq,Th.nme,Th.q,Th.me,Th.areas,Test(kk).lambda,Test(kk).mu,Num);
114  u=Test(kk).u;
115  v=Test(kk).v;
116  switch Num
117  case {0,2}
118  U=[u{1}(qf(1,:),qf(2,:));u{2}(qf(1,:),qf(2,:))];
119  U=U(:);
120  V=[v{1}(qf(1,:),qf(2,:));v{2}(qf(1,:),qf(2,:))];
121  V=V(:);
122  case {1,3}
123  U=[u{1}(qf(1,:),qf(2,:)) u{2}(qf(1,:),qf(2,:))]';
124  V=[v{1}(qf(1,:),qf(2,:)) v{2}(qf(1,:),qf(2,:))]';
125  end
126  %whos
127  Test(kk).error=abs(Test(kk).Stiff-U'*KK*V);
128  fprintf(' functions %d : u(x,y)=(%s,%s), v(x,y)=(%s,%s),\n -> Stiff error=%e\n', ...
129  kk,Test(kk).cu{1},Test(kk).cu{2},Test(kk).cv{1},Test(kk).cv{2},Test(kk).error);
130  end
131  checkTest2(Test)
132 
133 % TEST 3
134  disp('--------------------------------')
135  disp(' Test 3: Validations by order ')
136  disp('--------------------------------')
137  n=length(Test);
138  u=Test(n).u;
139  v=Test(n).v;
140  lambda=Test(n).lambda;
141  mu=Test(n).mu;
142  ExSol=Test(n).Stiff;
143 
144  for k=1:10
145  Th=SquareMesh(20*k+50);
146  qf=Th.q;
147  fprintf(' Matrix size : %d\n',Th.nq);
148  h(k)=GetMaxLengthEdges(Th.q,Th.me);
149  tic();
150  M=StiffElasAssemblingP1OptV2(Th.nq,Th.nme,Th.q,Th.me,Th.areas,lambda,mu,Num);
151  TT(k)=toc();
152  switch Num
153  case {0,2}
154  U=[u{1}(qf(1,:),qf(2,:));u{2}(qf(1,:),qf(2,:))];
155  U=U(:);
156  V=[v{1}(qf(1,:),qf(2,:));v{2}(qf(1,:),qf(2,:))];
157  V=V(:);
158  case {1,3}
159  U=[u{1}(qf(1,:),qf(2,:)) u{2}(qf(1,:),qf(2,:))]';
160  V=[v{1}(qf(1,:),qf(2,:)) v{2}(qf(1,:),qf(2,:))]';
161  end
162  Error(k)=abs(ExSol-U'*M*V);
163  fprintf(' StiffElasAssemblingP1OptV2 CPU times : %3.3f(s)\n',TT(k));
164  fprintf(' Error : %e\n',Error(k));
165  end
166 
167  figure(1)
168  loglog(h,Error,'+-r',h,h*1.1*Error(1)/h(1),'-sm',h,1.1*Error(1)*(h/h(1)).^2,'-db')
169  legend('Error','O(h)','O(h^2)')
170  xlabel('h')
171  title(sprintf('Test 3 : Stiff Elas. Matrix (Num=%d)',Num))
172  checkTest3(h,Error)
173 
174  %figure(3)
175  %spy(M)
176  %title(sprintf('Matrix sparsity for %s numbering',s))
177 end
178 
179 function checkTest1(Test)
180  I=find(Test.error>1e-14);
181  if isempty(I)
182  disp('------------------------')
183  disp(' Test 1 (results): OK')
184  disp('------------------------')
185  else
186  disp('----------------------------')
187  disp(' Test 1 (results): FAILED')
188  disp('----------------------------')
189  end
190 end
191 
192 function checkTest2(Test)
193  N=length(Test);
194  cntFalse=0;
195  for k=1:N
196  if (Test(k).degree<=1)
197  if (Test(k).error>1e-12)
198  cntFalse=cntFalse+1;
199  end
200  end
201  end
202  if (cntFalse==0)
203  disp('------------------------')
204  disp(' Test 2 (results): OK')
205  disp('------------------------')
206  else
207  disp('----------------------------')
208  disp(' Test 2 (results): FAILED')
209  disp('----------------------------')
210  end
211 end
212 
213 function checkTest3(h,error)
214  % order 2
215  P=polyfit(log(h),log(error),1);
216  if abs(P(1)-2)<1e-3
217  disp('------------------------')
218  disp(' Test 3 (results): OK')
219  fprintf(' -> found numerical order %f. Must be 2\n',P(1))
220  disp('------------------------')
221  else
222  disp('----------------------------')
223  disp(' Test 3 (results): FAILED')
224  fprintf(' -> found numerical order %f. Must be 2\n',P(1))
225  disp('----------------------------')
226  end
227 end