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

Go to the source code of this file.

Functions/Subroutines

subroutine massflow_percent (node1, node2, nodem, nelem, lakon, kon, ipkon, nactdog, identity, ielprop, prop, iflag, v, xflow, f, nodef, idirf, df, cp, r, physcon, dvi, numf, set, shcon, nshcon, rhcon, nrhcon, ntmat_, co, vold, mi, ttime, time, iaxial)
 

Function/Subroutine Documentation

◆ massflow_percent()

subroutine massflow_percent ( integer  node1,
integer  node2,
integer  nodem,
integer  nelem,
character*8, dimension(*)  lakon,
integer, dimension(*)  kon,
integer, dimension(*)  ipkon,
integer, dimension(0:3,*)  nactdog,
logical  identity,
integer, dimension(*)  ielprop,
real*8, dimension(*)  prop,
integer  iflag,
real*8, dimension(0:mi(2),*)  v,
real*8  xflow,
real*8  f,
integer, dimension(*)  nodef,
integer, dimension(*)  idirf,
real*8, dimension(*)  df,
real*8  cp,
real*8  r,
real*8, dimension(*)  physcon,
real*8  dvi,
integer  numf,
character*81, dimension(*)  set,
real*8, dimension(0:3,ntmat_,*)  shcon,
integer, dimension(*)  nshcon,
real*8, dimension(0:1,ntmat_,*)  rhcon,
integer, dimension(*)  nrhcon,
integer  ntmat_,
real*8, dimension(3,*)  co,
real*8, dimension(0:mi(2),*)  vold,
integer, dimension(*)  mi,
real*8  ttime,
real*8  time,
integer  iaxial 
)
23 !
24 ! partial massflow element
25 !
26 ! author: Yannick Muller
27 !
28  implicit none
29 !
30  logical identity
31  character*8 lakon(*)
32  character*81 set(*)
33 !
34  integer nelem,nactdog(0:3,*),node1,node2,nodem,numf,
35  & ielprop(*),nodef(*),idirf(*),index,iflag,
36  & inv,ipkon(*),kon(*),number,kgas,iaxial,
37  & nodea,nodeb,mi(*),i,itype,nodemup,
38  & nrhcon(*),ntmat_,nshcon(*)
39 !
40  real*8 prop(*),v(0:mi(2),*),xflow,f,df(*),kappa,r,a,d,xl,
41  & p1,p2,t1,physcon(*),pi,xflow_oil,t2,co(3,*),vold(0:mi(2),*),
42  & xflow_sum,percent_xflow,cp,dvi,pt1,pt2,tt1,tt2,ttime,time,
43  & shcon(0:3,ntmat_,*),rhcon(0:1,ntmat_,*)
44 !
45  pi=4.d0*datan(1.d0)
46  index=ielprop(nelem)
47 !
48  if(iflag.eq.0) then
49  identity=.true.
50 !
51  if(nactdog(2,node1).ne.0)then
52  identity=.false.
53  elseif(nactdog(2,node2).ne.0)then
54  identity=.false.
55  elseif(nactdog(1,nodem).ne.0)then
56  identity=.false.
57  endif
58 !
59  elseif(iflag.eq.1)then
60 !
61  percent_xflow=prop(index+1)
62  xflow_sum=0
63 !
64  do i=2,10
65  if(nint(prop(index+i)).ne.0) then
66  nodemup=kon(ipkon(nint(prop(index+i)))+2)
67  if(v(1,nodemup).gt.0)then
68  xflow_sum=xflow_sum+v(1,nodemup)*iaxial
69  endif
70  endif
71  enddo
72 !
73  if(xflow_sum.eq.0d0) then
74  xflow_sum=0.001d0
75  endif
76 !
77  xflow=xflow_sum*percent_xflow
78 !
79  elseif((iflag.eq.2).or.(iflag.eq.3))then
80 !
81  percent_xflow=prop(index+1)
82  xflow_sum=0
83  do i=2,10
84  if(nint(prop(index+i)).ne.0) then
85  nodemup=kon(ipkon(nint(prop(index+i)))+2)
86  if(v(1,nodemup).gt.0)then
87  xflow_sum=xflow_sum+v(1,nodemup)*iaxial
88  endif
89  endif
90  enddo
91 
92  if(xflow_sum.eq.0.d0) then
93  xflow_sum=1e-5
94  endif
95 !
96  inv=1
97 !
98  pt1=v(2,node1)
99  pt2=v(2,node2)
100  xflow=v(1,nodem)*iaxial
101  tt1=v(0,node1)-physcon(1)
102  tt2=v(0,node2)-physcon(1)
103 !
104  nodef(1)=node1
105  nodef(2)=node1
106  nodef(3)=nodem
107  nodef(4)=node2
108 !
109  idirf(1)=2
110  idirf(2)=0
111  idirf(3)=1
112  idirf(4)=2
113 !
114  if(iflag.eq.2) then
115  numf=4
116 !
117  f=xflow/xflow_sum-percent_xflow
118 !
119  df(1)=0
120  df(2)=0
121  df(3)=1/xflow_sum
122  df(4)=0
123 !
124 ! output
125 !
126  elseif(iflag.eq.3)then
127 !
128  xflow_oil=0
129 !
130  write(1,*) ''
131  write(1,55) ' from node ',node1,
132  & ' to node ', node2,' : air massflow rate = '
133  & ,inv*xflow,
134  & ', oil massflow rate = ',xflow_oil
135  55 format(1x,a,i6,a,i6,a,e11.4,a,a,e11.4,a)
136 !
137  write(1,56)' Inlet node ',node1,' : Tt1 = ',tt1,
138  & ' , Ts1 = ',tt1,' , Pt1 = ',pt1
139 !
140  write(1,*)' Element ',nelem,lakon(nelem)
141  write(1,57)' Massflow upstream = ',xflow_sum,
142  & ' [kg/s]'
143  write(1,58)' Massflow fraction = ', percent_xflow
144  write(1,56)' Outlet node ',node2,': Tt2=',tt2,
145  & ', Ts2=',tt2,', Pt2=',pt2
146 !
147  endif
148  endif
149 !
150  56 format(1x,a,i6,a,e11.4,a,e11.4,a,e11.4,a)
151  57 format(1x,a,e11.4,a)
152  58 format(1x,a,e11.4)
153 !
154  xflow=xflow/iaxial
155  df(3)=df(3)*iaxial
156 !
157  return
subroutine df(x, u, uprime, rpar, nev)
Definition: subspace.f:133
static double * e11
Definition: radflowload.c:42
Hosted by OpenAircraft.com, (Michigan UAV, LLC)