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

Go to the source code of this file.

Functions/Subroutines

subroutine catmoscnst1
 
subroutine catmoscnst2
 

Function/Subroutine Documentation

◆ catmoscnst1()

subroutine catmoscnst1 ( )

Definition at line 2 of file catmosCnst1.f.

References c, d, d0, fd1, and parameter().

Referenced by __atmosd.f__(), __atmosd2.f__(), __atmosd3.f__(), __atmosd4.f__(), and ciniatmos().

2 
3 ! catmosCnst1: compute some basic constants
4 ! a, b, d0, cum, H
5 ! where
6 ! a, b; const shown below
7 ! d0: see the formula below.
8 ! cumd: amount of atmospher above the node.
9 ! H: scale height at the node
10 ! all in mks unit.
11 !
12 ! The scale height is approximated by a
13 ! number of stright lines as a function of height. The data in
14 ! stdatmos1.d gives height, temperatur, etc at each nodal point.
15 !
16 ! The scale height, H, is expressed by H = H0 + a(z-z0)
17 ! = kT/mg
18 ! in each region.
19 ! We neglect height dependence of gravitational accelleration g,
20 ! and the average mass of air molecules, m.
21 ! Since the data table gives T(z)= T0 + b(z-z0) at the nodal points,
22 ! we can first get b,
23 ! and then a by a = dH/dz = k/mg * b. H at a nodal point, z,
24 ! is obtained as H(z) =kT(z)/mg.
25 !
26 ! The density off a nodal point is given by
27 !
28 ! rho = rho0 * (1+ a(z-z0)/H(z0))**(-1-1/a) (a != 0)
29 ! = rho0 * exp(- (z-z0)/H) (a =0; hence H is const)
30 !
31 ! (We employ H(z0) as the scale height in the segment)
32 ! The amount of air between given heights, z1 and z2 is by
33 !
34 ! d = d0 *(fd(z1) - fd(z2)) where
35 !
36 ! fd(z) = (1+ a(z-z0)/H(z0))**(-1/a) (a != 0)
37 !
38 ! = exp(-(z-z0)/H ) (a = 0)
39 !
40 ! where d0 = rho0*H(z0)
41 !
42 ! If z1=z0, d becomes
43 !
44 ! d= d0 ( 1 - fd(z2))
45 !
46 !
47 
48 !
49 !
50  implicit none
51 
52 #include "Zglobalc.h"
53 #include "Zatmos.h"
54 
55  integer i
56 
57 
58 
59  real*8 m, g, k, nuc, c2, kbymg, dtdz
60 ! m: avrage molecule mass number of air
61 ! g: gravitational acceleration.
62 ! k: Boltzman's const.
63 ! c2: c**2
64  parameter( m=14.5d0 *2, nuc =938.3d6, g=9.80665d0, c2=c*c,
65  * k= 8.617d-5, kbymg = k*c2/(m*nuc*g))
66 
67  integer nodes
68 
69 #include "Zstdatmosf.h"
70 
71  nodes = atmos%nodes
72  do i = 1, nodes-1
73  dtdz = (atmos%T(i+1) - atmos%T(i))/ (atmos%z(i+1) - atmos%z(i))
74  atmos%H(i) = kbymg * atmos%T(i) ! at atmos%z(i)
75  atmos%b(i) = dtdz
76  atmos%a(i) = kbymg * atmos%b(i)
77  atmos%d0(i) = atmos%rho(i) * atmos%H(i)
78 !///////////
79 ! write(0, *)
80 ! * atmos.z(i), atmos.d0(i), atmos.a(i), atmos.b(i), atmos.H(i)
81 !////////
82  enddo
83 
84  atmos%cumd(nodes) = atmos%rho(nodes) * hinf ! put very small amout
85 
86  do i = nodes-1, 1, -1
87  if(atmos%a(i) .eq. 0.) then
88 ! exponential atmosphere
89  atmos%cumd(i) = atmos%cumd(i+1) +
90  * atmos%d0(i)*
91 ! * (1.0- exp(-(atmos.z(i+1) - atmos.z(i))/atmos.H(i)))
92  * (1.0- fd0(atmos%z(i+1), atmos%z(i), atmos%H(i)))
93  else
94  atmos%cumd(i) = atmos%cumd(i+1) +
95  * atmos%d0(i)*
96  * (1.0 -
97  * fd1(atmos%z(i+1), atmos%a(i),atmos%z(i), atmos%H(i) ) )
98 
99 ! * (1.0+atmos.a(i)*(atmos.z(i+1) -atmos.z(i))/atmos.H(i))
100 ! * **(-1./atmos.a(i))
101 ! * )
102  endif
103  enddo
104  do i = 1, nodes
105 ! atmos.logP(i) = log(atmos.P(i))
106  atmos%logrho(i) = log(atmos%rho(i))
107  atmos%logcumd(i) = log(atmos%cumd(i))
108  atmos%logcumdi(atmos%nodes-i+1) = atmos%logcumd(i)
109  atmos%logrhoi(atmos%nodes-i+1) = atmos%logrho(i)
110  atmos%zi(atmos%nodes-i+1) = atmos%z(i)
111  enddo
integer npitbl real *nx parameter(n=101, npitbl=46, nx=n-1) real *8 uconst
nodes i
block data cblkEvhnp ! currently usable models data RegMdls ad *special data *Cekaon d0
Definition: cblkEvhnp.h:5
dE dx *! Nuc Int sampling table d
Definition: cblkMuInt.h:130
block data include Zlatfit h c fitting region data data data data data d0 data data d0 data data m
Definition: ZlatfitBD.h:35
dE dx *! Nuc Int sampling table g
Definition: cblkMuInt.h:130
! to be included just before the execution code ! density as a function of height real * fd1
Definition: Zstdatmosf.h:3
dE dx *! Nuc Int sampling table c
Definition: cblkMuInt.h:130
Here is the call graph for this function:
Here is the caller graph for this function:

◆ catmoscnst2()

subroutine catmoscnst2 ( )

Definition at line 115 of file catmosCnst1.f.

References kcsplcoef().

Referenced by __atmosd2.f__(), __atmosd3.f__(), __atmosd4.f__(), and ciniatmos().

115 ! compute c-spline coef. for later use
116  implicit none
117 #include "Zatmos.h"
118 
119  integer nodes
120 
121  nodes = atmos%nodes
122 ! height--> rho
123  call kcsplcoef(atmos%z, atmos%logrho, nodes, atmos%coefh2r,
124  * maxnodes)
125 ! height--> log(depth)
126  call kcsplcoef(atmos%z, atmos%logcumd, nodes, atmos%coefh2d,
127  * maxnodes)
128 ! h--> log(P)
129 ! call kcsplCoef(atmos.z, atmos.logP, nodes, atmos.coefh2P,
130 ! * maxnodes)
131 ! h--> scale H
132  call kcsplcoef(atmos%z, atmos%H, nodes-1, atmos%coefh2H,
133  * maxnodes)
134 ! h--> T
135  call kcsplcoef(atmos%z, atmos%T, nodes, atmos%coefh2T,
136  * maxnodes)
137 ! log(rho)--> h
138 ! call kcsplCoef(atmos.logrho, atmos.z, nodes, atmos.coefr2h,
139 ! * maxnodes)
140 ! log(depth)--> h
141  call kcsplcoef(atmos%logcumdi, atmos%zi, nodes, atmos%coefd2h,
142  * maxnodes)
143 
144 ! log(detph) --> log(rho)
145  call kcsplcoef(atmos%logcumdi, atmos%logrhoi, nodes,
146  * atmos%coefd2r, maxnodes)
147 ! log(P) --> h
148 ! call kcsplCoef(atmos.logP, atmos.z, nodes, atmos.coefP2h,
149 ! * maxnodes)
150 
subroutine kcsplcoef(x, y, n, coef, nc)
Definition: kcsplCoef.f:2
Here is the call graph for this function:
Here is the caller graph for this function: