OptFEM2DP1 Toolbox  V1.2b3
Matlab/Octave Optimized P1-Lagrange Finite Element Method in 2D
 All Files Functions Pages
StiffElasAssemblingP1OptV1.m
Go to the documentation of this file.
1 function K=StiffElasAssemblingP1OptV1(nq,nme,q,me,areas,lambda,mu,Num)
2 % function K=StiffElasAssemblingP1OptV1(nq,nme,q,me,areas,lambda,mu,Num)
3 % Assembly of the Stiffness Elasticity Matrix by `P_1`-Lagrange finite elements
4 % using "OptV1" version (see report).
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 `P_1`-Lagrange vector 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 % Parameters:
13 % nq: total number of nodes in the mesh, also denoted by `\nq`.
14 % nme: total number of triangles, also denoted by `\nme`.
15 % q: Array of vertices coordinates, `2\times\nq` array. <br/>
16 % `{\q}(\il,j)` is the
17 % `\il`-th coordinate of the `j`-th vertex, `\il\in\{1,2\}` and
18 % `j\in\ENS{1}{\nq}`
19 % me: Connectivity array, `3\times\nme` array. <br/>
20 % `\me(\jl,k)` is the storage index of the
21 % `\jl`-th vertex of the `k`-th triangle in the array `\q`, `\jl\in\{1,2,3\}` and
22 % `k\in{\ENS{1}{\nme}}`.
23 % areas: Array of areas, `1\times\nme` array. areas(k) is the area of the `k`-th triangle.
24 % lambda: the first Lame coefficient in Hooke's law
25 % mu: the second Lame coefficient in Hooke's law
26 % Num:
27 % - 0 global alternate numbering with local alternate numbering (classical method),
28 % - 1 global block numbering with local alternate numbering,
29 % - 2 global alternate numbering with local block numbering,
30 % - 3 global block numbering with local block numbering.
31 %
32 % Return values:
33 % K: `2\nq\times 2\nq` stiffness elasticity sparse matrix
34 %
35 % Example:
36 % @verbatim
37 % Th=SquareMesh(10);
38 % lambda=1; mu=1;
39 % Num = 0;
40 % KK=StiffElasAssemblingP1OptV1(Th.nq,Th.nme,Th.q,Th.me,Th.areas,lambda,mu,Num);@endverbatim
41 %
42 % See also:
43 % #BuildIkFunc, #BuildElemStiffElasMatOptV0Func
44 % Copyright:
45 % See \ref license
46 ElemStiffElasMatOptV0=BuildElemStiffElasMatOptV0Func(Num);
47 C=[lambda+2*mu,lambda,0;lambda,lambda+2*mu,0;0,0,mu];
48 GetI=BuildIkFunc(Num,nq);
49 Ig=zeros(36*nme,1);Jg=zeros(36*nme,1);Kg=zeros(36*nme,1);
50 kk=1:36;
51 for k=1:nme
52  Me=ElemStiffElasMatOptV0(q(:,me(:,k)),areas(k),C);
53  I=GetI(me,k);
54  jA=ones(6,1)*I;
55  iA=jA';
56 
57  Ig(kk)=iA(:);
58  Jg(kk)=jA(:);
59  Kg(kk)=Me(:);
60  kk=kk+36;
61 end
62 K=sparse(Ig,Jg,Kg,2*nq,2*nq);