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

Go to the source code of this file.

Functions/Subroutines

subroutine printoutelem (prlab, ipkon, lakon, kon, co, ener, mi, ii, nelem, energytot, volumetot, enerkintot, ne, stx, nodes, thicke, ielmat, ielem, iface, mortar, ielprop, prop, sideload, nload, nelemload, xload, bhetot)
 

Function/Subroutine Documentation

◆ printoutelem()

subroutine printoutelem ( character*6, dimension(*)  prlab,
integer, dimension(*)  ipkon,
character*8, dimension(*)  lakon,
integer, dimension(*)  kon,
real*8, dimension(3,*)  co,
real*8, dimension(mi(1),*)  ener,
integer, dimension(*)  mi,
integer  ii,
integer  nelem,
real*8  energytot,
real*8  volumetot,
real*8  enerkintot,
integer  ne,
real*8, dimension(6,mi(1),*)  stx,
integer  nodes,
real*8, dimension(mi(3),*)  thicke,
integer, dimension(mi(3),*)  ielmat,
integer  ielem,
integer  iface,
integer  mortar,
integer, dimension(*)  ielprop,
real*8, dimension(*)  prop,
character*20, dimension(*)  sideload,
integer  nload,
integer, dimension(2,*)  nelemload,
real*8, dimension(2,*)  xload,
real*8  bhetot 
)
23 !
24 ! stores whole element results for element "nelem" in the .dat file
25 !
26  implicit none
27 !
28  character*6 prlab(*)
29  character*8 lakon(*)
30  character*20 sideload(*)
31 !
32  integer ipkon(*),nelem,ii,kon(*),mi(*),nope,indexe,i,j,k,
33  & konl(20),iface,mortar,ielem,ielprop(*),nvol,nbhe,
34  & mint3d,jj,nener,iflag,nkin,ne,nodes,ki,kl,ilayer,nlayer,kk,
35  & nopes,ielmat(mi(3),*),mint2d,null,id,nload,nelemload(2,*)
36 !
37  real*8 ener(mi(1),*),energytot,volumetot,energy,volume,co(3,*),
38  & xl(3,20),xi,et,ze,xsj,shp(4,20),weight,enerkintot,enerkin,
39  & stx(6,mi(1),*),a,gs(8,4),dlayer(4),tlayer(4),thickness,
40  & thicke(mi(3),*),xlayer(mi(3),4),shp2(7,8),xs2(3,7),xsj2(3),
41  & xl2(3,8),prop(*),dflux,xload(2,*),bhe,bhetot
42 !
43  include "gauss.f"
44 !
45  data iflag /2/
46 !
47  if(ipkon(nelem).lt.0) return
48  indexe=ipkon(nelem)
49  null=0
50 !
51  nener=0
52  nkin=0
53  nvol=0
54  nbhe=0
55 !
56  if((prlab(ii)(1:4).eq.'ELSE').or.(prlab(ii)(1:4).eq.'CELS')) then
57  nener=1
58  elseif(prlab(ii)(1:4).eq.'ELKE') then
59  nkin=1
60  elseif(prlab(ii)(1:4).eq.'EVOL') then
61  nvol=1
62  elseif(prlab(ii)(1:4).eq.'EBHE') then
63  nbhe=1
64  endif
65 !
66 ! for contact displacements and stresses no integration has
67 ! to be performed
68 !
69  if((prlab(ii)(1:5).eq.'CDIS ').or.
70  & (prlab(ii)(1:5).eq.'CDIST')) then
71 !
72 ! contact displacements
73 !
74  if(mortar.eq.0) then
75  write(5,'(i10,1p,1x,e13.6,1p,1x,e13.6,1p,1x,e13.6)') nodes,
76  & stx(1,1,nelem),stx(2,1,nelem),stx(3,1,nelem)
77  elseif(mortar.eq.1) then
78  write(5,'(i10,1x,i10,1p,1x,e13.6,1p,1x,e13.6,1p,1x,e13.6)')
79  & ielem,iface,
80  & stx(1,1,nelem),stx(2,1,nelem),stx(3,1,nelem)
81  endif
82  return
83  elseif((prlab(ii)(1:5).eq.'CSTR ').or.
84  & (prlab(ii)(1:5).eq.'CSTRT')) then
85 !
86 ! contact stresses
87 !
88  if(mortar.eq.0) then
89  write(5,'(i10,1p,1x,e13.6,1p,1x,e13.6,1p,1x,e13.6)') nodes,
90  & stx(4,1,nelem),stx(5,1,nelem),stx(6,1,nelem)
91  elseif(mortar.eq.1) then
92  write(5,'(i10,1x,i10,1p,1x,e13.6,1p,1x,e13.6,1p,1x,e13.6)')
93  & ielem,iface,
94  & stx(4,1,nelem),stx(5,1,nelem),stx(6,1,nelem)
95  endif
96  return
97  elseif(prlab(ii)(1:4).eq.'EBHE ') then
98 !
99 ! body heating: check whether there is any body heating in
100 ! this elements
101 !
102  dflux=0.d0
103  if(nload.gt.0) then
104  call nident2(nelemload,nelem,nload,id)
105  do
106  if((id.eq.0).or.(nelemload(1,id).ne.nelem)) exit
107  if(sideload(id)(1:2).ne.'BF') then
108  id=id-1
109  cycle
110  endif
111  dflux=xload(1,id)
112  exit
113  enddo
114  endif
115 !
116 ! if no body heat flux: print and leave
117 !
118  if(dflux.eq.0.d0) then
119  if((prlab(ii)(1:5).eq.'EBHE ').or.
120  & (prlab(ii)(1:5).eq.'EBHET')) then
121  write(5,'(i10,1p,1x,e13.6)') nelem,dflux
122  endif
123  return
124  endif
125  endif
126 !
127  if(lakon(nelem)(1:5).eq.'C3D8I') then
128  nope=11
129  elseif(lakon(nelem)(4:4).eq.'2') then
130  nope=20
131  elseif(lakon(nelem)(4:4).eq.'8') then
132  nope=8
133  elseif(lakon(nelem)(4:5).eq.'10') then
134  nope=10
135  elseif(lakon(nelem)(4:4).eq.'4') then
136  nope=4
137  elseif(lakon(nelem)(4:5).eq.'15') then
138  nope=15
139  elseif(lakon(nelem)(4:5).eq.'6') then
140  nope=6
141  else
142  nope=0
143  endif
144 !
145 ! composite materials
146 !
147  if(lakon(nelem)(7:8).eq.'LC') then
148 !
149 ! determining the number of layers
150 !
151  nlayer=0
152  do k=1,mi(3)
153  if(ielmat(k,nelem).ne.0) then
154  nlayer=nlayer+1
155  endif
156  enddo
157 !
158  if(lakon(nelem)(4:4).eq.'2') then
159 !
160  mint2d=4
161  nopes=8
162 !
163 ! determining the layer thickness and global thickness
164 ! at the shell integration points
165 !
166  iflag=1
167  indexe=ipkon(nelem)
168  do kk=1,mint2d
169  xi=gauss3d2(1,kk)
170  et=gauss3d2(2,kk)
171  call shape8q(xi,et,xl2,xsj2,xs2,shp2,iflag)
172  tlayer(kk)=0.d0
173  do i=1,nlayer
174  thickness=0.d0
175  do j=1,nopes
176  thickness=thickness+thicke(i,indexe+j)*shp2(4,j)
177  enddo
178  tlayer(kk)=tlayer(kk)+thickness
179  xlayer(i,kk)=thickness
180  enddo
181  enddo
182  iflag=2
183  !
184  ilayer=0
185  do i=1,4
186  dlayer(i)=0.d0
187  enddo
188 
189  elseif(lakon(nelem)(4:5).eq.'15') then
190  !
191  mint2d=3
192  nopes=6
193  !
194  ! determining the layer thickness and global thickness
195  ! at the shell integration points
196  !
197  iflag=1
198  indexe=ipkon(nelem)
199  do kk=1,mint2d
200  xi=gauss3d10(1,kk)
201  et=gauss3d10(2,kk)
202  call shape6tri(xi,et,xl2,xsj2,xs2,shp2,iflag)
203  tlayer(kk)=0.d0
204  do i=1,nlayer
205  thickness=0.d0
206  do j=1,nopes
207  thickness=thickness+thicke(i,indexe+j)*shp2(4,j)
208  enddo
209  tlayer(kk)=tlayer(kk)+thickness
210  xlayer(i,kk)=thickness
211  enddo
212  enddo
213  iflag=2
214 !
215  ilayer=0
216  do i=1,3
217  dlayer(i)=0.d0
218  enddo
219  endif
220 !
221  endif
222 !
223  do j=1,nope
224  konl(j)=kon(indexe+j)
225  do k=1,3
226  xl(k,j)=co(k,konl(j))
227  enddo
228  enddo
229 !
230  energy=0.d0
231  volume=0.d0
232  bhe=0.d0
233  enerkin=0.d0
234 !
235  if(lakon(nelem)(4:5).eq.'8R') then
236  mint3d=1
237  elseif(lakon(nelem)(4:7).eq.'20RB') then
238  if((lakon(nelem)(8:8).eq.'R').or.
239  & (lakon(nelem)(8:8).eq.'C')) then
240  mint3d=50
241  else
242  call beamintscheme(lakon(nelem),mint3d,ielprop(nelem),prop,
243  & null,xi,et,ze,weight)
244  endif
245  elseif((lakon(nelem)(4:4).eq.'8').or.
246  & (lakon(nelem)(4:6).eq.'20R')) then
247  if(lakon(nelem)(7:8).eq.'LC') then
248  mint3d=8*nlayer
249  else
250  mint3d=8
251  endif
252  elseif(lakon(nelem)(4:4).eq.'2') then
253  mint3d=27
254  elseif(lakon(nelem)(4:5).eq.'10') then
255  mint3d=4
256  elseif(lakon(nelem)(4:4).eq.'4') then
257  mint3d=1
258  elseif(lakon(nelem)(4:5).eq.'15') then
259  if(lakon(nelem)(7:8).eq.'LC') then
260  mint3d=6*nlayer
261  else
262  mint3d=9
263  endif
264  elseif(lakon(nelem)(4:5).eq.'6') then
265  mint3d=2
266  else
267  if(nener.eq.1)then
268  energy=ener(1,nelem)
269  endif
270  mint3d=0
271  endif
272 !
273  do jj=1,mint3d
274  if(lakon(nelem)(4:5).eq.'8R') then
275  xi=gauss3d1(1,jj)
276  et=gauss3d1(2,jj)
277  ze=gauss3d1(3,jj)
278  weight=weight3d1(jj)
279  elseif(lakon(nelem)(4:7).eq.'20RB') then
280  if((lakon(nelem)(8:8).eq.'R').or.
281  & (lakon(nelem)(8:8).eq.'C')) then
282  xi=gauss3d13(1,kk)
283  et=gauss3d13(2,kk)
284  ze=gauss3d13(3,kk)
285  weight=weight3d13(kk)
286  else
287  call beamintscheme(lakon(nelem),mint3d,ielprop(nelem),
288  & prop,kk,xi,et,ze,weight)
289  endif
290  elseif((lakon(nelem)(4:4).eq.'8').or.
291  & (lakon(nelem)(4:6).eq.'20R'))
292  & then
293  if(lakon(nelem)(7:8).ne.'LC') then
294  xi=gauss3d2(1,jj)
295  et=gauss3d2(2,jj)
296  ze=gauss3d2(3,jj)
297  weight=weight3d2(jj)
298  else
299  kl=mod(jj,8)
300  if(kl.eq.0) kl=8
301 !
302  xi=gauss3d2(1,kl)
303  et=gauss3d2(2,kl)
304  ze=gauss3d2(3,kl)
305  weight=weight3d2(kl)
306 !
307  ki=mod(jj,4)
308  if(ki.eq.0) ki=4
309 !
310  if(kl.eq.1) then
311  ilayer=ilayer+1
312  if(ilayer.gt.1) then
313  do i=1,4
314  dlayer(i)=dlayer(i)+xlayer(ilayer-1,i)
315  enddo
316  endif
317  endif
318  ze=2.d0*(dlayer(ki)+(ze+1.d0)/2.d0*xlayer(ilayer,ki))/
319  & tlayer(ki)-1.d0
320  weight=weight*xlayer(ilayer,ki)/tlayer(ki)
321  endif
322  elseif(lakon(nelem)(4:4).eq.'2') then
323  xi=gauss3d3(1,jj)
324  et=gauss3d3(2,jj)
325  ze=gauss3d3(3,jj)
326  weight=weight3d3(jj)
327  elseif(lakon(nelem)(4:5).eq.'10') then
328  xi=gauss3d5(1,jj)
329  et=gauss3d5(2,jj)
330  ze=gauss3d5(3,jj)
331  weight=weight3d5(jj)
332  elseif(lakon(nelem)(4:4).eq.'4') then
333  xi=gauss3d4(1,jj)
334  et=gauss3d4(2,jj)
335  ze=gauss3d4(3,jj)
336  weight=weight3d4(jj)
337  elseif(lakon(nelem)(4:5).eq.'15') then
338  if(lakon(nelem)(7:8).ne.'LC') then
339  xi=gauss3d8(1,jj)
340  et=gauss3d8(2,jj)
341  ze=gauss3d8(3,jj)
342  weight=weight3d8(jj)
343  else
344  kl=mod(jj,6)
345  if(kl.eq.0) kl=6
346 !
347  xi=gauss3d10(1,kl)
348  et=gauss3d10(2,kl)
349  ze=gauss3d10(3,kl)
350  weight=weight3d10(kl)
351 !
352  ki=mod(jj,3)
353  if(ki.eq.0) ki=3
354 !
355  if(kl.eq.1) then
356  ilayer=ilayer+1
357  if(ilayer.gt.1) then
358  do i=1,3
359  dlayer(i)=dlayer(i)+xlayer(ilayer-1,i)
360  enddo
361  endif
362  endif
363  ze=2.d0*(dlayer(ki)+(ze+1.d0)/2.d0*xlayer(ilayer,ki))/
364  & tlayer(ki)-1.d0
365  weight=weight*xlayer(ilayer,ki)/tlayer(ki)
366  endif
367  else
368  xi=gauss3d7(1,jj)
369  et=gauss3d7(2,jj)
370  ze=gauss3d7(3,jj)
371  weight=weight3d7(jj)
372  endif
373 !
374  if(lakon(nelem)(1:5).eq.'C3D8R') then
375  call shape8hr(xl,xsj,shp,gs,a)
376  elseif(lakon(nelem)(1:5).eq.'C3D8I') then
377  call shape8hu(xi,et,ze,xl,xsj,shp,iflag)
378  elseif(nope.eq.20) then
379  call shape20h(xi,et,ze,xl,xsj,shp,iflag)
380  elseif(nope.eq.8) then
381  call shape8h(xi,et,ze,xl,xsj,shp,iflag)
382  elseif(nope.eq.10) then
383  call shape10tet(xi,et,ze,xl,xsj,shp,iflag)
384  elseif(nope.eq.4) then
385  call shape4tet(xi,et,ze,xl,xsj,shp,iflag)
386  elseif(nope.eq.15) then
387  call shape15w(xi,et,ze,xl,xsj,shp,iflag)
388  else
389  call shape6w(xi,et,ze,xl,xsj,shp,iflag)
390  endif
391 !
392  if(nener.eq.1) then
393  energy=energy+weight*xsj*ener(jj,nelem)
394  elseif(nkin.eq.1) then
395  enerkin=enerkin+weight*xsj*ener(jj,nelem+ne)
396  elseif(nvol.eq.1) then
397  volume=volume+weight*xsj
398  elseif(nbhe.eq.1) then
399  bhe=bhe+dflux*weight*xsj
400  endif
401  enddo
402 !
403  if(nener.eq.1) then
404  energytot=energytot+energy
405  elseif(nkin.eq.1) then
406  enerkintot=enerkintot+enerkin
407  elseif(nvol.eq.1) then
408  volumetot=volumetot+volume
409  elseif(nbhe.eq.1) then
410  bhetot=bhetot+bhe
411  endif
412 !
413 ! writing to file
414 !
415  if((prlab(ii)(1:5).eq.'ELSE ').or.
416  & (prlab(ii)(1:5).eq.'ELSET')) then
417  write(5,'(i10,1p,1x,e13.6)') nelem,energy
418  elseif((prlab(ii)(1:5).eq.'CELS ').or.
419  & (prlab(ii)(1:5).eq.'CELST')) then
420  if(mortar.eq.0) then
421  write(5,'(i10,1p,1x,e13.6)') nodes,energy
422  elseif(mortar.eq.1) then
423  write(5,'(i10,1x,i10,1p,1x,e13.6)') ielem,iface,energy
424  endif
425  elseif((prlab(ii)(1:5).eq.'EVOL ').or.
426  & (prlab(ii)(1:5).eq.'EVOLT')) then
427  write(5,'(i10,1p,1x,e13.6)') nelem,volume
428  elseif((prlab(ii)(1:5).eq.'ELKE ').or.
429  & (prlab(ii)(1:5).eq.'ELKET')) then
430  write(5,'(i10,1p,1x,e13.6)') nelem,enerkin
431  elseif((prlab(ii)(1:5).eq.'EBHE ').or.
432  & (prlab(ii)(1:5).eq.'EBHET')) then
433  write(5,'(i10,1p,1x,e13.6)') nelem,bhe
434  endif
435 !
436  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 nident2(x, px, n, id)
Definition: nident2.f:27
subroutine shape15w(xi, et, ze, xl, xsj, shp, iflag)
Definition: shape15w.f:20
subroutine shape20h(xi, et, ze, xl, xsj, shp, iflag)
Definition: shape20h.f:20
subroutine shape8hr(xl, xsj, shp, gs, a)
Definition: shape8hr.f:20
subroutine beamintscheme(lakonl, mint3d, npropstart, prop, kk, xi, et, ze, weight)
Definition: beamintscheme.f:21
subroutine nodes(inpc, textpart, co, nk, nk_, set, istartset, iendset, ialset, nset, nset_, nalset, nalset_, istep, istat, n, iline, ipol, inl, ipoinp, inp, ipoinpc)
Definition: nodes.f:22
subroutine shape4tet(xi, et, ze, xl, xsj, shp, iflag)
Definition: shape4tet.f:20
subroutine shape8hu(xi, et, ze, xl, xsj, shp, iflag)
Definition: shape8hu.f:20
subroutine dflux(flux, sol, kstep, kinc, time, noel, npt, coords, jltyp, temp, press, loadtype, area, vold, co, lakonl, konl, ipompc, nodempc, coefmpc, nmpc, ikmpc, ilmpc, iscale, mi, sti, xstateini, xstate, nstate_, dtime)
Definition: dflux.f:23
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)