COSMOS v7.655  COSMOSv7655
(AirShowerMC)
reanal.f
Go to the documentation of this file.
1  implicit none
2 #include "Zprivate.h"
3 #include "Ztrack.h"
4 
5 
6  integer num, i, nr, firstev, lastev
7  type(track):: Zfirst
8  character*160 str
9  character*80 skelin, stro(3)
10 
11  real*8 temp
12 
13  integer klena, nlow
14 
15  integer cumnum, irevent(2)
16  common /cselev/ cumnum, irevent
17 
18  rdev = 70
19 
20  str = ' '
21  read(*,'(a)') str
22 !
23  do i = 1, 3
24  stro(i) = ' '
25  enddo
26 !
27  call ksplit(str, 80, 3, stro, nr)
28  if(nr .lt. 3) stop 'must give 1 files and 2 event #'
29 
30  skelin = stro(1)
31  read( stro(2), * ) firstev
32  read( stro(3), * ) lastev
33 
34  open(rdev, file=skelin(1:klena(skelin)), form='unformatted',
35  * status='old')
36 
37  do while(.true.)
38  read(rdev, end= 1000) cumnum, num, irevent, zfirst
39 
40  if(cumnum .gt. lastev) goto 1000
41 
42  call cgethes(rdev)
43  nlow = 1
44  do while (nlow .ne. -1)
45  read(rdev) nlow, p
46  do i = 1, nlow
47  read(rdev) c
48  enddo
49  enddo
50  enddo
51  1000 continue
52  if(lastev .gt. cumnum) then
53  call cerrormsg(
54  * 'EOF before the specified last event reached; ok?', 1)
55  else
56  call cerrormsg(
57  * ' all events processed successfully', 1)
58  endif
59 
60  end
61 
62 
63  subroutine cgethes(from)
64  implicit none
65 #include "Zprivate.h"
66  integer from
67 
68  integer cumnum, irevent(2)
69  common /cselev/ cumnum, irevent
70 
71 !
72  integer i
73  read(from) np
74  do i = 1, np
75  read(from) o(i)
76  enddo
77 ! ==========================================
78 ! do your analysis and if you want to
79 ! memorize the it for flesh, make
80 ! Copy = .true. Then, the event number
81 ! (or seed and event number) is written on stdout.
82 ! You can use
83 ! o(i).xxx where xxx is
84 !
85 ! where: obsv. location
86 ! code, subcode, charge: code=9 is Nucleus, A=subcode,Z=charge
87 ! atime: arrival time ns
88 ! erg: energy GeV
89 ! mass: mass in GeV
90 ! x,y: in m
91 ! wx,wy,wx: direction cos.
92 ! zenith: cos of zenith angle.
93 !
94 
95  copy = .true.
96 
97 
98 ! ==========================================
99  if(copy) then
100  write(*,*) cumnum
101 ! or
102 ! write(*,*) irevent, cumnum
103  endif
104  end
subroutine ksplit(a, m, n, b, nr)
Definition: ksplit.f:2
subroutine cerrormsg(msg, needrtn)
Definition: cerrorMsg.f:4
nodes i
Definition: Ztrack.h:44
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
********************block data cblkHeavy ********************integer j data *HeavyG2symbol p
Definition: cblkHeavy.h:7
subroutine cgethes(from)
Definition: chookFlesh.f:322
struct ob o[NpMax]
Definition: Zprivate.h:34
integer function klena(cha)
Definition: klena.f:20
dE dx *! Nuc Int sampling table c
Definition: cblkMuInt.h:130