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

Go to the source code of this file.

Functions/Subroutines

subroutine resultsprint (co, nk, kon, ipkon, lakon, ne, v, stn, inum, stx, ielorien, norien, orab, t1, ithermal, filab, een, iperturb, fn, nactdof, iout, vold, nodeboun, ndirboun, nboun, nmethod, ttime, xstate, epn, mi, nstate_, ener, enern, xstaten, eei, set, nset, istartset, iendset, ialset, nprint, prlab, prset, qfx, qfn, trab, inotr, ntrans, nelemload, nload, ikin, ielmat, thicke, eme, emn, rhcon, nrhcon, shcon, nshcon, cocon, ncocon, ntmat_, sideload, icfd, inomat, pslavsurf, islavact, cdn, mortar, islavnode, nslavnode, ntie, islavsurf, time, ielprop, prop, veold, ne0, nmpc, ipompc, nodempc, labmpc, energyini, energy, orname, xload)
 

Function/Subroutine Documentation

◆ resultsprint()

subroutine resultsprint ( real*8, dimension(3,*)  co,
integer  nk,
integer, dimension(*)  kon,
integer, dimension(*)  ipkon,
character*8, dimension(*)  lakon,
integer  ne,
real*8, dimension(0:mi(2),*)  v,
real*8, dimension(6,*)  stn,
integer, dimension(*)  inum,
real*8, dimension(6,mi(1),*)  stx,
integer, dimension(mi(3),*)  ielorien,
integer  norien,
real*8, dimension(7,*)  orab,
real*8, dimension(*)  t1,
integer, dimension(2)  ithermal,
character*87, dimension(*)  filab,
real*8, dimension(6,*)  een,
integer, dimension(*)  iperturb,
real*8, dimension(0:mi(2),*)  fn,
integer, dimension(0:mi(2),*)  nactdof,
integer  iout,
real*8, dimension(0:mi(2),*)  vold,
integer, dimension(*)  nodeboun,
integer, dimension(*)  ndirboun,
integer  nboun,
integer  nmethod,
real*8  ttime,
real*8, dimension(nstate_,mi(1),*)  xstate,
real*8, dimension(*)  epn,
integer, dimension(*)  mi,
integer  nstate_,
real*8, dimension(mi(1),*)  ener,
real*8, dimension(*)  enern,
real*8, dimension(nstate_,*)  xstaten,
real*8, dimension(6,mi(1),*)  eei,
character*81, dimension(*)  set,
integer  nset,
integer, dimension(*)  istartset,
integer, dimension(*)  iendset,
integer, dimension(*)  ialset,
integer  nprint,
character*6, dimension(*)  prlab,
character*81, dimension(*)  prset,
real*8, dimension(3,mi(1),*)  qfx,
real*8, dimension(3,*)  qfn,
real*8, dimension(7,*)  trab,
integer, dimension(2,*)  inotr,
integer  ntrans,
integer, dimension(2,*)  nelemload,
integer  nload,
integer  ikin,
integer, dimension(mi(3),*)  ielmat,
real*8, dimension(mi(3),*)  thicke,
real*8, dimension(6,mi(1),*)  eme,
real*8, dimension(6,*)  emn,
real*8, dimension(0:1,ntmat_,*)  rhcon,
integer, dimension(*)  nrhcon,
real*8, dimension(0:3,ntmat_,*)  shcon,
integer, dimension(*)  nshcon,
real*8, dimension(0:6,ntmat_,*)  cocon,
integer, dimension(2,*)  ncocon,
integer  ntmat_,
character*20, dimension(*)  sideload,
integer  icfd,
integer, dimension(*)  inomat,
real*8, dimension(3,*)  pslavsurf,
integer, dimension(*)  islavact,
real*8, dimension(6,*)  cdn,
integer  mortar,
integer, dimension(*)  islavnode,
integer, dimension(*)  nslavnode,
integer  ntie,
integer, dimension(2,*)  islavsurf,
real*8  time,
integer, dimension(*)  ielprop,
real*8, dimension(*)  prop,
real*8, dimension(0:mi(2),*)  veold,
integer  ne0,
integer  nmpc,
integer, dimension(*)  ipompc,
integer, dimension(3,*)  nodempc,
character*20, dimension(*)  labmpc,
real*8, dimension(*)  energyini,
real*8, dimension(*)  energy,
character*80, dimension(*)  orname,
real*8, dimension(2,*)  xload 
)
29 !
30 ! - stores the results in the .dat file, if requested
31 ! - nodal quantities at the nodes
32 ! - element quantities at the integration points
33 ! - calculates the extrapolation of element quantities to
34 ! the nodes (if requested for .frd output)
35 ! - calculates 1d/2d results for 1d/2d elements by
36 ! interpolation
37 !
38  implicit none
39 !
40  logical force,rfprint
41 !
42  character*1 cflag
43  character*6 prlab(*)
44  character*8 lakon(*)
45  character*20 sideload(*),labmpc(*)
46  character*80 orname(*)
47  character*81 set(*),prset(*)
48  character*87 filab(*)
49 !
50  integer kon(*),inum(*),iperm(20),mi(*),ielorien(mi(3),*),
51  & ipkon(*),nactdof(0:mi(2),*),nodeboun(*),icompressible,
52  & nelemload(2,*),ndirboun(*),ielmat(mi(3),*),nrhcon(*),
53  & inotr(2,*),iorienloc,iflag,nload,mt,nk,ne,ithermal(2),i,
54  & norien,iperturb(*),iout,nboun,nmethod,node,nshcon(*),
55  & nfield,ndim,nstate_,nset,istartset(*),iendset(*),ialset(*),
56  & nprint,ntrans,ikin,ncocon(2,*),ntmat_,icfd,inomat(*),mortar,
57  & islavact(*),islavnode(*),nslavnode(*),ntie,islavsurf(2,*),
58  & ielprop(*),ne0,index,nmpc,ipompc(*),nodempc(3,*),nactdoh,
59  & iextrapolate
60 !
61  real*8 co(3,*),v(0:mi(2),*),stx(6,mi(1),*),stn(6,*),cdn(6,*),
62  & qfx(3,mi(1),*),qfn(3,*),orab(7,*),fn(0:mi(2),*),pslavsurf(3,*),
63  & t1(*),een(6,*),vold(0:mi(2),*),epn(*),thicke(mi(3),*),time,
64  & ener(mi(1),*),enern(*),eei(6,mi(1),*),rhcon(0:1,ntmat_,*),
65  & ttime,xstate(nstate_,mi(1),*),trab(7,*),xstaten(nstate_,*),
66  & eme(6,mi(1),*),emn(6,*),shcon(0:3,ntmat_,*),cocon(0:6,ntmat_,*),
67  & prop(*),veold(0:mi(2),*),energy(*),energyini(*),xload(2,*)
68 !
69  data iflag /3/
70  data iperm /5,6,7,8,1,2,3,4,13,14,15,16,9,10,11,12,17,18,19,20/
71 !
72  mt=mi(2)+1
73  iextrapolate=0
74 !
75 ! no print requests
76 !
77  if(iout.le.0) then
78 !
79 ! 2d basic dof results (displacements, temperature) are
80 ! calculated in each iteration, so that they are available
81 ! in the user subroutines
82 !
83  if(filab(1)(5:5).ne.' ') then
84  nfield=mt
85  call map3dto1d2d_v(v,ipkon,inum,kon,lakon,nfield,nk,
86  & ne,nactdof)
87  endif
88 !
89 ! the total energy should not be calculated:
90 ! - for non-dynamical calculations (nmethod!=4)
91 ! - for modal dynamics (iperturb(1)<=1)
92 ! - for thermal and thermomechanical calculations (ithermal(1)>1)
93 ! - for electromagnetic calculations (mi(2)=5)
94 !
95  if((nmethod.eq.4).and.(iperturb(1).gt.1).and.
96  & (ithermal(1).le.1).and.(mi(2).ne.5)) then
97  call calcenergy(ipkon,lakon,kon,co,ener,mi,ne,thicke,
98  & ielmat,energyini,energy,ielprop,prop)
99  endif
100 !
101  return
102  endif
103 !
104 ! output in dat file (with *NODE PRINT or *EL PRINT)
105 !
106  call printout(set,nset,istartset,iendset,ialset,nprint,
107  & prlab,prset,v,t1,fn,ipkon,lakon,stx,eei,xstate,ener,
108  & mi(1),nstate_,ithermal,co,kon,qfx,ttime,trab,inotr,ntrans,
109  & orab,ielorien,norien,nk,ne,inum,filab,vold,ikin,ielmat,thicke,
110  & eme,islavsurf,mortar,time,ielprop,prop,veold,orname,
111  & nelemload,nload,sideload,xload)
112 !
113 ! for facial information (*section print): if forces and/or
114 ! moments in sections are requested, the stresses have to be
115 ! extrapolated from the integration points to the nodes first
116 !
117  do i=1,nprint
118  if(prlab(i)(1:3).eq.'SOF') then
119  nfield=6
120  ndim=6
121  if((norien.gt.0).and.(filab(3)(6:6).eq.'L')) then
122  iorienloc=1
123  else
124  iorienloc=0
125  endif
126  cflag=filab(3)(5:5)
127  force=.false.
128 !
129  call extrapolate(stx,stn,ipkon,inum,kon,lakon,nfield,nk,
130  & ne,mi(1),ndim,orab,ielorien,co,iorienloc,cflag,
131  & vold,force,ielmat,thicke,ielprop,prop)
132  iextrapolate=1
133  exit
134  endif
135  enddo
136 !
137  icompressible=0
138  call printoutface(co,rhcon,nrhcon,ntmat_,vold,shcon,nshcon,
139  & cocon,ncocon,icompressible,istartset,iendset,ipkon,lakon,kon,
140  & ialset,prset,ttime,nset,set,nprint,prlab,ielmat,mi,
141  & ithermal,nactdoh,icfd,time,stn)
142 !
143 ! interpolation in the original nodes of 1d and 2d elements
144 ! this operation has to be performed in any case since
145 ! the interpolated values may be needed as boundary conditions
146 ! in the next step (e.g. the temperature in a heat transfer
147 ! calculation as boundary condition in a subsequent static
148 ! step)
149 !
150  if(filab(1)(5:5).ne.' ') then
151  nfield=mt
152  cflag=filab(1)(5:5)
153  force=.false.
154  call map3dto1d2d(v,ipkon,inum,kon,lakon,nfield,nk,
155  & ne,cflag,co,vold,force,mi)
156  endif
157 !
158  if((filab(2)(1:4).eq.'NT ').and.(ithermal(1).le.1)) then
159  if(filab(2)(5:5).eq.'I') then
160  nfield=1
161  cflag=filab(2)(5:5)
162  force=.false.
163  call map3dto1d2d(t1,ipkon,inum,kon,lakon,nfield,nk,
164  & ne,cflag,co,vold,force,mi)
165  endif
166  endif
167 !
168 ! check whether forces are requested in the frd-file. If so, but
169 ! none are requested in the .dat file, and output=2d,
170 ! map3dto1d2d has to be called
171 !
172  if(filab(5)(1:2).eq.'RF') then
173  if(filab(5)(5:5).eq.'I') then
174  rfprint=.false.
175  do i=1,nprint
176  if(prlab(i)(1:2).eq.'RF') then
177  rfprint=.true.
178  exit
179  endif
180  enddo
181  if(.not.rfprint) then
182  nfield=mt
183  cflag=' '
184  force=.true.
185  call map3dto1d2d(fn,ipkon,inum,kon,lakon,nfield,nk,
186  & ne,cflag,co,vold,force,mi)
187  endif
188  endif
189  endif
190 !
191 ! for composites:
192 ! interpolation of the displacements and temperatures
193 ! from the expanded nodes to the layer nodes
194 !
195  if(mi(3).gt.1) then
196  if((filab(1)(1:3).eq.'U ').or.
197  & ((filab(2)(1:4).eq.'NT ').and.(ithermal(1).gt.1))) then
198  nfield=mt
199  call map3dtolayer(v,ipkon,kon,lakon,nfield,
200  & ne,co,ielmat,mi)
201  endif
202  if((filab(2)(1:4).eq.'NT ').and.(ithermal(1).le.1)) then
203  nfield=1
204  call map3dtolayer(t1,ipkon,kon,lakon,nfield,
205  & ne,co,ielmat,mi)
206  endif
207  endif
208 !
209 ! determining the contact differential displacements and stresses
210 ! in the contact nodes for output in frd format (only for face-
211 ! to-face penalty; for node-to-face penalty these quantities are
212 ! determined in the slave nodes and no extrapolation is necessary)
213 !
214 ! This block must precede all calls to extrapolate, since the
215 ! field inum from extrapolatecontact.f is not correct; by a
216 ! subsequent call to extrapolate inum is corrected.
217 !
218  if((filab(26)(1:4).eq.'CONT').or.(filab(46)(1:4).eq.'PCON')) then
219  if(mortar.eq.1) then
220  nfield=6
221  ndim=6
222  cflag=filab(3)(5:5)
223  force=.false.
224  call extrapolatecontact(stx,cdn,ipkon,inum,kon,lakon,nfield,
225  & nk,ne,mi(1),ndim,co,cflag,vold,force,pslavsurf,
226  & islavact,islavnode,nslavnode,ntie,islavsurf,ielprop,prop,
227  & ielmat,ne0)
228  endif
229  endif
230 !
231 ! determining the stresses in the nodes for output in frd format
232 !
233  if((filab(3)(1:4).eq.'S ').or.(filab(18)(1:4).eq.'PHS ').or.
234  & (filab(20)(1:4).eq.'MAXS').or.
235  & (((filab(44)(1:4).eq.'EMFE').or.(filab(45)(1:4).eq.'EMFB'))
236  & .and.(ithermal(1).ne.2))) then
237  nfield=6
238  ndim=6
239  if((norien.gt.0).and.(filab(3)(6:6).eq.'L')) then
240  iorienloc=1
241  else
242  iorienloc=0
243  endif
244  cflag=filab(3)(5:5)
245  force=.false.
246 !
247  call extrapolate(stx,stn,ipkon,inum,kon,lakon,nfield,nk,
248  & ne,mi(1),ndim,orab,ielorien,co,iorienloc,cflag,
249  & vold,force,ielmat,thicke,ielprop,prop)
250  iextrapolate=1
251 !
252  endif
253 !
254 ! determining the total strains in the nodes for output in frd format
255 !
256  if((filab(4)(1:4).eq.'E ').or.(filab(30)(1:4).eq.'MAXE')) then
257  nfield=6
258  ndim=6
259  if((norien.gt.0).and.(filab(4)(6:6).eq.'L')) then
260  iorienloc=1
261  else
262  iorienloc=0
263  endif
264  cflag=filab(4)(5:5)
265  force=.false.
266  call extrapolate(eei,een,ipkon,inum,kon,lakon,nfield,nk,
267  & ne,mi(1),ndim,orab,ielorien,co,iorienloc,cflag,
268  & vold,force,ielmat,thicke,ielprop,prop)
269  iextrapolate=1
270  endif
271 !
272 ! determining the mechanical strains in the nodes for output in
273 ! frd format
274 !
275  if(filab(32)(1:4).eq.'ME ') then
276  nfield=6
277  ndim=6
278  if((norien.gt.0).and.(filab(4)(6:6).eq.'L')) then
279  iorienloc=1
280  else
281  iorienloc=0
282  endif
283  cflag=filab(4)(5:5)
284  force=.false.
285  call extrapolate(eme,emn,ipkon,inum,kon,lakon,nfield,nk,
286  & ne,mi(1),ndim,orab,ielorien,co,iorienloc,cflag,
287  & vold,force,ielmat,thicke,ielprop,prop)
288  iextrapolate=1
289  endif
290 !
291 ! determining the plastic equivalent strain in the nodes
292 ! for output in frd format
293 !
294  if(filab(6)(1:4).eq.'PEEQ') then
295  nfield=1
296  ndim=nstate_
297  iorienloc=0
298  cflag=filab(6)(5:5)
299  force=.false.
300  call extrapolate(xstate,epn,ipkon,inum,kon,lakon,nfield,nk,
301  & ne,mi(1),ndim,orab,ielorien,co,iorienloc,cflag,
302  & vold,force,ielmat,thicke,ielprop,prop)
303  iextrapolate=1
304  endif
305 !
306 ! determining the total energy in the nodes
307 ! for output in frd format
308 !
309  if(filab(7)(1:4).eq.'ENER') then
310  nfield=1
311  ndim=1
312  iorienloc=0
313  cflag=filab(7)(5:5)
314  force=.false.
315  call extrapolate(ener,enern,ipkon,inum,kon,lakon,nfield,nk,
316  & ne,mi(1),ndim,orab,ielorien,co,iorienloc,cflag,
317  & vold,force,ielmat,thicke,ielprop,prop)
318  iextrapolate=1
319  endif
320 !
321 ! determining the internal state variables in the nodes
322 ! for output in frd format
323 !
324  if(filab(8)(1:4).eq.'SDV ') then
325  nfield=nstate_
326  ndim=nstate_
327  if((norien.gt.0).and.(filab(9)(6:6).eq.'L')) then
328  write(*,*) '*WARNING in results: SDV variables cannot'
329  write(*,*) ' be stored in a local frame;'
330  write(*,*) ' the global frame will be used'
331  endif
332  iorienloc=0
333  cflag=filab(8)(5:5)
334  force=.false.
335  call extrapolate(xstate,xstaten,ipkon,inum,kon,lakon,nfield,nk,
336  & ne,mi(1),ndim,orab,ielorien,co,iorienloc,cflag,
337  & vold,force,ielmat,thicke,ielprop,prop)
338  iextrapolate=1
339  endif
340 !
341 ! determining the heat flux in the nodes for output in frd format
342 !
343  if(((filab(9)(1:4).eq.'HFL ').and.(ithermal(1).gt.1)).or.
344  & ((filab(42)(1:3).eq.'ECD').and.(ithermal(1).eq.2))) then
345  nfield=3
346  ndim=3
347  if((norien.gt.0).and.(filab(9)(6:6).eq.'L')) then
348  iorienloc=1
349  else
350  iorienloc=0
351  endif
352  cflag=filab(9)(5:5)
353  force=.false.
354  call extrapolate(qfx,qfn,ipkon,inum,kon,lakon,nfield,nk,
355  & ne,mi(1),ndim,orab,ielorien,co,iorienloc,cflag,
356  & vold,force,ielmat,thicke,ielprop,prop)
357  iextrapolate=1
358  endif
359 !
360 ! if no element quantities requested in the nodes: calculate
361 ! inum if nodal quantities are requested: used in subroutine frd
362 ! to determine which nodes are active in the model
363 !
364 c if((filab(3)(1:4).ne.'S ').and.(filab(4)(1:4).ne.'E ').and.
365 c & (filab(6)(1:4).ne.'PEEQ').and.(filab(7)(1:4).ne.'ENER').and.
366 c & (filab(8)(1:4).ne.'SDV ').and.(filab(9)(1:4).ne.'HFL ').and.
367 c & (filab(42)(1:3).ne.'ECD').and.(filab(32)(1:4).ne.'ME ').and.
368 c & ((nmethod.ne.4).or.(iperturb(1).ge.2))) then
369 c
370  if((iextrapolate.eq.0).and.
371  & ((nmethod.ne.4).or.(iperturb(1).ge.2))) then
372 !
373  nfield=0
374  ndim=0
375  iorienloc=0
376  cflag=filab(1)(5:5)
377  call createinum(ipkon,inum,kon,lakon,nk,ne,cflag,nelemload,
378  & nload,nodeboun,nboun,ndirboun,ithermal,co,vold,mi,ielmat)
379  endif
380 !
381 c if(ithermal(1).gt.1) then
382  if(ithermal(2).gt.1) then
383 !
384 ! next section is executed if at least one step is thermal
385 ! or thermomechanical
386 !
387 ! extrapolation for the network
388 ! -interpolation for the total pressure and temperature
389 ! in the middle nodes
390 ! -extrapolation for the mass flow in the end nodes
391 !
392  call networkextrapolate(v,ipkon,inum,kon,lakon,ne,mi)
393 !
394 ! printing values for environmental film and
395 ! pressure nodes (these nodes are considered to be network
396 ! nodes)
397 !
398  do i=1,nload
399  if((sideload(i)(3:4).ne.'FC').and.
400  & (sideload(i)(3:4).ne.'NP')) cycle
401  node=nelemload(2,i)
402  if(icfd.eq.1) then
403  if(node.gt.0) then
404  if(inomat(node).ne.0) cycle
405  endif
406  endif
407  if((node.gt.0).and.(sideload(i)(1:1).ne.' ')) then
408  if(inum(node).lt.0) cycle
409  inum(node)=-1
410  endif
411  enddo
412 !
413 ! printing values radiation
414 ! (these nodes are considered to be network nodes, unless
415 ! they were already assigned to the structure)
416 !
417  do i=1,nload
418  if((sideload(i)(3:4).ne.'CR')) cycle
419  node=nelemload(2,i)
420  if(icfd.eq.1) then
421  if(node.gt.0) then
422  if(inomat(node).ne.0) cycle
423  endif
424  endif
425  if((node.gt.0).and.(sideload(i)(1:1).ne.' ')) then
426  if(inum(node).ne.0) cycle
427  inum(node)=-1
428  endif
429  enddo
430 !
431 ! printing values for nodes belonging to network MPC's
432 ! (these nodes are considered to be network nodes)
433 !
434  do i=1,nmpc
435  if(labmpc(i)(1:7).eq.'NETWORK') then
436  index=ipompc(i)
437  do
438  node=nodempc(1,index)
439  if(inum(node).ge.0) inum(node)=-1
440  index=nodempc(3,index)
441  if(index.eq.0) exit
442  enddo
443  endif
444  enddo
445 !
446 ! printing values of prescribed boundary conditions (these
447 ! nodes are considered to be network nodes)
448 !
449  do i=1,nboun
450  node=nodeboun(i)
451  if(inum(node).ne.0) cycle
452  if(icfd.eq.1) then
453  if(inomat(node).ne.0) cycle
454  endif
455  if((cflag.ne.' ').and.(ndirboun(i).eq.3)) cycle
456  inum(node)=-1
457  enddo
458  endif
459 !
460  return
subroutine networkextrapolate(v, ipkon, inum, kon, lakon, ne, mi)
Definition: networkextrapolate.f:21
subroutine extrapolate(yi, yn, ipkon, inum, kon, lakon, nfield, nk, ne, mi, ndim, orab, ielorien, co, iorienloc, cflag, vold, force, ielmat, thicke, ielprop, prop)
Definition: extrapolate.f:22
subroutine extrapolatecontact(yi, yn, ipkon, inum, kon, lakon, nfield, nk, ne, mi, ndim, co, cflag, vold, force, pslavsurf, islavact, islavnode, nslavnode, ntie, islavsurf, ielprop, prop, ielmat, ne0)
Definition: extrapolatecontact.f:22
subroutine map3dtolayer(yn, ipkon, kon, lakon, nfield, ne, co, ielmat, mi)
Definition: map3dtolayer.f:21
subroutine printout(set, nset, istartset, iendset, ialset, nprint, prlab, prset, v, t1, fn, ipkon, lakon, stx, eei, xstate, ener, mi, nstate_, ithermal, co, kon, qfx, ttime, trab, inotr, ntrans, orab, ielorien, norien, nk, ne, inum, filab, vold, ikin, ielmat, thicke, eme, islavsurf, mortar, time, ielprop, prop, veold, orname, nelemload, nload, sideload, xload)
Definition: printout.f:25
subroutine map3dto1d2d(yn, ipkon, inum, kon, lakon, nfield, nk, ne, cflag, co, vold, force, mi)
Definition: map3dto1d2d.f:21
subroutine printoutface(co, rhcon, nrhcon, ntmat_, vold, shcon, nshcon, cocon, ncocon, icompressible, istartset, iendset, ipkonf, lakonf, konf, ialset, prset, ttime, nset, set, nprint, prlab, ielmat, mi, ithermal, nactdoh, icfd, time, stn)
Definition: printoutface.f:23
subroutine createinum(ipkon, inum, kon, lakon, nk, ne, cflag, nelemload, nload, nodeboun, nboun, ndirboun, ithermal, co, vold, mi, ielmat)
Definition: createinum.f:21
subroutine map3dto1d2d_v(yn, ipkon, inum, kon, lakon, nfield, nk, ne, nactdof)
Definition: map3dto1d2d_v.f:21
subroutine calcenergy(ipkon, lakon, kon, co, ener, mi, ne, thicke, ielmat, energyini, energy, ielprop, prop)
Definition: calcenergy.f:21
Hosted by OpenAircraft.com, (Michigan UAV, LLC)