COSMOS v7.655  COSMOSv7655
(AirShowerMC)
cl2tTblA.f File Reference

Go to the source code of this file.

Functions/Subroutines

real *8 function clen2thickta (z, leng)
 
subroutine cl2tintp (x, y, xx, ans)
 
real *8 function ct2lta (z, t)
 

Function/Subroutine Documentation

◆ cl2tintp()

subroutine cl2tintp ( real*8, dimension(*)  x,
real*8, dimension(*)  y,
real*8  xx,
real*8  ans 
)

Definition at line 34 of file cl2tTblA.f.

References kdwhereis(), kpolintp(), and parameter().

Referenced by clen2thickta(), and ct2lta().

34  implicit none
35 #include "Zatmos.h"
36 
37  real*8 x(*) ! input.
38  real*8 y(*) ! input.
39  real*8 xx ! input. some value inside x(*).
40  real*8 ans ! output. interpolated value of y at xx.
41 
42  real*8 error
43  integer m
44  integer loca, k
45  parameter(m = 3) ! old value was 5
46 ! where is xx in x
47  call kdwhereis(xx, numstep, x, 1, loca)
48  if(loca .ge. numstep) then
49 ! give arbitray large value.
50  ans = y(numstep)
51  else
52 ! use max of m points around loca for interpolation
53  k = min(max(loca - (m-1)/2,1), numstep+1-m)
54  ! max of m points from k
55  call kpolintp(x(k), 1, y(k), 1, m, xx, ans, error)
56  endif
integer npitbl real *nx parameter(n=101, npitbl=46, nx=n-1) real *8 uconst
atmos%rho(atmos%nodes) **exp(-(z-atmos%z(atmos%nodes))/Hinf) elseif(z .lt. atmos%z(1)) then ans=atmos%rho(1) **exp((atmos%z(1) -z)/atmos%H(1)) else call kdwhereis(z, atmos%nodes, atmos%z, 1, i) a=atmos%a(i) if(a .ne. 0.d0) then ans=atmos%rho(i) **(1+a *(z-atmos%z(i))/atmos%H(i)) **(-1.0d0-1.d0/a) else ans=*atmos%rho(i) *exp(-(z-atmos%z(i))/atmos%H(i)) endif endif ! zsave=z ! endif cvh2den=ans end ! ---------------------------------- real *8 function cvh2temp(z) implicit none ! vettical height to temperatur(Kelvin) real *8 z ! input. vertical height in m ! output is temperature of the atmospher in Kelvin real *8 ans integer i if(z .gt. atmos%z(atmos%nodes)) then ans=atmos%T(atmos%nodes) elseif(z .lt. atmos%z(1)) then ans=atmos%T(1)+atmos%b(1) *(z - atmos%z(1)) else call kdwhereis(z, atmos%nodes, atmos%z, 1, i) ans=atmos%T(i)+atmos%b(i) *(z-atmos%z(i)) endif cvh2temp=ans end !--------------------------------------------- real *8 function cthick2h(t) implicit none real *8 t ! input. air thickness in kg/m^2 real *8 logt, ans integer i real *8 dod0, fd, a logt=log(t) if(t .ge. atmos%cumd(1)) then ans=atmos%z(1) - *(logt - atmos%logcumd(1)) *atmos%H(1) elseif(t .le. atmos%cumd(atmos%nodes)) then ans=atmos%z(atmos%nodes) - *Hinf *log(t/atmos%cumd(atmos%nodes)) else call kdwhereis(t, atmos%nodes, atmos%cumd, 1, i) ! i is such that X(i) > x >=x(i+1) ans
subroutine kpolintp(xa, xstep, ya, ystep, n, x, y, error)
Definition: kpolintp.f:163
subroutine kdwhereis(x, in, a, step, loc)
Definition: kdwhereis.f:27
latitude latitude this system is used *****************************************************************! type coord sequence union map real y
Definition: Zcoord.h:25
block data include Zlatfit h c fitting region data data data data data d0 data data d0 data data m
Definition: ZlatfitBD.h:35
! structure defining a particle at production ! Basic idea of what is to be contained in ! the particle structue is that dynamical ones should be included those derivable from the particle code ! is not included ******************************************************type fmom momentum sequence union map real e endmap map real * x
Definition: Zptcl.h:21
Here is the call graph for this function:
Here is the caller graph for this function:

◆ clen2thickta()

real*8 function clen2thickta ( real*8  z,
real*8  leng 
)

Definition at line 5 of file cl2tTblA.f.

References cl2tintp().

5 !
6 ! z: real*8. input. vertical height in m.
7 ! leng: real*8. input. length along cosz direction.
8 ! function value. thickness of air in kg/m2. for leng
9  implicit none
10 #include "Zatmos.h"
11  real*8 z, leng
12  real*8 s1, t1, s2, t2
13 
14 ! get slant length along cosz for z: s1
15 ! get slant length ; s2 = s1 + leng
16 ! get thickness <-- by T vs slant length
17  call cl2tintp(heighttbl, lentbl, z, s1)
18  call cl2tintp(lentbl, thicktbl, s1, t1)
19  s2 = s1 + leng
20  call cl2tintp(lentbl, thicktbl, s2, t2)
21  clen2thickta = t2 - t1
22  if(clen2thickta .le. 0.) then
23 ! actually case of < 0 dose not happen.
24 ! when particle comes to the last segment of
25 ! the table, this happens. and the
26 ! particle does not move at all. To avoid this
27 ! give some thikness. (actually this is not
28 ! used but this is for comfort only)
29  clen2thickta = 10.
30  endif
nodes z
integer leng
Definition: interface2.h:1
subroutine cl2tintp(x, y, xx, ans)
Definition: cl2tTblA.f:34
real *8 function clen2thickta(z, leng)
Definition: cl2tTblA.f:5
Here is the call graph for this function:

◆ ct2lta()

real*8 function ct2lta ( real*8  z,
real*8  t 
)

Definition at line 61 of file cl2tTblA.f.

References cl2tintp().

61  implicit none
62 #include "Zatmos.h"
63 
64  real*8 z ! input vertical height in m
65  real*8 t ! input. thickness of air in kg/m^2 along cosz
66 !
67  real*8 t1, t2, s1, s2
68 ! z ---> thick
69  call cl2tintp(heighttbl, thicktbl, z, t1)
70  call cl2tintp(heighttbl, lentbl, z, s1)
71  t2 = t1 + t
72  call cl2tintp(thicktbl, lentbl, t2, s2)
73  ct2lta = s2- s1
74  if(ct2lta .le. 0.) then
75 ! actually case of < 0 dose not happen.
76 ! when a particle comes to an end segment of
77 ! a table, this happnes and the particle
78 ! cannot move. To avoid this, give some
79 ! finite value to move it.
80 ! clen2thickTA a thickness corresponding to this length
81 ! but is not used since the particle reaches
82 ! the bottom observation depth.
83  ct2lta = 10.
84  endif
nodes z
subroutine cl2tintp(x, y, xx, ans)
Definition: cl2tTblA.f:34
nodes t
real *8 function ct2lta(z, t)
Definition: cl2tTblA.f:61
Here is the call graph for this function: