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

Go to the source code of this file.

Functions/Subroutines

subroutine printoutfluid (set, nset, istartset, iendset, ialset, nprint, prlab, prset, ipkonf, lakonf, sti, eei, xstate, ener, mi, nstate_, co, konf, qfx, ttime, trab, inotr, ntrans, orab, ielorienf, norien, vold, ielmatf, thicke, eme, xturb, physcon, nactdoh, ielpropf, prop, xkappa, xmach, ithermal, orname)
 

Function/Subroutine Documentation

◆ printoutfluid()

subroutine printoutfluid ( 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,
integer, dimension(*)  ipkonf,
character*8, dimension(*)  lakonf,
real*8, dimension(6,mi(1),*)  sti,
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_,
real*8, dimension(3,*)  co,
integer, dimension(*)  konf,
real*8, dimension(3,mi(1),*)  qfx,
real*8  ttime,
real*8, dimension(7,*)  trab,
integer, dimension(2,*)  inotr,
integer  ntrans,
real*8, dimension(7,*)  orab,
integer, dimension(mi(3),*)  ielorienf,
integer  norien,
real*8, dimension(0:mi(2),*)  vold,
integer, dimension(mi(3),*)  ielmatf,
real*8, dimension(mi(3),*)  thicke,
real*8, dimension(6,mi(1),*)  eme,
real*8, dimension(2,*)  xturb,
real*8, dimension(*)  physcon,
integer, dimension(*)  nactdoh,
integer, dimension(*)  ielpropf,
real*8, dimension(*)  prop,
real*8, dimension(*)  xkappa,
real*8, dimension(*)  xmach,
integer, dimension(2)  ithermal,
character*80, dimension(*)  orname 
)
25 !
26 ! stores results in the .dat file
27 !
28  implicit none
29 !
30  character*6 prlab(*)
31  character*8 lakonf(*)
32  character*80 noset,elset,orname(*)
33  character*81 set(*),prset(*)
34 !
35  integer nset,istartset(*),iendset(*),ialset(*),nprint,ipkonf(*),
36  & mi(*),nstate_,ii,jj,iset,l,limit,node,ipos,nelel,ithermal(2),
37  & nelem,konf(*),inotr(2,*),ntrans,ielorienf(mi(3),*),norien,
38  & mt,ielmatf(mi(3),*),nactdoh(*),
39  & ielpropf(*)
40 !
41  real*8 sti(6,mi(1),*),xkappa(*),xmach(*),
42  & eei(6,mi(1),*),xstate(nstate_,mi(1),*),ener(mi(1),*),
43  & co(3,*),qfx(3,mi(1),*),ttime,
44  & trab(7,*),orab(7,*),vold(0:mi(2),*),thicke(mi(3),*),
45  & eme(6,mi(1),*),xturb(2,*),physcon(*),prop(*)
46 !
47  mt=mi(2)+1
48 !
49  do ii=1,nprint
50 !
51 ! nodal values
52 !
53  if((prlab(ii)(1:4).eq.'VF ').or.(prlab(ii)(1:4).eq.'PSF ').or.
54  & (prlab(ii)(1:4).eq.'TSF ').or.(prlab(ii)(1:4).eq.'PTF ').or.
55  & (prlab(ii)(1:4).eq.'TTF ').or.(prlab(ii)(1:4).eq.'CP ').or.
56  & (prlab(ii)(1:4).eq.'TURB').or.(prlab(ii)(1:4).eq.'MACH'))
57  & then
58 !
59  ipos=index(prset(ii),' ')
60  noset=' '
61  noset(1:ipos-1)=prset(ii)(1:ipos-1)
62 !
63 ! printing the header
64 !
65  if(prlab(ii)(1:4).eq.'VF ') then
66  write(5,*)
67  write(5,100) noset(1:ipos-2),ttime
68  100 format(' velocities (vx,vy,vz) for set ',a,
69  & ' and time ',e14.7)
70  write(5,*)
71  elseif(prlab(ii)(1:4).eq.'PSF ') then
72  write(5,*)
73  write(5,101) noset(1:ipos-2),ttime
74  101 format(' static pressures for set ',a,' and time ',e14.7)
75  write(5,*)
76  elseif(prlab(ii)(1:5).eq.'TSF ') then
77  write(5,*)
78  write(5,102) noset(1:ipos-2),ttime
79  102 format(' static temperatures for set ',a,
80  & ' and time ',e14.7)
81  write(5,*)
82  elseif(prlab(ii)(1:5).eq.'PTF ') then
83  write(5,*)
84  write(5,103) noset(1:ipos-2),ttime
85  103 format(' total pressures for set ',a,' and time ',e14.7)
86  write(5,*)
87  elseif(prlab(ii)(1:4).eq.'TTF ') then
88  write(5,*)
89  write(5,115) noset(1:ipos-2),ttime
90  115 format(' total temperatures for set ',a,' and time ',
91  & e14.7)
92  write(5,*)
93  elseif(prlab(ii)(1:4).eq.'CP ') then
94  write(5,*)
95  write(5,117) noset(1:ipos-2),ttime
96  117 format(' pressure coefficients for set ',
97  &a,' and time ',e14.7)
98  write(5,*)
99  elseif(prlab(ii)(1:4).eq.'TURB') then
100  write(5,*)
101  write(5,118) noset(1:ipos-2),ttime
102  118 format(' turbulence variables for set ',a,
103  & ' and time ',e14.7)
104  write(5,*)
105  elseif(prlab(ii)(1:4).eq.'MACH') then
106  write(5,*)
107  write(5,119) noset(1:ipos-2),ttime
108  119 format(' Mach numbers for set ',a,
109  & ' and time ',e14.7)
110  write(5,*)
111  endif
112 !
113 ! printing the data
114 !
115  do iset=1,nset
116  if(set(iset).eq.prset(ii)) exit
117  enddo
118  do jj=istartset(iset),iendset(iset)
119  if(ialset(jj).lt.0) cycle
120  if(jj.eq.iendset(iset)) then
121  node=ialset(jj)
122  call printoutnodefluid(prlab,vold,xturb,physcon,
123  & ii,node,trab,inotr,ntrans,co,mi,xkappa,xmach)
124  elseif(ialset(jj+1).gt.0) then
125  node=ialset(jj)
126  call printoutnodefluid(prlab,vold,xturb,physcon,
127  & ii,node,trab,inotr,ntrans,co,mi,xkappa,xmach)
128  else
129  do node=ialset(jj-1)-ialset(jj+1),ialset(jj),
130  & -ialset(jj+1)
131  call printoutnodefluid(prlab,vold,xturb,physcon,
132  & ii,node,trab,inotr,ntrans,co,mi,xkappa,xmach)
133  enddo
134  endif
135  enddo
136 !
137 ! integration point values
138 !
139  elseif((prlab(ii)(1:4).eq.'SVF ').or.
140  & (prlab(ii)(1:4).eq.'HFLF')) then
141 !
142  ipos=index(prset(ii),' ')
143  elset=' '
144  elset(1:ipos-1)=prset(ii)(1:ipos-1)
145 !
146  limit=1
147 !
148  do l=1,limit
149 !
150 ! printing the header
151 !
152  if(prlab(ii)(1:4).eq.'SVF ') then
153  write(5,*)
154  write(5,106) elset(1:ipos-2),ttime
155  106 format(' viscous stresses (elem, integ.pnt.,sxx,syy,sz
156  &z,sxy,sxz,syz) for set ',a,' and time ',e14.7)
157  write(5,*)
158  elseif(prlab(ii)(1:4).eq.'HFLF') then
159  write(5,*)
160  write(5,112) elset(1:ipos-2),ttime
161  112 format(' heat flux (elem, integ.pnt.,qx,qy,qz) for set
162  & ',a,' and time ',e14.7)
163  write(5,*)
164  endif
165 !
166 ! printing the data
167 !
168  do iset=1,nset
169  if(set(iset).eq.prset(ii)) exit
170  enddo
171  do jj=istartset(iset),iendset(iset)
172  if(ialset(jj).lt.0) cycle
173  if(jj.eq.iendset(iset)) then
174  nelem=ialset(jj)
175  nelel=nactdoh(nelem)
176  call printoutint(prlab,ipkonf,lakonf,sti,eei,
177  & xstate,ener,mi(1),nstate_,ii,nelem,qfx,
178  & orab,ielorienf,norien,co,konf,ielmatf,thicke,
179  & eme,ielpropf,prop,nelel,ithermal,
180  & orname)
181  elseif(ialset(jj+1).gt.0) then
182  nelem=ialset(jj)
183  nelel=nactdoh(nelem)
184  call printoutint(prlab,ipkonf,lakonf,sti,eei,
185  & xstate,ener,mi(1),nstate_,ii,nelem,qfx,orab,
186  & ielorienf,norien,co,konf,ielmatf,thicke,eme,
187  & ielpropf,prop,nelel,ithermal,
188  & orname)
189  else
190  do nelem=ialset(jj-1)-ialset(jj+1),ialset(jj),
191  & -ialset(jj+1)
192  nelel=nactdoh(nelem)
193  call printoutint(prlab,ipkonf,lakonf,sti,eei,
194  & xstate,ener,mi(1),nstate_,ii,nelem,
195  & qfx,orab,ielorienf,norien,co,konf,ielmatf,
196  & thicke,eme,ielpropf,prop,nelel,ithermal,
197  & orname)
198  enddo
199  endif
200  enddo
201 !
202  enddo
203  endif
204  enddo
205 !
206  return
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)
Definition: printoutint.f:22
subroutine printoutnodefluid(prlab, vold, xturb, physcon, ii, node, trab, inotr, ntrans, co, mi, xkappa, xmach)
Definition: printoutnodefluid.f:21
Hosted by OpenAircraft.com, (Michigan UAV, LLC)