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

Go to the source code of this file.

Functions/Subroutines

subroutine tt_calc (xflow, Tt, Pt, kappa, r, A, Ts, icase)
 

Function/Subroutine Documentation

◆ tt_calc()

subroutine tt_calc ( real*8  xflow,
real*8  Tt,
real*8  Pt,
real*8  kappa,
real*8  r,
real*8  A,
real*8  Ts,
integer  icase 
)
20 !
21 ! this subroutine solves the implicit equation
22 ! f=xflow*dsqrt(Tt)/(a*Pt)-C*(TtdT)**expon*(Ttdt-1)**0.5d0
23 ! in order to find Tt when Ts , xflow, pt and a are given
24 !
25 ! author: Yannick Muller
26 !
27  implicit none
28 !
29  integer icase,i
30 !
31  real*8 xflow,tt,pt,ts,kappa,r,f,df,a,expon,tt_old,c,ttzts,
32  & deltatt,ttzts_crit,qred,h1,h2,h3
33 !
34  expon=-0.5d0*(kappa+1.d0)/(kappa-1.d0)
35 !
36  c=dsqrt(2.d0/r*kappa/(kappa-1.d0))
37 !
38 ! f=xflow*dsqrt(Tt)/(A*Pt)-C*(Tt/Ts)**expon*(Tt/Ts-1)**0.5d0
39 !
40 ! df=-C*Ttdt**expon*(expon/Ts*(Tt/Ts-1)**0.5d0
41 ! & -0.5d0*Tt/Ts/Ts*(Tt/Ts-1.d0)**(-0.5d0))
42 !
43  if(dabs(xflow).le.1e-9) then
44  tt=ts
45  return
46  endif
47 !
48 ! initial guess
49 !
50  qred=abs(xflow)*dsqrt(ts)/(a*pt)
51  tt=ts*(1+(qred**2/c**2))
52 !
53 ! adiabatic
54 !
55  if(icase.eq.0) then
56 !
57  ttzts_crit=(kappa+1.d0)/2.d0
58 !
59 ! isothermal
60 !
61  else
62 !
63  ttzts_crit=(1d0+(kappa-1.d0)/(2.d0*kappa))
64 !
65  endif
66 !
67  if(tt/ts.gt.ttzts_crit) then
68  tt=ts*(ttzts_crit+1.d0)/2.d0
69  endif
70 !
71  i=0
72  tt_old=tt
73 !
74  do
75  i=i+1
76  ttzts=tt/ts
77  h1=ttzts-1.d0
78  h2=dsqrt(h1)
79  h3=ttzts**expon
80 !
81  f=c*h2*h3
82 !
83  df=0.5*dabs(xflow)/(a*pt*dsqrt(tt))
84  & -c*h2*h3*(expon/tt+1.d0/(2.d0*h1*ts))
85 !
86  qred=abs(xflow)*dsqrt(tt)/(a*pt)
87 !
88  f=qred-f
89  deltatt=-f/df
90 !
91  tt=tt+deltatt
92 !
93  if((((dabs(tt-tt_old)/tt_old).le.1.e-8))
94  & .or.((dabs(tt-tt_old)).le.1.e-10)
95  & .or.((dabs(f).le.1e-5).and.(deltatt.lt.1e-3))) then
96 !
97  if(tt/ts.gt.ttzts_crit) then
98  tt=ts*ttzts_crit
99  endif
100  exit
101  else if((i.gt.40)) then
102  tt=0.99*ts*ttzts_crit
103  exit
104  endif
105  tt_old=tt
106  enddo
107 !
108  return
subroutine df(x, u, uprime, rpar, nev)
Definition: subspace.f:133
Hosted by OpenAircraft.com, (Michigan UAV, LLC)