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

Go to the source code of this file.

Functions/Subroutines

subroutine printoutint (prlab, ipkon, lakon, stx, eei, xstate, ener, mi, nstate_, ii, nelem, qfx, orab, ielorien, norien, co, konf, ielmat, thicke, eme, ielprop, prop, nelel, ithermal, orname)
 

Function/Subroutine Documentation

◆ printoutint()

subroutine printoutint ( character*6, dimension(*)  prlab,
integer, dimension(*)  ipkon,
character*8, dimension(*)  lakon,
real*8, dimension(6,mi(1),*)  stx,
real*8, dimension(6,mi(1),*)  eei,
real*8, dimension(nstate_,mi(1),*)  xstate,
real*8, dimension(mi(1),*)  ener,
integer, dimension(*)  mi,
integer  nstate_,
integer  ii,
integer  nelem,
real*8, dimension(3,mi(1),*)  qfx,
real*8, dimension(7,*)  orab,
integer, dimension(mi(3),*)  ielorien,
integer  norien,
real*8, dimension(3,*)  co,
integer, dimension(*)  konf,
integer, dimension(mi(3),*)  ielmat,
real*8, dimension(mi(3),*)  thicke,
real*8, dimension(6,mi(1),*)  eme,
integer, dimension(*)  ielprop,
real*8, dimension(*)  prop,
integer  nelel,
integer, dimension(2)  ithermal,
character*80, dimension(*)  orname 
)
22 !
23 ! stores integration point results for element "nelem" in the .dat file
24 !
25 ! nelem is the element number from the input deck
26 ! nelel is a local number, created e.g. by renumbering; thus
27 ! far only applicable for CFD
28 !
29  implicit none
30 !
31  character*6 prlab(*)
32  character*8 lakon(*)
33  character*80 orname(*)
34 !
35  integer ipkon(*),mi(*),nstate_,nelem,l,ii,mint3d,j,k,nope,
36  & ielorien(mi(3),*),norien,konf(*),konl,indexe,m,iorien,iflag,
37  & ielmat(mi(3),*),nopes,mint2d,kk,ki,kl,nlayer,ilayer,
38  & null,ielprop(*),nelel,ithermal(2)
39 !
40  real*8 stx(6,mi(1),*),eei(6,mi(1),*),xstate(nstate_,mi(1),*),
41  & ener(mi(1),*),qfx(3,mi(1),*),xi,et,ze,xl(3,20),xsj,shp(4,20),
42  & coords(3,27),weight,orab(7,*),co(3,*),a(3,3),b(3,3),c(3,3),
43  & qfxl(3),thicke(mi(3),*),xsj2(3),shp2(7,8),xl2(3,8),xs2(3,7),
44  & thickness,tlayer(4),dlayer(4),xlayer(mi(3),4),eme(6,mi(1),*),
45  & prop(*)
46 !
47  include "gauss.f"
48 !
49  data iflag /1/
50  data null /0/
51 !
52  if(ipkon(nelel).lt.0) then
53  return
54  else
55  indexe=ipkon(nelel)
56  endif
57 !
58 ! check whether transformation is necessary (if orientation
59 ! is applied and output in local system is requested)
60 !
61  if(lakon(nelel)(7:8).ne.'LC') then
62  if((norien.eq.0).or.(prlab(ii)(6:6).eq.'G')) then
63  iorien=0
64  else
65  iorien=ielorien(1,nelel)
66  endif
67  elseif(lakon(nelel)(4:5).eq.'20') then
68 !
69 ! composite materials
70 !
71  if((norien.eq.0).or.(prlab(ii)(6:6).eq.'G')) then
72  iorien=0
73  else
74 !
75 ! check whether at least one layer has a transformation
76 !
77  iorien=0
78  do k=1,mi(3)
79  if(ielorien(k,nelel).ne.0) then
80  iorien=ielorien(k,nelel)
81  exit
82  endif
83  enddo
84  endif
85 !
86 ! determining the number of layers
87 !
88  mint2d=4
89  nopes=8
90 !
91  nlayer=0
92  do k=1,mi(3)
93  if(ielmat(k,nelel).ne.0) then
94  nlayer=nlayer+1
95  endif
96  enddo
97 !
98 ! determining the layer thickness and global thickness
99 ! at the shell integration points
100 !
101  iflag=1
102  do kk=1,mint2d
103  xi=gauss3d2(1,kk)
104  et=gauss3d2(2,kk)
105  call shape8q(xi,et,xl2,xsj2,xs2,shp2,iflag)
106  tlayer(kk)=0.d0
107  do k=1,nlayer
108  thickness=0.d0
109  do j=1,nopes
110  thickness=thickness+thicke(k,indexe+j)*shp2(4,j)
111  enddo
112  tlayer(kk)=tlayer(kk)+thickness
113  xlayer(k,kk)=thickness
114  enddo
115  enddo
116  iflag=3
117 !
118  ilayer=0
119  do k=1,4
120  dlayer(k)=0.d0
121  enddo
122 !
123 !
124  elseif(lakon(nelel)(4:5).eq.'15') then
125 !
126 ! composite materials
127 !
128  if((norien.eq.0).or.(prlab(ii)(6:6).eq.'G')) then
129  iorien=0
130  else
131 !
132 ! check whether at least one layer has a transformation
133 !
134  iorien=0
135  do k=1,mi(3)
136  if(ielorien(k,nelel).ne.0) then
137  iorien=ielorien(k,nelel)
138  exit
139  endif
140  enddo
141  endif
142 !
143 ! determining the number of layers
144 !
145  mint2d=3
146  nopes=6
147 !
148  nlayer=0
149  do k=1,mi(3)
150  if(ielmat(k,nelel).ne.0) then
151  nlayer=nlayer+1
152  endif
153  enddo
154 !
155 ! determining the layer thickness and global thickness
156 ! at the shell integration points
157 !
158  iflag=1
159  do kk=1,mint2d
160  xi=gauss3d10(1,kk)
161  et=gauss3d10(2,kk)
162  call shape6tri(xi,et,xl2,xsj2,xs2,shp2,iflag)
163  tlayer(kk)=0.d0
164  do k=1,nlayer
165  thickness=0.d0
166  do j=1,nopes
167  thickness=thickness+thicke(k,indexe+j)*shp2(4,j)
168  enddo
169  tlayer(kk)=tlayer(kk)+thickness
170  xlayer(k,kk)=thickness
171  enddo
172  enddo
173  iflag=3
174 !
175  ilayer=0
176  do k=1,3
177  dlayer(k)=0.d0
178  enddo
179 !
180  endif
181 !
182 ! number of integration points
183 !
184  if((lakon(nelel)(4:5).eq.'8R').or.
185  & (lakon(nelel)(1:1).eq.'F')) then
186  mint3d=1
187  elseif(lakon(nelel)(4:7).eq.'20RB') then
188  if((lakon(nelel)(8:8).eq.'R').or.
189  & (lakon(nelel)(8:8).eq.'C')) then
190  mint3d=50
191  else
192  call beamintscheme(lakon(nelel),mint3d,ielprop(nelel),prop,
193  & null,xi,et,ze,weight)
194  endif
195  elseif((lakon(nelel)(4:4).eq.'8').or.
196  & (lakon(nelel)(4:6).eq.'20R')) then
197  if(lakon(nelel)(7:8).eq.'LC') then
198  mint3d=8*nlayer
199  else
200  mint3d=8
201  endif
202  elseif(lakon(nelel)(4:4).eq.'2') then
203  mint3d=27
204  elseif(lakon(nelel)(4:5).eq.'10') then
205  mint3d=4
206  elseif(lakon(nelel)(4:4).eq.'4') then
207  mint3d=1
208  elseif(lakon(nelel)(4:5).eq.'15') then
209  if(lakon(nelel)(7:8).eq.'LC') then
210  mint3d=6*nlayer
211  else
212  mint3d=9
213  endif
214  elseif(lakon(nelel)(4:4).eq.'6') then
215  mint3d=2
216  elseif(lakon(nelel)(1:1).eq.'U') then
217 c if(lakon(nelel)(4:4).eq.' ') then
218 c mint3d=ichar(lakon(nelel)(5:5))-48
219 c else
220 c mint3d=10*(ichar(lakon(nelel)(4:4))-48)
221 c & +ichar(lakon(nelel)(5:5))-48
222 c endif
223  mint3d=ichar(lakon(nelel)(6:6))
224  else
225  return
226  endif
227 !
228 ! calculation of the integration point coordinates for
229 ! output in the local system (if needed)
230 !
231  if(iorien.ne.0) then
232  if(lakon(nelel)(4:4).eq.'2') then
233  nope=20
234  elseif(lakon(nelel)(4:4).eq.'8') then
235  nope=8
236  elseif(lakon(nelel)(4:5).eq.'10') then
237  nope=10
238  elseif(lakon(nelel)(4:4).eq.'4') then
239  nope=4
240  elseif(lakon(nelel)(4:5).eq.'15') then
241  nope=15
242  elseif(lakon(nelel)(4:4).eq.'6') then
243  nope=6
244  endif
245 !
246 c indexe=ipkon(nelel)
247  do j=1,nope
248  konl=konf(indexe+j)
249  do k=1,3
250  xl(k,j)=co(k,konl)
251  enddo
252  enddo
253 !
254  do j=1,mint3d
255  if((lakon(nelel)(4:5).eq.'8R').or.
256  & (lakon(nelel)(1:4).eq.'F3D8')) then
257  xi=gauss3d1(1,j)
258  et=gauss3d1(2,j)
259  ze=gauss3d1(3,j)
260  weight=weight3d1(j)
261  elseif(lakon(nelel)(4:8).eq.'20RB') then
262  if((lakon(nelel)(8:8).eq.'B').or.
263  & (lakon(nelel)(8:8).eq.'C')) then
264  xi=gauss3d13(1,j)
265  et=gauss3d13(2,j)
266  ze=gauss3d13(3,j)
267  weight=weight3d13(j)
268  else
269  call beamintscheme(lakon(nelel),mint3d,ielprop(nelel),
270  & prop,j,xi,et,ze,weight)
271  endif
272  elseif((lakon(nelel)(4:4).eq.'8').or.
273  & (lakon(nelel)(4:6).eq.'20R'))
274  & then
275  if(lakon(nelel)(7:8).ne.'LC') then
276  xi=gauss3d2(1,j)
277  et=gauss3d2(2,j)
278  ze=gauss3d2(3,j)
279  weight=weight3d2(j)
280  else
281  kl=mod(j,8)
282  if(kl.eq.0) kl=8
283 !
284  xi=gauss3d2(1,kl)
285  et=gauss3d2(2,kl)
286  ze=gauss3d2(3,kl)
287  weight=weight3d2(kl)
288 !
289  ki=mod(j,4)
290  if(ki.eq.0) ki=4
291 !
292  if(kl.eq.1) then
293  ilayer=ilayer+1
294  if(ilayer.gt.1) then
295  do k=1,4
296  dlayer(k)=dlayer(k)+xlayer(ilayer-1,k)
297  enddo
298  endif
299  endif
300  ze=2.d0*(dlayer(ki)+(ze+1.d0)/2.d0*xlayer(ilayer,ki))/
301  & tlayer(ki)-1.d0
302  weight=weight*xlayer(ilayer,ki)/tlayer(ki)
303  endif
304  elseif(lakon(nelel)(4:4).eq.'2') then
305  xi=gauss3d3(1,j)
306  et=gauss3d3(2,j)
307  ze=gauss3d3(3,j)
308  weight=weight3d3(j)
309  elseif(lakon(nelel)(4:5).eq.'10') then
310  xi=gauss3d5(1,j)
311  et=gauss3d5(2,j)
312  ze=gauss3d5(3,j)
313  weight=weight3d5(j)
314  elseif(lakon(nelel)(4:4).eq.'4') then
315  xi=gauss3d4(1,j)
316  et=gauss3d4(2,j)
317  ze=gauss3d4(3,j)
318  weight=weight3d4(j)
319  elseif(lakon(nelel)(4:5).eq.'15') then
320  if(lakon(nelel)(7:8).ne.'LC') then
321  xi=gauss3d8(1,j)
322  et=gauss3d8(2,j)
323  ze=gauss3d8(3,j)
324  weight=weight3d8(j)
325  else
326  kl=mod(j,6)
327  if(kl.eq.0) kl=6
328 !
329  xi=gauss3d10(1,kl)
330  et=gauss3d10(2,kl)
331  ze=gauss3d10(3,kl)
332  weight=weight3d10(kl)
333 !
334  ki=mod(j,3)
335  if(ki.eq.0) ki=3
336 !
337  if(kl.eq.1) then
338  ilayer=ilayer+1
339  if(ilayer.gt.1) then
340  do k=1,3
341  dlayer(k)=dlayer(k)+xlayer(ilayer-1,k)
342  enddo
343  endif
344  endif
345  ze=2.d0*(dlayer(ki)+(ze+1.d0)/2.d0*xlayer(ilayer,ki))/
346  & tlayer(ki)-1.d0
347  weight=weight*xlayer(ilayer,ki)/tlayer(ki)
348  endif
349  elseif(lakon(nelel)(1:4).eq.'C3D6') then
350  xi=gauss3d7(1,j)
351  et=gauss3d7(2,j)
352  ze=gauss3d7(3,j)
353  weight=weight3d7(j)
354  elseif(lakon(nelel)(1:4).eq.'F3D6') then
355  xi=gauss3d14(1,j)
356  et=gauss3d14(2,j)
357  ze=gauss3d14(3,j)
358  weight=weight3d14(j)
359  endif
360 !
361  if(nope.eq.20) then
362  call shape20h(xi,et,ze,xl,xsj,shp,iflag)
363  elseif(nope.eq.8) then
364  call shape8h(xi,et,ze,xl,xsj,shp,iflag)
365  elseif(nope.eq.10) then
366  call shape10tet(xi,et,ze,xl,xsj,shp,iflag)
367  elseif(nope.eq.4) then
368  call shape4tet(xi,et,ze,xl,xsj,shp,iflag)
369  elseif(nope.eq.15) then
370  call shape15w(xi,et,ze,xl,xsj,shp,iflag)
371  else
372  call shape6w(xi,et,ze,xl,xsj,shp,iflag)
373  endif
374 !
375  do k=1,3
376  coords(k,j)=0.d0
377  do l=1,nope
378  coords(k,j)=coords(k,j)+xl(k,l)*shp(4,l)
379  enddo
380  enddo
381  enddo
382  endif
383 !
384  if((prlab(ii)(1:4).eq.'S ').or.(prlab(ii)(1:4).eq.'SVF ')) then
385  do j=1,mint3d
386 !
387 ! composite materials
388 !
389  if(lakon(nelel)(7:8).eq.'LC') then
390  if((norien.eq.0).or.(prlab(ii)(6:6).eq.'G')) then
391  iorien=0
392  elseif(lakon(nelel)(4:5).eq.'20') then
393  kl=mod(j,8)
394  if(kl.eq.0) kl=8
395  ilayer=(j-kl)/8+1
396  iorien=ielorien(ilayer,nelel)
397  elseif(lakon(nelel)(4:5).eq.'15') then
398  kl=mod(j,6)
399  if(kl.eq.0) kl=6
400  ilayer=(j-kl)/6+1
401  iorien=ielorien(ilayer,nelel)
402  endif
403  endif
404 !
405  if(iorien.eq.0) then
406  write(5,'(i10,1x,i3,1p,6(1x,e13.6))') nelem,j,
407  & (stx(k,j,nelel),k=1,6)
408  else
409  call transformatrix(orab(1,iorien),coords(1,j),a)
410  b(1,1)=stx(1,j,nelel)
411  b(2,2)=stx(2,j,nelel)
412  b(3,3)=stx(3,j,nelel)
413  b(1,2)=stx(4,j,nelel)
414  b(1,3)=stx(5,j,nelel)
415  b(2,3)=stx(6,j,nelel)
416  b(2,1)=b(1,2)
417  b(3,1)=b(1,3)
418  b(3,2)=b(2,3)
419  do k=1,3
420  do l=1,3
421  c(k,l)=0.d0
422  do m=1,3
423  c(k,l)=c(k,l)+b(k,m)*a(m,l)
424  enddo
425  enddo
426  enddo
427  do k=1,3
428  do l=k,3
429  b(k,l)=0.d0
430  do m=1,3
431  b(k,l)=b(k,l)+a(m,k)*c(m,l)
432  enddo
433  enddo
434  enddo
435  write(5,'(i10,1x,i3,1p,6(1x,e13.6),1x,a20)') nelem,j,
436  & b(1,1),b(2,2),b(3,3),b(1,2),b(1,3),b(2,3),
437  & orname(iorien)(1:20)
438  endif
439  enddo
440  elseif(prlab(ii)(1:4).eq.'E ') then
441  do j=1,mint3d
442 !
443 ! composite materials
444 !
445  if(lakon(nelel)(7:8).eq.'LC') then
446  if((norien.eq.0).or.(prlab(ii)(6:6).eq.'G')) then
447  iorien=0
448  elseif(lakon(nelel)(4:5).eq.'20') then
449  kl=mod(j,8)
450  if(kl.eq.0) kl=8
451  ilayer=(j-kl)/8+1
452  iorien=ielorien(ilayer,nelel)
453  elseif(lakon(nelel)(4:5).eq.'15') then
454  kl=mod(j,6)
455  if(kl.eq.0) kl=6
456  ilayer=(j-kl)/6+1
457  iorien=ielorien(ilayer,nelel)
458  endif
459  endif
460 !
461  if(iorien.eq.0) then
462  write(5,'(i10,1x,i3,1p,6(1x,e13.6))') nelem,j,
463  & (eei(k,j,nelel),k=1,6)
464  else
465  call transformatrix(orab(1,iorien),coords(1,j),a)
466  b(1,1)=eei(1,j,nelel)
467  b(2,2)=eei(2,j,nelel)
468  b(3,3)=eei(3,j,nelel)
469  b(1,2)=eei(4,j,nelel)
470  b(1,3)=eei(5,j,nelel)
471  b(2,3)=eei(6,j,nelel)
472  b(2,1)=b(1,2)
473  b(3,1)=b(1,3)
474  b(3,2)=b(2,3)
475  do k=1,3
476  do l=1,3
477  c(k,l)=0.d0
478  do m=1,3
479  c(k,l)=b(k,m)*a(m,l)
480  enddo
481  enddo
482  enddo
483  do k=1,3
484  do l=k,3
485  b(k,l)=0.d0
486  do m=1,3
487  b(k,l)=a(m,k)*c(m,l)
488  enddo
489  enddo
490  enddo
491  write(5,'(i10,1x,i3,1p,6(1x,e13.6),1x,a20)') nelem,j,
492  & b(1,1),b(2,2),b(3,3),b(1,2),b(1,3),b(2,3),
493  & orname(iorien)(1:20)
494  endif
495  enddo
496  elseif(prlab(ii)(1:4).eq.'ME ') then
497  do j=1,mint3d
498 !
499 ! composite materials
500 !
501  if(lakon(nelel)(7:8).eq.'LC') then
502  if((norien.eq.0).or.(prlab(ii)(6:6).eq.'G')) then
503  iorien=0
504  elseif(lakon(nelel)(4:5).eq.'20') then
505  kl=mod(j,8)
506  if(kl.eq.0) kl=8
507  ilayer=(j-kl)/8+1
508  iorien=ielorien(ilayer,nelel)
509  elseif(lakon(nelel)(4:5).eq.'15') then
510  kl=mod(j,6)
511  if(kl.eq.0) kl=6
512  ilayer=(j-kl)/6+1
513  iorien=ielorien(ilayer,nelel)
514  endif
515  endif
516 !
517  if(iorien.eq.0) then
518  write(5,'(i10,1x,i3,1p,6(1x,e13.6))') nelem,j,
519  & (eme(k,j,nelel),k=1,6)
520  else
521  call transformatrix(orab(1,iorien),coords(1,j),a)
522  b(1,1)=eme(1,j,nelel)
523  b(2,2)=eme(2,j,nelel)
524  b(3,3)=eme(3,j,nelel)
525  b(1,2)=eme(4,j,nelel)
526  b(1,3)=eme(5,j,nelel)
527  b(2,3)=eme(6,j,nelel)
528  b(2,1)=b(1,2)
529  b(3,1)=b(1,3)
530  b(3,2)=b(2,3)
531  do k=1,3
532  do l=1,3
533  do m=1,3
534  c(k,l)=b(k,m)*a(m,l)
535  enddo
536  enddo
537  enddo
538  do k=1,3
539  do l=k,3
540  do m=1,3
541  b(k,l)=a(m,k)*c(m,l)
542  enddo
543  enddo
544  enddo
545  write(5,'(i10,1x,i3,1p,6(1x,e13.6),1x,a20)') nelem,j,
546  & b(1,1),b(2,2),b(3,3),b(1,2),b(1,3),b(2,3),
547  & orname(iorien)(1:20)
548  endif
549  enddo
550  elseif(prlab(ii)(1:4).eq.'PEEQ') then
551  do j=1,mint3d
552  write(5,'(i10,1x,i3,1p,6(1x,e13.6))') nelem,j,
553  & xstate(1,j,nelel)
554  enddo
555  elseif(prlab(ii)(1:4).eq.'ENER') then
556  do j=1,mint3d
557  write(5,'(i10,1x,i3,1p,6(1x,e13.6))') nelem,j,
558  & ener(j,nelel)
559  enddo
560  elseif(prlab(ii)(1:4).eq.'SDV ') then
561  do j=1,mint3d
562 !
563 ! composite materials
564 !
565  if(lakon(nelel)(7:8).eq.'LC') then
566  if((norien.eq.0).or.(prlab(ii)(6:6).eq.'G')) then
567  iorien=0
568  elseif(lakon(nelel)(4:5).eq.'20') then
569  kl=mod(j,8)
570  if(kl.eq.0) kl=8
571  ilayer=(j-kl)/8+1
572  iorien=ielorien(ilayer,nelel)
573  elseif(lakon(nelel)(4:5).eq.'15') then
574  kl=mod(j,6)
575  if(kl.eq.0) kl=6
576  ilayer=(j-kl)/6+1
577  iorien=ielorien(ilayer,nelel)
578  endif
579  endif
580 !
581  if(iorien.ne.0) then
582  write(*,*) '*WARNING in printoutint: SDV cannot be'
583  write(*,*) ' printed in the local system'
584  write(*,*) ' results are in the global system'
585  endif
586  write(5,'(i10,1x,i3,1p,99(1x,e13.6))') nelem,j,
587  & (xstate(k,j,nelel),k=1,nstate_)
588  enddo
589  elseif(((prlab(ii)(1:4).eq.'HFL ').or.(prlab(ii)(1:4).eq.'HFLF'))
590  & .and.(ithermal(1).gt.1)) then
591  do j=1,mint3d
592 !
593 ! composite materials
594 !
595  if(lakon(nelel)(7:8).eq.'LC') then
596  if((norien.eq.0).or.(prlab(ii)(6:6).eq.'G')) then
597  iorien=0
598  elseif(lakon(nelel)(4:5).eq.'20') then
599  kl=mod(j,8)
600  if(kl.eq.0) kl=8
601  ilayer=(j-kl)/8+1
602  iorien=ielorien(ilayer,nelel)
603  elseif(lakon(nelel)(4:5).eq.'15') then
604  kl=mod(j,6)
605  if(kl.eq.0) kl=6
606  ilayer=(j-kl)/6+1
607  iorien=ielorien(ilayer,nelel)
608  endif
609  endif
610 !
611  if(iorien.eq.0) then
612  write(5,'(i10,1x,i3,1p,3(1x,e13.6))') nelem,j,
613  & (qfx(k,j,nelel),k=1,3)
614  else
615  do k=1,3
616  qfxl(k)=qfx(k,j,nelel)
617  enddo
618  call transformatrix(orab(1,iorien),coords(1,j),a)
619  write(5,'(i10,1x,i3,1p,3(1x,e13.6),1x,a20)') nelem,j,
620  & qfxl(1)*a(1,1)+qfxl(2)*a(2,1)+qfxl(3)*a(3,1),
621  & qfxl(1)*a(1,2)+qfxl(2)*a(2,2)+qfxl(3)*a(3,2),
622  & qfxl(1)*a(1,3)+qfxl(2)*a(2,3)+qfxl(3)*a(3,3),
623  & orname(iorien)(1:20)
624  endif
625  enddo
626  endif
627 !
628  return
subroutine shape6w(xi, et, ze, xl, xsj, shp, iflag)
Definition: shape6w.f:20
subroutine shape8q(xi, et, xl, xsj, xs, shp, iflag)
Definition: shape8q.f:20
subroutine shape10tet(xi, et, ze, xl, xsj, shp, iflag)
Definition: shape10tet.f:20
subroutine shape8h(xi, et, ze, xl, xsj, shp, iflag)
Definition: shape8h.f:20
subroutine shape15w(xi, et, ze, xl, xsj, shp, iflag)
Definition: shape15w.f:20
subroutine transformatrix(xab, p, a)
Definition: transformatrix.f:20
subroutine shape20h(xi, et, ze, xl, xsj, shp, iflag)
Definition: shape20h.f:20
subroutine beamintscheme(lakonl, mint3d, npropstart, prop, kk, xi, et, ze, weight)
Definition: beamintscheme.f:21
subroutine shape4tet(xi, et, ze, xl, xsj, shp, iflag)
Definition: shape4tet.f:20
subroutine shape6tri(xi, et, xl, xsj, xs, shp, iflag)
Definition: shape6tri.f:20
subroutine thickness(dgdx, nobject, nodedesiboun, ndesiboun, objectset, xo, yo, zo, x, y, z, nx, ny, nz, co, ifree, ndesia, ndesib, iobject, ndesi, dgdxglob, nk)
Definition: thickness.f:22
Hosted by OpenAircraft.com, (Michigan UAV, LLC)