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

Go to the source code of this file.

Functions/Subroutines

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

Function/Subroutine Documentation

◆ ts_calc()

subroutine ts_calc ( real*8, intent(in)  xflow,
real*8, intent(in)  Tt,
real*8, intent(in)  Pt,
real*8, intent(in)  kappa,
real*8, intent(in)  r,
real*8, intent(in)  A,
real*8, intent(inout)  Ts,
integer, intent(in)  icase 
)
20 !
21 ! Calculation of the static temperature Ts from the total
22 ! temperature Tt, the total pressure Pt and the mass flow xflow
23 !
24 ! this subroutine solves the implicit equation
25 ! xflow*dsqrt(Tt)/(A*Pt)-C*(Tt/Ts)**expon*(Tt/Ts-1)**0.5d0=0.d0
26 !
27 ! expon=-0.5d0*(kappa+1.d0)/(kappa-1.d0)
28 ! C=dsqrt(2.d0/r*kappa/(kappa-1.d0))
29 !
30 ! xflow is the mass flow for 360 degrees, sign is irrelevant
31 ! for this routine
32 !
33 ! author: Yannick Muller
34 !
35  implicit none
36 !
37  integer icase,i
38 !
39  real*8 xflow,tt,pt,ts,kappa,r,f,df,a,expon,ts_old,c,ttzts,
40  & deltats,ttzts_crit, qred_crit,qred,h1,h2,h3
41 !
42  intent(in) xflow,tt,pt,kappa,r,a,icase
43 !
44  intent(inout) ts
45 !
46  expon=-0.5d0*(kappa+1.d0)/(kappa-1.d0)
47 !
48  c=dsqrt(2.d0/r*kappa/(kappa-1.d0))
49 !
50 ! f=xflow*dsqrt(Tt)/(A*Pt)-C*(Tt/Ts)**expon*(Tt/Ts-1)**0.5d0
51 !
52 ! df=-C*Tt/Ts**expon*(expon/Ts*(Tt/Ts-1)**0.5d0
53 ! & -0.5d0*Tt/Ts/Ts*(Tt/Ts-1.d0)**(-0.5d0))
54 !
55  ts_old=tt
56 !
57  if(dabs(xflow).le.1e-9) then
58  ts=tt
59  return
60  endif
61 !
62  qred=abs(xflow)*dsqrt(tt)/(a*pt)
63 !
64 ! optimised estimate of T static
65 !
66  ts=tt/(1+(qred**2/c**2))
67 !
68 ! adiabatic
69 !
70  if(icase.eq.0) then
71  ttzts_crit=(kappa+1.d0)/2.d0
72 !
73 ! isothermal
74 !
75  else
76  ttzts_crit=(1d0+(kappa-1.d0)/(2.d0*kappa))
77  endif
78 !
79  qred_crit=c*(ttzts_crit)**expon*(ttzts_crit-1.d0)**0.5d0
80 !
81  if(qred.ge.qred_crit) then
82  ts=tt/ttzts_crit
83  return
84  endif
85 !
86  i=0
87 !
88 ! start of the Newton-Raphson-Procedure to solve the nonlinear
89 ! equation
90 !
91  do
92  i=i+1
93  ttzts=tt/ts
94  h1=ttzts-1.d0
95  h2=dsqrt(h1)
96  h3=ttzts**expon
97 !
98  f=c*h2*h3
99 !
100  df=f*(expon+0.5d0*ttzts/h1)/ts
101 !
102  f=qred-f
103  deltats=-f/df
104 !
105  ts=ts+deltats
106 !
107  if((((dabs(ts-ts_old)/ts_old).le.1.e-8))
108  & .or.((dabs(ts-ts_old)).le.1.e-10)) then
109  exit
110  else if(i.gt.20) then
111  ts=0.9*tt
112  exit
113  endif
114  ts_old=ts
115  enddo
116 !
117  return
subroutine df(x, u, uprime, rpar, nev)
Definition: subspace.f:133
Hosted by OpenAircraft.com, (Michigan UAV, LLC)