-
Notifications
You must be signed in to change notification settings - Fork 2
/
Copy pathTSDT_5dof.m
218 lines (190 loc) · 8.7 KB
/
TSDT_5dof.m
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
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
function [tangent,stiffness,force] = TSDT_5dof(Parameter,Load,Lx,Ly,h,...
nel,nnel,ndof,nnode,edof,sdof,Coordinate,IEN,...
B1Matrix,B2Matrix,B3Matrix,B4Matrix,...
S1Matrix,S2Matrix,S3Matrix,S4Matrix,displacement)
%--------------------------------------------------------------------------
% Order of Gauss Quadrature
%--------------------------------------------------------------------------
% For bending stiffness
nglb = 3; % Full Gauss-Legendre quadrature for bending
[pointb,weightb]=GaussQuadrature(nglb); % sampling points & weights
% For shear stiffness
ngls = 2; % Reduced Gauss-Legendre quadrature for shear
[points,weights] = GaussQuadrature(ngls); % sampling points & weights
%--------------------------------------------------------------------------
% Initialization of matrices and vectors
%--------------------------------------------------------------------------
tangent = zeros(sdof); % system tangent stiffness matrix
stiffness = zeros(sdof); % system stiffness matrix
force = zeros(sdof,1); % System Force Vector
for iel=1:nel % loop for the total number of elements
tSigmab_e = zeros(edof,edof); % initialization of initial stress matrix
kLLb_e = zeros(edof,edof); % initialization of linear bending matrix
kLNLb_e = zeros(edof,edof);
kNLLb_e = zeros(edof,edof);
kNLNLb_e = zeros(edof,edof);
tSigmas_e = zeros(edof,edof);
kLLs_e = zeros(edof,edof); % initialization of linear shear matrix
kLNLs_e = zeros(edof,edof);
kNLLs_e = zeros(edof,edof);
kNLNLs_e = zeros(edof,edof);
kp_e = zeros(edof,edof);
fm_e = zeros(edof,1); % initialization of force vector
% keyboard
IEN_e = IEN(iel,:);
x_e = Coordinate(IEN_e,1)';
y_e = Coordinate(IEN_e,2)';
u_e = displacement(IEN_e*ndof-(ndof-1))';
v_e = displacement(IEN_e*ndof-(ndof-2))';
w_e = displacement(IEN_e*ndof-(ndof-3))';
phix_e = displacement(IEN_e*ndof-(ndof-4))';
phiy_e = displacement(IEN_e*ndof-(ndof-5))';
thetax_e = displacement(IEN_e*ndof-(ndof-6))';
thetay_e = displacement(IEN_e*ndof-(ndof-7))';
%--------------------------------------------------------------------------
% Numerical integration for bending term
%--------------------------------------------------------------------------
for intx=1:nglb
xi=pointb(intx,1); % sampling point in x-axis
wtx=weightb(intx,1); % weight in x-axis
for inty=1:nglb
eta=pointb(inty,2); % sampling point in y-axis
wty=weightb(inty,2) ; % weight in y-axis
% keyboard
% Compute Shape Function
[R,dRdxi,dRdeta] = ShapeFunctionQ9(xi,eta);
% compute Jacobian
[detJ22,invJ22] = Jacobian(dRdxi,dRdeta,x_e,y_e);
% derivatives w.r.t. physical coordinate
[dRdx,dRdy] = ShapefunctionDerivatives(dRdxi,dRdeta,invJ22);
%--------------------------------------------------------------------------
% compute bending element matrix
%--------------------------------------------------------------------------
% linear bending kinematic matrix
BLb = LinearBending(ndof,nnel,B1Matrix,dRdx,dRdy);
kLLb_e = kLLb_e + BLb'*B1Matrix*BLb*wtx*wty*detJ22;
% % non-linear bending kinematic matrix
% [Ab,Gb] = NonLinearBending(ndof,nnel,dRdx,dRdy,u_e,v_e,w_e,phix_e,phiy_e,thetax_e,thetay_e);
% if strcmp(Parameter.Type,'Linear')
% Gb(:,:) = 0;
% end
% if strcmp(Parameter.Type,'vonKarman')
% Gb([1:4,7:end],:) = 0;
% end
% BNLb = Ab*Gb;
% kLNLb_e = kLNLb_e + BLb'*B2Matrix*BNLb*wtx*wty*detJ22;
% kNLLb_e = kNLLb_e + BNLb'*B3Matrix*BLb*wtx*wty*detJ22;
% kNLNLb_e = kNLNLb_e + BNLb'*B4Matrix*BNLb*wtx*wty*detJ22;
% Sb = InitialStressb_e(Parameter,ndof,nnel,R,dRdx,dRdy,...
% B1Matrix,B2Matrix,B3Matrix,B4Matrix,u_e,v_e,w_e,phix_e,phiy_e,thetax_e,thetay_e);
% tSigmab_e = tSigmab_e + Gb'*Sb*Gb*wtx*wty*detJ22;
%--------------------------------------------------------------------------
% compute force element vector
%--------------------------------------------------------------------------
if strcmp(Parameter.Load,'SSL')
fm_e = fm_e + Load*wtx*wty*detJ22*Force(ndof,nnel,R)*cos(pi/Lx*(x_e*R))*cos(pi/Ly*(y_e*R));
elseif strcmp(Parameter.Load,'UDL')
fm_e = fm_e + Load*wtx*wty*detJ22*Force(ndof,nnel,R);
else
error('Spartial Load is wrong\n');
end
end
end % end of numerical integration loop for bending term
%--------------------------------------------------------------------------
% numerical integration for shear term
%--------------------------------------------------------------------------
for intx=1:ngls
xi=points(intx,1); % sampling point in x-axis
wtx=weights(intx,1); % weight in x-axis
for inty=1:ngls
eta=points(inty,2); % sampling point in y-axis
wty=weights(inty,2) ; % weight in y-axis
% Compute Shape Function
[R,dRdxi,dRdeta] = ShapeFunctionQ9(xi,eta);
% compute Jacobian
[detJ22,invJ22] = Jacobian(dRdxi,dRdeta,x_e,y_e);
% derivatives w.r.t. physical coordinate
[dRdx,dRdy] = ShapefunctionDerivatives(dRdxi,dRdeta,invJ22);
%%
%--------------------------------------------------------------------------
% compute non-linear bending element matrix
%--------------------------------------------------------------------------
% linear bending kinematic matrix
BLb = LinearBending(ndof,nnel,B1Matrix,dRdx,dRdy);
% non-linear bending kinematic matrix
[Ab,Gb] = NonLinearBending(ndof,nnel,dRdx,dRdy,u_e,v_e,w_e,phix_e,phiy_e,thetax_e,thetay_e);
if strcmp(Parameter.Type,'Linear')
Gb(:,:) = 0;
end
if strcmp(Parameter.Type,'vonKarman')
Gb([1:4,7:end],:) = 0;
end
BNLb = Ab*Gb;
kLNLb_e = kLNLb_e + BLb'*B2Matrix*BNLb*wtx*wty*detJ22;
kNLLb_e = kNLLb_e + BNLb'*B3Matrix*BLb*wtx*wty*detJ22;
kNLNLb_e = kNLNLb_e + BNLb'*B4Matrix*BNLb*wtx*wty*detJ22;
%--------------------------------------------------------------------------
% compute initial stress tanget element matrix
%--------------------------------------------------------------------------
% non-linear initial inplane stress matrix
Sb = InitialStressb_e(Parameter,ndof,nnel,R,dRdx,dRdy,...
B1Matrix,B2Matrix,B3Matrix,B4Matrix,u_e,v_e,w_e,phix_e,phiy_e,thetax_e,thetay_e);
tSigmab_e = tSigmab_e + Gb'*Sb*Gb*wtx*wty*detJ22;
%%
%--------------------------------------------------------------------------
% compute linear and non-linear shear element matrix
%--------------------------------------------------------------------------
% Linear shear kinematic matrix
BLs = LinearShear(ndof,nnel,S1Matrix,R,dRdx,dRdy); % shear kinematic matrix
% non-linear shear kinematic matrix
[As,Gs] = NonLinearShear(ndof,nnel,R,dRdx,dRdy,u_e,v_e,w_e,phix_e,phiy_e,thetax_e,thetay_e);
if strcmp(Parameter.Type,'Linear')
Gs(:,:) = 0;
end
if strcmp(Parameter.Type,'vonKarman')
Gs(:,:) = 0;
end
BNLs = As*Gs;
kLLs_e = kLLs_e + BLs'*S1Matrix*BLs*wtx*wty*detJ22;
kLNLs_e = kLNLs_e + BLs'*S2Matrix*BNLs*wtx*wty*detJ22;
kNLLs_e = kNLLs_e + BNLs'*S3Matrix*BLs*wtx*wty*detJ22;
kNLNLs_e = kNLNLs_e + BNLs'*S4Matrix*BNLs*wtx*wty*detJ22;
%--------------------------------------------------------------------------
% compute initial stress tanget element matrix
%--------------------------------------------------------------------------
% non-linear initial shear stress matrix
Ss = InitialStresss_e(ndof,nnel,R,dRdx,dRdy,...
S1Matrix,S2Matrix,S3Matrix,S4Matrix,u_e,v_e,w_e,phix_e,phiy_e,thetax_e,thetay_e);
tSigmas_e = tSigmas_e + Gs'*Ss*Gs*wtx*wty*detJ22;
%%
%--------------------------------------------------------------------------
% compute penalty stiffness matrix
%--------------------------------------------------------------------------
gamma = Parameter.gamma;
[Px,Py] = PenaltyStiffness(ndof,nnel,R,dRdx,dRdy);
kp_e = kp_e + gamma*(Px'*Px+Py'*Py)*wtx*wty*detJ22*h;
end
end % end of numerical integration loop for shear term
% keyboard
%--------------------------------------------------------------------------
% compute element matrix and vector
%--------------------------------------------------------------------------
if strcmp(Parameter.Type,'Linear')
NL = 0;
else
NL = 1;
end
te = kp_e +...
kLLb_e + (kLNLb_e + kNLLb_e + kNLNLb_e + tSigmab_e)*NL +...
kLLs_e + (kLNLs_e + kNLLs_e + kNLNLs_e + tSigmas_e)*NL;
ke = kp_e +...
kLLb_e + (0.5*kLNLb_e + kNLLb_e + 0.5*kNLNLb_e)*NL +...
kLLs_e + (0.5*kLNLs_e + kNLLs_e + 0.5*kNLNLs_e)*NL;
fe = fm_e;
%--------------------------------------------------------------------------
% Assemble
%--------------------------------------------------------------------------
[tangent,stiffness,force] = Assemble(...
tangent,te,stiffness,ke,force,fe,IEN_e,ndof);
end
end