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

Go to the source code of this file.

Functions/Subroutines

subroutine csampaf0 (iowk, filen, sampInfo)
 
subroutine csampaf (sampInfo, xs)
 
subroutine csampafintp (sampInfo, xv, ans)
 
subroutine csampafmax (sampInfo, xmax, fmax)
 

Function/Subroutine Documentation

◆ csampaf()

subroutine csampaf ( type(sampaf)  sampInfo,
real*8  xs 
)

Definition at line 33 of file csampAF.f.

References ksampaf(), n, and x.

33  implicit none
34 #include "ZsampAF.h"
35  type(sampaf):: sampinfo
36  real*8 xs
37  call ksampaf(sampinfo.x, sampinfo.yi, sampinfo.n,
38  * sampinfo.coef2, sampinfo.n, xs)
subroutine ksampaf(x, yi, n, coef2, nc, xs)
Definition: ksampAF.f:120
! constants thru Cosmos real ! if multiplied to deg radian Torad ! light velocity m sec ! infinty ! kg m2 *Togpcm2 g cm2 ! g cm2 *Tokgpm2 kg m2 ! cm *Tom m ! m *Tocm cm ! g cm3 *Tokgpm3 kg m3 ! kg m3 *Togpcm3 g cm3 ! sec *Tonsec nsec ! Tesla m ! Avogadro *A2deninv ! mfp *n * xs
Definition: Zglobalc.h:18
integer n
Definition: Zcinippxc.h:1
! 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:

◆ csampaf0()

subroutine csampaf0 ( integer  iowk,
character*(*)  filen,
type(sampaf)  sampInfo 
)

Definition at line 2 of file csampAF.f.

References copenf(), cskipcomment(), ksampaf0(), n, true, x, and y.

2  implicit none
3 #include "ZsampAF.h"
4  integer iowk ! input file logical number temporarily used
5  character*(*) filen ! input. file name which contains (x,dn/dx)
6 
7  type(sampaf):: sampinfo
8 
9  integer icon
10  integer i
11  call copenf(iowk, filen, icon)
12  if(icon .ne. 0 ) then
13  write(0,*) ' error '
14  write(0,*) 'file: ',filen
15  write(0,*) ' could not be opened'
16  stop
17  endif
18  i = 0
19  call cskipcomment(iowk, icon)
20  if(icon .ne. 0 ) stop
21  do while ( .true. )
22  read(iowk, *, end=100) sampinfo.x(i+1), sampinfo.y(i+1)
23  i = i + 1
24  enddo
25  100 continue
26  close(iowk)
27  sampinfo.n = i
28  call ksampaf0(sampinfo.x, sampinfo.y, sampinfo.n,
29  * sampinfo.coef, sampinfo.n, sampinfo.yi,
30  * sampinfo.sum, sampinfo.coef2)
nodes i
block data cblkElemag data *AnihiE ! Eposi< 1 TeV, anihilation considered *X0/365.667/, ! radiation length of air in kg/m2 *Ecrit/81.e-3/, ! critical energy of air in GeV *MaxComptonE/1./, ! compton is considered below 1 GeV *MaxPhotoE/1.e-3/, ! above this, PhotoElectric effect neg. *MinPhotoProdE/153.e-3/, ! below 153 MeV, no gp --> hadrons ! scattering const not MeV *Knockon true
Definition: cblkElemag.h:7
subroutine ksampaf0(x, y, n, coef, nc, yi, total, coef2)
Definition: ksampAF.f:55
subroutine copenf(io, fnin, icon)
Definition: copenf.f:8
latitude latitude this system is used *****************************************************************! type coord sequence union map real y
Definition: Zcoord.h:25
subroutine cskipcomment(io, icon)
Definition: cskipComment.f:19
integer n
Definition: Zcinippxc.h:1
! 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:

◆ csampafintp()

subroutine csampafintp ( type(sampaf)  sampInfo,
real*8  xv,
real*8  ans 
)

Definition at line 41 of file csampAF.f.

References kcsplintp(), n, x, and y.

Referenced by csampafmax().

41 ! this is not for sampling but simply
42 ! get value of y at xf
43 
44  implicit none
45 #include "ZsampAF.h"
46  type(sampaf):: sampinfo ! input obtained by csampAF0
47  real*8 xv ! input
48  real*8 ans ! output y at xv
49  call kcsplintp(sampinfo.x, sampinfo.y, sampinfo.n,
50  * sampinfo.coef, sampinfo.n, xv, ans)
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
latitude latitude this system is used *****************************************************************! type coord sequence union map real y
Definition: Zcoord.h:25
subroutine kcsplintp(x, y, n, coef, nc, v, f)
Definition: kcsplIntp.f:2
integer n
Definition: Zcinippxc.h:1
! 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:

◆ csampafmax()

subroutine csampafmax ( type(sampaf)  sampInfo,
real*8  xmax,
real*8  fmax 
)

Definition at line 53 of file csampAF.f.

References csampafintp(), and n.

53 ! find max position and value of given function
54 ! (approx value)
55  implicit none
56 #include "ZsampAF.h"
57  type(sampaf):: sampinfo ! input obtained by csampAF0
58  real*8 xmax ! ouput. max position in (x1,x2) ; approx value
59  real*8 fmax ! outpu. max function value
60 
61  real*8 x1, x2
62  real*8 x, dx, temp
63  integer i
64 
65  x1 = sampinfo.x(1)
66  x2 = sampinfo.x(sampinfo.n)
67  x = x1
68  dx = (x2-x1)/sampinfo.n/10.
69  xmax = x1
70  call csampafintp(sampinfo, xmax, fmax)
71  call csampafintp(sampinfo, x2, temp)
72  if( fmax .lt. temp ) then
73  xmax = x2
74  fmax = temp
75  endif
76  x = x + dx
77  do while (x .lt. x2-dx/2 )
78  call csampafintp(sampinfo, x, temp)
79  if( fmax .lt. temp ) then
80  xmax = x
81  fmax = temp
82  endif
83  x = x + dx
84  enddo
subroutine csampafintp(sampInfo, xv, ans)
Definition: csampAF.f:41
block data include Zlatfit h c fitting region data x1(1)/0.03/
nodes i
integer npitbl real *nx dx real dx
Definition: Zcinippxc.h:10
real(4), dimension(:), allocatable, save temp
Definition: cNRLAtmos.f:29
block data include Zlatfit h c fitting region data x2(1)/0.5/data x1(2)/0.3/
integer n
Definition: Zcinippxc.h:1
! 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: