COSMOS v7.655  COSMOSv7655
(AirShowerMC)
bin2ascii.f
Go to the documentation of this file.
1 ! **********************
2  module modbin2ascii
3  integer,parameter:: ndepth=50
4  integer fn1
5  real*8 asdep(ndepth), munit(ndepth)
6  real*8 esize0(ndepth),
8  * seloss0(ndepth),
10  * cog0
11  integer evno0
12 
13  contains
14  subroutine get1hyb( rew )
15  implicit none
16  logical rew
17  character*256 input0
18  integer i, klena
19 
20 
21 ! data format
22 ! 1 xxx
23 ! 1 xxx
24 ! 1 xxx
25 ! 1 xxx
26 !
27 ! 2 xxx
28 ! 2 xxx
29 ! 2 xxx
30 ! 2 xxx
31 !
32 ! 3 xxx
33 ! 3 xxx
34 ! 3 xxx
35 ! 3 xxx
36  if(rew) rewind fn1
37  input0 = "x"
38  do while (input0(1:10) .ne. " ")
39  input0=" "
40  read( fn1 ,'(a)') input0
41  if(input0(1:10) .ne. " ") then
42  read(input0(1:klena(input0)), *)
43  * evno0, i, asdep(i), esize0(i), age0(i),
44  * cogdep0(i), seloss0(i),
45  * munit(i), ng0(i), ne0(i), nmu0(i), cog0
46  endif
47  enddo
48  end subroutine get1hyb
49 ! ***********************
50  subroutine mergehyb1(h1)
51  use modhistogram1
52  use modhistogram2
53  use modhistogram3
54  implicit none
55  type(histogram1) h1
56  type(histogram2) h2
57  type(histogram3) h3
58 
59 
60  integer klena
61  integer j
62 
63  do while (h1%c%eventno .ne. evno0)
64  call get1hyb( h1%c%eventno .lt. evno0)
65  enddo
66 !
67 ! available variables
68 ! * idx, ASdep(idx), Esize0(idx), age0(idx),
69 ! * SEloss0(idx), munit(idx),
70 ! * Ng0(idx), Ne0(idx), Nmu0(idx),
71 ! * ASdep(idx)/cog0, cog0
72 !
73 ! this part must be consistent with
74 ! FleshHist/interface.f output for evid
75  read(h1%c%id, '(i3)') j
76  write(h1%c%id,
77  * '(i3, i5, f5.2, f5.2,
78  * i5, i4)')
79  * j, int( asdep(j) ),
80  * age0(j), asdep(j)/cog0,
81  * int(munit(j)), int(cog0)
82 
83  return
84 ! *******************
85  entry mergehyb2(h2)
86 ! *******************
87 
88 
89  do while (h2%c%eventno .ne. evno0)
90  call get1hyb( h2%c%eventno .lt. evno0)
91  enddo
92  read(h2%c%id, '(i3)') j
93 
94  write(h2%c%id,
95  * '(i3, i5, f5.2, f5.2,
96  * i5, i4)')
97  * j, int( asdep(j) ),
98  * age0(j), asdep(j)/cog0,
99  * int(munit(j)), int(cog0)
100 
101  return
102 ! *****************
103  entry mergehyb3(h3)
104 ! ****************
105 
106 
107  do while (h3%c%eventno .ne. evno0)
108  call get1hyb( h3%c%eventno .lt. evno0)
109  enddo
110 
111  read(h3%c%id, '(i3)') j
112  write(h3%c%id,
113  * '(i3, i5, f5.2, f5.2,
114  * i5, i4)')
115  * j, int( asdep(j) ),
116  * age0(j), asdep(j)/cog0,
117  * int(munit(j)), int(cog0)
118 
119 
120  end subroutine mergehyb1
121 
122  subroutine openhyb(icon)
123  implicit none
124  integer icon ! output. 1--> hybrid must be read
125  ! 0--> hybrid need not be used
126  integer leng
127 
128  character*256 hyb0
129  integer kgetenv2
130 
131 
132  fn1= 3
133  leng = kgetenv2("HYBFILE0", hyb0)
134  call copenfw2(fn1, hyb0, 1, icon)
135  if(icon .ne. 1) then
136  write(0,*)
137  * '*************** caution ************'
138  write(0,*)
139  * "You haven't given env. var. HYBFILE0"
140  write(0,*)
141  * "or File specified by HYBFILE0"
142  if( icon .eq. 0) then
143  write(0,*) 'not exists'
144  else
145  write(0,*) ' cannot be opened '
146  endif
147  write(0,*)
148  * "It's ok if you don't merge hybrid data file"
149  icon = 0
150  else
151  write(0,*) hyb0(1:leng), ' opened'
152  icon = 1
153  endif
154 
155  evno0 =0
156 
157  end subroutine openhyb
158 
159  end module modbin2ascii
160 
161  program main
162  use modhistogram1
163  use modhistogram2
164  use modhistogram3
165 
166  use modbin2ascii
167 
168  implicit none
169 ! binary hist file: binary output file made by programs
170 ! in UserHook/Hist
171 ! Generic routine for converting a binary hist file
172 ! to an ascii hist file.
173 ! Usage: compile: make -f bin2ascii.mk
174 ! execution:
175 ! set environmental variable HISTFILE0 to be
176 ! a binary hist file.
177 ! Optionally you may give another env. var.
178 ! HYBFILE0 to be a hybrid data output file
179 ! made by a program in UserHook/DisPara/FleshHist
180 ! or UserHook/ForTA. If this is not given,
181 ! no special treatment of mergining hybrid result.
182 ! (If the event is generated by using DisPara,
183 ! you need to merge hybrid data, otherwise
184 ! the hybrid data in the binary hist is not
185 ! correct. )
186 ! bin2ascii$ARCH > xxxx.hist
187 !
188  type(histogram1) h10
189  type(histogram2) h20
190  type(histogram3) h30
191 
192 
193  integer kgetenv2
194  integer leng, i, fn0
195  integer icon0, iconhyb
196  character*256 hist0
197  character*6 histid0, oldhist
198  real normf
199  data normf/-1.0/
200 
201  fn0 = 2
202  leng = kgetenv2("HISTFILE0", hist0)
203  call copenfw2(fn0, hist0, 2, icon0)
204  if(icon0 .ne. 1) then
205  write(0,*) "File specified by HISTFILE0 "
206  write(0,*) hist0
207  if( icon0 .eq. 0) then
208  write(0,*) 'not exists'
209  else
210  write(0,*) ' cannot be opened '
211  endif
212  write(0,*) ' icon=',icon0
213  stop 9999
214  else
215  write(0,*) hist0(1:leng), ' opened'
216  endif
217  leng = kgetenv2("OLDHIST", oldhist)
218  if( oldhist .eq. "yes") then
219  write(0,*) ' old histogram format is assumed'
220  call kwhistfmt(.true.) ! for old fomat
221  else
222  call kwhistfmt(.false.)
223  endif
224  call openhyb(iconhyb)
225  do while(.true.)
226  read( fn0, end=1000 ) histid0
227 !/////////////
228 ! write(0,*) ' id=', histid0
229 !/////////////////////
230  100 continue
231  if( histid0 .eq. '#hist1' ) then
232 !//////////////
233 ! write(0,*) ' reading hist'
234 !////////////
235  call kwhistr(h10, fn0, icon0)
236 ! if(iconhyb .eq. 1) then
237  if(icon0 .eq. 1) then
238  call mergehyb1(h10)
239  endif
240  call kwhists(h10, normf)
241  call kwhistpr(h10, -6) ! to stdout
242  call kwhistd(h10)
243  elseif(histid0 .eq. '#hist2' ) then
244  call kwhistr2(h20, fn0, icon0)
245  if(icon0 .eq. 1) then
246  call mergehyb2(h20)
247  endif
248  call kwhists2(h20, normf)
249  call kwhistpr2(h20, -6) ! to stdout
250  call kwhistd2(h20)
251  elseif(histid0 .eq. '#hist3' ) then
252  call kwhistr3(h30, fn0, icon0)
253  if(icon0 .eq. 1) then
254  call mergehyb3(h30)
255  endif
256  call kwhists3(h30, normf)
257  call kwhistpr3(h30, -6)
258  call kwhistd3(h30)
259  else
260  write(0,*) 'histid=', histid0, ' invalid'
261  stop 9000
262  endif
263  enddo
264  1000 continue
265  write(0,*) 'all events processed '
266  end program main
267 
subroutine get1hyb
real *8, dimension(ndepth) nmu0
Definition: bin2ascii.f:6
real *8, dimension(ndepth) esize0
Definition: bin2ascii.f:6
void kwhistd(struct histogram1 *h)
real *8, dimension(ndepth) ng0
Definition: bin2ascii.f:6
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 openhyb(icon)
Definition: bin2ascii.f:123
subroutine copenfw2(io, fnin, form, icon)
Definition: copenf.f:205
real *8 cog0
Definition: bin2ascii.f:6
subroutine mergehyb1(h1)
Definition: bin2ascii.f:51
void kwhistr(struct histogram1 *h, FILE *bfnor, int icon)
integer fn1
Definition: bin2ascii.f:4
real *8, dimension(ndepth) cogdep0
Definition: bin2ascii.f:6
program main
Definition: ascii2bin.f:1
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 ! knockon is considered Obsolete *PhotoProd false
Definition: cblkElemag.h:7
real *8, dimension(ndepth) ne0
Definition: bin2ascii.f:6
void kwhists(struct histogram1 *h, float inorm)
void kwhistpr(struct histogram1 *h, FILE *fno)
real *8, dimension(ndepth) munit
Definition: bin2ascii.f:5
real *8, dimension(ndepth) seloss0
Definition: bin2ascii.f:6
real *8, dimension(ndepth) age0
Definition: bin2ascii.f:6
integer evno0
Definition: bin2ascii.f:11
real *8, dimension(ndepth) asdep
Definition: bin2ascii.f:5
integer, parameter ndepth
Definition: bin2ascii.f:3