CalculiX  2.13
A Free Software Three-Dimensional Structural Finite Element Program
orthonl.f File Reference

Go to the source code of this file.

Functions/Subroutines

subroutine orthonl (w, vo, elas, s, ii1, jj1, weight)
 

Function/Subroutine Documentation

◆ orthonl()

subroutine orthonl ( real*8, dimension(3,3), intent(in)  w,
real*8, dimension(3,3), intent(in)  vo,
real*8, dimension(21), intent(in)  elas,
real*8, dimension(60,60), intent(inout)  s,
integer, intent(in)  ii1,
integer, intent(in)  jj1,
real*8, intent(in)  weight 
)
20 !
21 ! This routine replaces the following lines in e_c3d.f for
22 ! an orthotropic material
23 !
24 ! do i1=1,3
25 ! iii1=ii1+i1-1
26 ! do j1=1,3
27 ! jjj1=jj1+j1-1
28 ! do k1=1,3
29 ! do l1=1,3
30 ! s(iii1,jjj1)=s(iii1,jjj1)
31 ! & +anisox(i1,k1,j1,l1)*w(k1,l1)
32 ! do m1=1,3
33 ! s(iii1,jjj1)=s(iii1,jjj1)
34 ! & +anisox(i1,k1,m1,l1)*w(k1,l1)
35 ! & *vo(j1,m1)
36 ! & +anisox(m1,k1,j1,l1)*w(k1,l1)
37 ! & *vo(i1,m1)
38 ! do n1=1,3
39 ! s(iii1,jjj1)=s(iii1,jjj1)
40 ! & +anisox(m1,k1,n1,l1)
41 ! & *w(k1,l1)*vo(i1,m1)*vo(j1,n1)
42 ! enddo
43 ! enddo
44 ! enddo
45 ! enddo
46 ! enddo
47 ! enddo
48 !
49  implicit none
50 !
51  integer ii1,jj1
52 !
53  real*8 w(3,3),vo(3,3),elas(21),s(60,60),weight
54 !
55  intent(in) w,vo,elas,ii1,jj1,weight
56 !
57  intent(inout) s
58 !
59  s(ii1,jj1)=s(ii1,jj1)+((elas( 1)+elas( 1)*vo(1,1)
60  &+(elas( 1)+elas( 1)*vo(1,1)
61  &)*vo(1,1)+(elas( 7)*vo(1,2))*vo(1,2)
62  &+(elas( 8)*vo(1,3))
63  &*vo(1,3))*w(1,1)
64  &+(elas( 2)*vo(1,2)+(elas( 2)*vo(1,2))*vo(1,1)+(elas( 7)
65  &+elas( 7)*vo(1,1))*vo(1,2)
66  &)*w(1,2)
67  &+(elas( 4)*vo(1,3)+(elas( 4)*vo(1,3))*vo(1,1)
68  &+(elas( 8)+elas( 8)*vo(1,1))
69  &*vo(1,3))*w(1,3)
70  &+(elas( 7)*vo(1,2)+(elas( 7)*vo(1,2))*vo(1,1)+(elas( 2)
71  &+elas( 2)*vo(1,1))*vo(1,2)
72  &)*w(2,1)
73  &+(elas( 7)+elas( 7)*vo(1,1)
74  &+(elas( 7)+elas( 7)*vo(1,1)
75  &)*vo(1,1)+(elas( 3)*vo(1,2))*vo(1,2)
76  &+(elas( 9)*vo(1,3))
77  &*vo(1,3))*w(2,2)
78  &+((elas( 5)*vo(1,3))*vo(1,2)
79  &+(elas( 9)*vo(1,2))
80  &*vo(1,3))*w(2,3)
81  &+(elas( 8)*vo(1,3)+(elas( 8)*vo(1,3))*vo(1,1)
82  &+(elas( 4)+elas( 4)*vo(1,1))
83  &*vo(1,3))*w(3,1)
84  &+((elas( 9)*vo(1,3))*vo(1,2)
85  &+(elas( 5)*vo(1,2))
86  &*vo(1,3))*w(3,2)
87  &+(elas( 8)+elas( 8)*vo(1,1)
88  &+(elas( 8)+elas( 8)*vo(1,1)
89  &)*vo(1,1)+(elas( 9)*vo(1,2))*vo(1,2)
90  &+(elas( 6)*vo(1,3))
91  &*vo(1,3))*w(3,3))*weight
92  s(ii1,jj1+1)=s(ii1,jj1+1)+((elas( 1)*vo(2,1)
93  &+(elas( 1)*vo(2,1)
94  &)*vo(1,1)+(elas( 7)
95  &+elas( 7)*vo(2,2))*vo(1,2)
96  &+(elas( 8)*vo(2,3))
97  &*vo(1,3))*w(1,1)
98  &+(elas( 2)
99  &+elas( 2)*vo(2,2)+(elas( 2)
100  &+elas( 2)*vo(2,2))*vo(1,1)+(elas( 7)*vo(2,1))*vo(1,2)
101  &)*w(1,2)
102  &+(elas( 4)*vo(2,3)+(elas( 4)*vo(2,3))*vo(1,1)
103  &+(elas( 8)*vo(2,1))
104  &*vo(1,3))*w(1,3)
105  &+(elas( 7)
106  &+elas( 7)*vo(2,2)+(elas( 7)
107  &+elas( 7)*vo(2,2))*vo(1,1)+(elas( 2)*vo(2,1))*vo(1,2)
108  &)*w(2,1)
109  &+(elas( 7)*vo(2,1)
110  &+(elas( 7)*vo(2,1)
111  &)*vo(1,1)+(elas( 3)
112  &+elas( 3)*vo(2,2))*vo(1,2)
113  &+(elas( 9)*vo(2,3))
114  &*vo(1,3))*w(2,2)
115  &+((elas( 5)*vo(2,3))*vo(1,2)
116  &+(elas( 9)+elas( 9)*vo(2,2))
117  &*vo(1,3))*w(2,3)
118  &+(elas( 8)*vo(2,3)+(elas( 8)*vo(2,3))*vo(1,1)
119  &+(elas( 4)*vo(2,1))
120  &*vo(1,3))*w(3,1)
121  &+((elas( 9)*vo(2,3))*vo(1,2)
122  &+(elas( 5)+elas( 5)*vo(2,2))
123  &*vo(1,3))*w(3,2)
124  &+(elas( 8)*vo(2,1)
125  &+(elas( 8)*vo(2,1)
126  &)*vo(1,1)+(elas( 9)
127  &+elas( 9)*vo(2,2))*vo(1,2)
128  &+(elas( 6)*vo(2,3))
129  &*vo(1,3))*w(3,3))*weight
130  s(ii1,jj1+2)=s(ii1,jj1+2)+((elas( 1)*vo(3,1)
131  &+(elas( 1)*vo(3,1)
132  &)*vo(1,1)+(elas( 7)*vo(3,2))*vo(1,2)
133  &+(elas( 8)+elas( 8)*vo(3,3))
134  &*vo(1,3))*w(1,1)
135  &+(elas( 2)*vo(3,2)
136  &+(elas( 2)*vo(3,2))*vo(1,1)+(elas( 7)*vo(3,1))*vo(1,2)
137  &)*w(1,2)
138  &+(elas( 4)
139  &+elas( 4)*vo(3,3)+(elas( 4)
140  &+elas( 4)*vo(3,3))*vo(1,1)
141  &+(elas( 8)*vo(3,1))
142  &*vo(1,3))*w(1,3)
143  &+(elas( 7)*vo(3,2)+(elas( 7)*vo(3,2))*vo(1,1)
144  &+(elas( 2)*vo(3,1))*vo(1,2)
145  &)*w(2,1)
146  &+(elas( 7)*vo(3,1)
147  &+(elas( 7)*vo(3,1)
148  &)*vo(1,1)+(elas( 3)*vo(3,2))*vo(1,2)
149  &+(elas( 9)+elas( 9)*vo(3,3))
150  &*vo(1,3))*w(2,2)
151  &+((elas( 5)
152  &+elas( 5)*vo(3,3))*vo(1,2)
153  &+(elas( 9)*vo(3,2))
154  &*vo(1,3))*w(2,3)
155  &+(elas( 8)
156  &+elas( 8)*vo(3,3)+(elas( 8)
157  &+elas( 8)*vo(3,3))*vo(1,1)
158  &+(elas( 4)*vo(3,1))
159  &*vo(1,3))*w(3,1)
160  &+((elas( 9)
161  &+elas( 9)*vo(3,3))*vo(1,2)
162  &+(elas( 5)*vo(3,2))
163  &*vo(1,3))*w(3,2)
164  &+(elas( 8)*vo(3,1)
165  &+(elas( 8)*vo(3,1)
166  &)*vo(1,1)+(elas( 9)*vo(3,2))*vo(1,2)
167  &+(elas( 6)+elas( 6)*vo(3,3))
168  &*vo(1,3))*w(3,3))*weight
169  s(ii1+1,jj1)=s(ii1+1,jj1)+((elas( 7)*vo(1,2)
170  &+(elas( 1)+elas( 1)*vo(1,1)
171  &)*vo(2,1)+(elas( 7)*vo(1,2))*vo(2,2)
172  &+(elas( 8)*vo(1,3))
173  &*vo(2,3))*w(1,1)
174  &+(elas( 7)+elas( 7)*vo(1,1)
175  &+(elas( 2)*vo(1,2))*vo(2,1)+(elas( 7)
176  &+elas( 7)*vo(1,1))*vo(2,2)
177  &)*w(1,2)
178  &+((elas( 4)*vo(1,3))*vo(2,1)
179  &+(elas( 8)+elas( 8)*vo(1,1))
180  &*vo(2,3))*w(1,3)
181  &+(elas( 2)+elas( 2)*vo(1,1)
182  &+(elas( 7)*vo(1,2))*vo(2,1)+(elas( 2)
183  &+elas( 2)*vo(1,1))*vo(2,2)
184  &)*w(2,1)
185  &+(elas( 3)*vo(1,2)+(elas( 7)+elas( 7)*vo(1,1)
186  &)*vo(2,1)+(elas( 3)*vo(1,2))*vo(2,2)
187  &+(elas( 9)*vo(1,3))
188  &*vo(2,3))*w(2,2)
189  &+(elas( 5)*vo(1,3)+(elas( 5)*vo(1,3))*vo(2,2)
190  &+(elas( 9)*vo(1,2))
191  &*vo(2,3))*w(2,3)
192  &+((elas( 8)*vo(1,3))*vo(2,1)
193  &+(elas( 4)+elas( 4)*vo(1,1))
194  &*vo(2,3))*w(3,1)
195  &+(elas( 9)*vo(1,3)+(elas( 9)*vo(1,3))*vo(2,2)
196  &+(elas( 5)*vo(1,2))
197  &*vo(2,3))*w(3,2)
198  &+(elas( 9)*vo(1,2)+(elas( 8)+elas( 8)*vo(1,1)
199  &)*vo(2,1)+(elas( 9)*vo(1,2))*vo(2,2)
200  &+(elas( 6)*vo(1,3))
201  &*vo(2,3))*w(3,3))*weight
202  s(ii1+1,jj1+1)=s(ii1+1,jj1+1)+((elas( 7)
203  &+elas( 7)*vo(2,2)+(elas( 1)*vo(2,1)
204  &)*vo(2,1)+(elas( 7)
205  &+elas( 7)*vo(2,2))*vo(2,2)
206  &+(elas( 8)*vo(2,3))
207  &*vo(2,3))*w(1,1)
208  &+(elas( 7)*vo(2,1)
209  &+(elas( 2)
210  &+elas( 2)*vo(2,2))*vo(2,1)+(elas( 7)*vo(2,1))*vo(2,2)
211  &)*w(1,2)
212  &+((elas( 4)*vo(2,3))*vo(2,1)
213  &+(elas( 8)*vo(2,1))
214  &*vo(2,3))*w(1,3)
215  &+(elas( 2)*vo(2,1)
216  &+(elas( 7)
217  &+elas( 7)*vo(2,2))*vo(2,1)+(elas( 2)*vo(2,1))*vo(2,2)
218  &)*w(2,1)
219  &+(elas( 3)
220  &+elas( 3)*vo(2,2)+(elas( 7)*vo(2,1)
221  &)*vo(2,1)+(elas( 3)
222  &+elas( 3)*vo(2,2))*vo(2,2)
223  &+(elas( 9)*vo(2,3))
224  &*vo(2,3))*w(2,2)
225  &+(elas( 5)*vo(2,3)+(elas( 5)*vo(2,3))*vo(2,2)
226  &+(elas( 9)+elas( 9)*vo(2,2))
227  &*vo(2,3))*w(2,3)
228  &+((elas( 8)*vo(2,3))*vo(2,1)
229  &+(elas( 4)*vo(2,1))
230  &*vo(2,3))*w(3,1)
231  &+(elas( 9)*vo(2,3)+(elas( 9)*vo(2,3))*vo(2,2)
232  &+(elas( 5)+elas( 5)*vo(2,2))
233  &*vo(2,3))*w(3,2)
234  &+(elas( 9)
235  &+elas( 9)*vo(2,2)+(elas( 8)*vo(2,1)
236  &)*vo(2,1)+(elas( 9)
237  &+elas( 9)*vo(2,2))*vo(2,2)
238  &+(elas( 6)*vo(2,3))
239  &*vo(2,3))*w(3,3))*weight
240  s(ii1+1,jj1+2)=s(ii1+1,jj1+2)+((elas( 7)*vo(3,2)+(elas( 1)*vo(3,1)
241  &)*vo(2,1)+(elas( 7)*vo(3,2))*vo(2,2)
242  &+(elas( 8)+elas( 8)*vo(3,3))
243  &*vo(2,3))*w(1,1)
244  &+(elas( 7)*vo(3,1)
245  &+(elas( 2)*vo(3,2))*vo(2,1)+(elas( 7)*vo(3,1))*vo(2,2)
246  &)*w(1,2)
247  &+((elas( 4)
248  &+elas( 4)*vo(3,3))*vo(2,1)
249  &+(elas( 8)*vo(3,1))
250  &*vo(2,3))*w(1,3)
251  &+(elas( 2)*vo(3,1)
252  &+(elas( 7)*vo(3,2))*vo(2,1)+(elas( 2)*vo(3,1))*vo(2,2)
253  &)*w(2,1)
254  &+(elas( 3)*vo(3,2)+(elas( 7)*vo(3,1)
255  &)*vo(2,1)+(elas( 3)*vo(3,2))*vo(2,2)
256  &+(elas( 9)+elas( 9)*vo(3,3))
257  &*vo(2,3))*w(2,2)
258  &+(elas( 5)
259  &+elas( 5)*vo(3,3)+(elas( 5)
260  &+elas( 5)*vo(3,3))*vo(2,2)
261  &+(elas( 9)*vo(3,2))
262  &*vo(2,3))*w(2,3)
263  &+((elas( 8)
264  &+elas( 8)*vo(3,3))*vo(2,1)
265  &+(elas( 4)*vo(3,1))
266  &*vo(2,3))*w(3,1)
267  &+(elas( 9)
268  &+elas( 9)*vo(3,3)+(elas( 9)
269  &+elas( 9)*vo(3,3))*vo(2,2)
270  &+(elas( 5)*vo(3,2))
271  &*vo(2,3))*w(3,2)
272  &+(elas( 9)*vo(3,2)+(elas( 8)*vo(3,1)
273  &)*vo(2,1)+(elas( 9)*vo(3,2))*vo(2,2)
274  &+(elas( 6)+elas( 6)*vo(3,3))
275  &*vo(2,3))*w(3,3))*weight
276  s(ii1+2,jj1)=s(ii1+2,jj1)+((elas( 8)*vo(1,3)
277  &+(elas( 1)+elas( 1)*vo(1,1)
278  &)*vo(3,1)+(elas( 7)*vo(1,2))*vo(3,2)
279  &+(elas( 8)*vo(1,3))
280  &*vo(3,3))*w(1,1)
281  &+((elas( 2)*vo(1,2))*vo(3,1)+(elas( 7)
282  &+elas( 7)*vo(1,1))*vo(3,2)
283  &)*w(1,2)
284  &+(elas( 8)+elas( 8)*vo(1,1)
285  &+(elas( 4)*vo(1,3))*vo(3,1)
286  &+(elas( 8)+elas( 8)*vo(1,1))
287  &*vo(3,3))*w(1,3)
288  &+((elas( 7)*vo(1,2))*vo(3,1)+(elas( 2)
289  &+elas( 2)*vo(1,1))*vo(3,2)
290  &)*w(2,1)
291  &+(elas( 9)*vo(1,3)+(elas( 7)+elas( 7)*vo(1,1)
292  &)*vo(3,1)+(elas( 3)*vo(1,2))*vo(3,2)
293  &+(elas( 9)*vo(1,3))
294  &*vo(3,3))*w(2,2)
295  &+(elas( 9)*vo(1,2)+(elas( 5)*vo(1,3))*vo(3,2)
296  &+(elas( 9)*vo(1,2))
297  &*vo(3,3))*w(2,3)
298  &+(elas( 4)+elas( 4)*vo(1,1)
299  &+(elas( 8)*vo(1,3))*vo(3,1)
300  &+(elas( 4)+elas( 4)*vo(1,1))
301  &*vo(3,3))*w(3,1)
302  &+(elas( 5)*vo(1,2)+(elas( 9)*vo(1,3))*vo(3,2)
303  &+(elas( 5)*vo(1,2))
304  &*vo(3,3))*w(3,2)
305  &+(elas( 6)*vo(1,3)+(elas( 8)+elas( 8)*vo(1,1)
306  &)*vo(3,1)+(elas( 9)*vo(1,2))*vo(3,2)
307  &+(elas( 6)*vo(1,3))
308  &*vo(3,3))*w(3,3))*weight
309  s(ii1+2,jj1+1)=s(ii1+2,jj1+1)+((elas( 8)*vo(2,3)
310  &+(elas( 1)*vo(2,1)
311  &)*vo(3,1)+(elas( 7)
312  &+elas( 7)*vo(2,2))*vo(3,2)
313  &+(elas( 8)*vo(2,3))
314  &*vo(3,3))*w(1,1)
315  &+((elas( 2)
316  &+elas( 2)*vo(2,2))*vo(3,1)+(elas( 7)*vo(2,1))*vo(3,2)
317  &)*w(1,2)
318  &+(elas( 8)*vo(2,1)
319  &+(elas( 4)*vo(2,3))*vo(3,1)
320  &+(elas( 8)*vo(2,1))
321  &*vo(3,3))*w(1,3)
322  &+((elas( 7)
323  &+elas( 7)*vo(2,2))*vo(3,1)+(elas( 2)*vo(2,1))*vo(3,2)
324  &)*w(2,1)
325  &+(elas( 9)*vo(2,3)+(elas( 7)*vo(2,1)
326  &)*vo(3,1)+(elas( 3)
327  &+elas( 3)*vo(2,2))*vo(3,2)
328  &+(elas( 9)*vo(2,3))
329  &*vo(3,3))*w(2,2)
330  &+(elas( 9)
331  &+elas( 9)*vo(2,2)+(elas( 5)*vo(2,3))*vo(3,2)
332  &+(elas( 9)+elas( 9)*vo(2,2))
333  &*vo(3,3))*w(2,3)
334  &+(elas( 4)*vo(2,1)
335  &+(elas( 8)*vo(2,3))*vo(3,1)
336  &+(elas( 4)*vo(2,1))
337  &*vo(3,3))*w(3,1)
338  &+(elas( 5)
339  &+elas( 5)*vo(2,2)+(elas( 9)*vo(2,3))*vo(3,2)
340  &+(elas( 5)+elas( 5)*vo(2,2))
341  &*vo(3,3))*w(3,2)
342  &+(elas( 6)*vo(2,3)+(elas( 8)*vo(2,1)
343  &)*vo(3,1)+(elas( 9)
344  &+elas( 9)*vo(2,2))*vo(3,2)
345  &+(elas( 6)*vo(2,3))
346  &*vo(3,3))*w(3,3))*weight
347  s(ii1+2,jj1+2)=s(ii1+2,jj1+2)+((elas( 8)
348  &+elas( 8)*vo(3,3)+(elas( 1)*vo(3,1)
349  &)*vo(3,1)+(elas( 7)*vo(3,2))*vo(3,2)
350  &+(elas( 8)+elas( 8)*vo(3,3))
351  &*vo(3,3))*w(1,1)
352  &+((elas( 2)*vo(3,2))*vo(3,1)+(elas( 7)*vo(3,1))*vo(3,2)
353  &)*w(1,2)
354  &+(elas( 8)*vo(3,1)
355  &+(elas( 4)
356  &+elas( 4)*vo(3,3))*vo(3,1)
357  &+(elas( 8)*vo(3,1))
358  &*vo(3,3))*w(1,3)
359  &+((elas( 7)*vo(3,2))*vo(3,1)+(elas( 2)*vo(3,1))*vo(3,2)
360  &)*w(2,1)
361  &+(elas( 9)
362  &+elas( 9)*vo(3,3)+(elas( 7)*vo(3,1)
363  &)*vo(3,1)+(elas( 3)*vo(3,2))*vo(3,2)
364  &+(elas( 9)+elas( 9)*vo(3,3))
365  &*vo(3,3))*w(2,2)
366  &+(elas( 9)*vo(3,2)+(elas( 5)
367  &+elas( 5)*vo(3,3))*vo(3,2)
368  &+(elas( 9)*vo(3,2))
369  &*vo(3,3))*w(2,3)
370  &+(elas( 4)*vo(3,1)
371  &+(elas( 8)
372  &+elas( 8)*vo(3,3))*vo(3,1)
373  &+(elas( 4)*vo(3,1))
374  &*vo(3,3))*w(3,1)
375  &+(elas( 5)*vo(3,2)+(elas( 9)
376  &+elas( 9)*vo(3,3))*vo(3,2)
377  &+(elas( 5)*vo(3,2))
378  &*vo(3,3))*w(3,2)
379  &+(elas( 6)
380  &+elas( 6)*vo(3,3)+(elas( 8)*vo(3,1)
381  &)*vo(3,1)+(elas( 9)*vo(3,2))*vo(3,2)
382  &+(elas( 6)+elas( 6)*vo(3,3))
383  &*vo(3,3))*w(3,3))*weight
384 !
385  return
Hosted by OpenAircraft.com, (Michigan UAV, LLC)