COSMOS v7.655  COSMOSv7655
(AirShowerMC)
chookFlesh.f
Go to the documentation of this file.
1 #include "cmain.f"
2 #include "../chookHybAS.f"
3 #include "ctemplCeren.f"
4 !
5 ! *************************************** hook for Beginning of a Run
6 ! * At this moment, all (system-level) initialization for this run
7 ! * has been ended. After this routine is executed, the system goes into the
8 ! * event creation loop.
9 ! *
10  subroutine chookbgrun
11  implicit none
12 #include "Zmanagerp.h"
13 #include "../../SkelFlesh/Zprivate.h"
14 
15  real*8 temp
16  character*100 msg
17  integer icon
18  integer i
19  eventno = 0
20 
21 ! namelist output
22  call cwriteparam(errorout, 0)
23 ! primary information
24  call cprintprim(errorout)
25 ! observation level information
26  call cprintobs(errorout)
27 
28  call cquhooki(1, mdev) ! get skeleton memo dev #
29  call cquhookc(1, msg) ! get file name for sekelton data
30  call cgetfname(msg, mskel) ! add host name etc if needed
31  call copenfw2(mdev, mskel, 2, icon)
32  if(icon .ne. 1) then
33  call cerrormsg(mskel,1)
34  call cerrormsg(' could not be opened',0)
35  endif
36 
37  call xbgrun
38  end
39 
40 ! *********************************** hook for Beginning of 1 event
41 ! * All system-level initialization for 1 event generation has been
42 ! * eneded at this moment.
43 ! * After this is executed, event generation starts.
44 ! *
45  subroutine chookbgevent
46  implicit none
47 #include "../../SkelFlesh/Zprivate.h"
48 
49 
50  integer nomore
51  call cbegin1ev( nomore )
52  if( nomore .eq. 1) then
53  call cerrormsg('all events have been fleshed', 1)
54  stop !!!!!!!!!!!!
55  endif
56  call cpushinci
57 
58  call xbgevent
59  return
60  end
61  subroutine cbegin1ev(nomore)
62  implicit none
63 #include "../../SkelFlesh/Zprivate.h"
64 #include "Ztrack.h"
65 #include "Ztrackv.h"
66 #include "Ztrackp.h"
67 #include "Zobs.h"
68 #include "Zobsp.h"
69 #include "Zobsv.h"
70 #include "Zcode.h"
71 #include "Zmanager.h"
72 #include "Zmanagerp.h"
73 
74  integer nomore ! output. 0 still there are showers
75  ! 1 no more skeleton showers to be fleshed
76 ! event number, primary
77 
78  type(track):: incident, zsave
79  type(coord):: angle
80 
81  integer i
82  integer seed(2)
83  integer cumnum, num, jeof, fin
84  read( mdev, end=1000, err=999 ) cumnum, num, seedsave, zfirst
85 
86  eventsintherun = eventsintherun + 1
87  eventno = eventno + 1
88 ! reset the seed.
89  call rnd1r(seedsave)
90 ! next incident; confirmed to be the same one as preserved one
91  call cmkincident(incident, fin)
92  if(fin .ne. 0 ) goto 1000
93  zsave = zfirst ! save; this is reset in next
94  call cinitracking( incident )
95 ! set first interaction pos
96  zfirst = zsave
97 !ccc call cresetTimer(Zfirst)
98 
99 
100 
101 ! do your own init for a one event here
102 ! ==========================================================
103 
104 
105 ! ==========================================================
106 !
107 
108  call cgethes(mdev) ! get high energy ptlcs
109  call cobshes ! imitate their observation
110  nomore = 0
111  return
112 
113  1000 continue
114  nomore = 1
115  return
116  999 continue
117  write(0,*) ' Mdev read err'
118  stop 1111
119  end
120 
121 ! ************************************ hook for observation
122 ! * One particel information is brought here by the system.
123 ! * All information of the particle is in aTrack
124 ! *
125  subroutine chookobs(aTrack, id)
126 !
127 ! Note that every real variable is in double precision so
128 ! that you may output it in sigle precision to save the memory.
129 ! In some cases it is essential to put it in sigle (say,
130 ! for gnuplot).
131 !
132  implicit none
133 #include "Zcode.h"
134 #include "Ztrack.h"
135 
136 
137  integer id ! input. 2 ==> reached at an observation level
138 ! 1 ==> aTrack is going out from
139 ! outer boundery.
140 ! 2 ==> reached at an observation level
141 ! 3 ==> reached at inner boundery.
142  type(track):: aTrack
143 !
144 !
145  integer n, i
146  real*8 eps, u
147 
148  if(atrack.wgt .gt. 1.) then
149  n=atrack.wgt
150  eps = atrack.wgt - n
151  call rndc(u)
152  if(u .lt. eps) then
153  n = n + 1
154  endif
155  else
156  n = 1
157  endif
158  do i = 1, n
159  call xobs(atrack, id)
160  enddo
161 
162  end
163 
164 ! *********************************** hook for end of 1 event
165 ! * At this moment, 1 event generation has been ended.
166 ! *
167  subroutine chookenevent
169  implicit none
170 #include "Ztrack.h"
171 #include "Ztrackv.h"
172 #include "Zobs.h"
173 #include "Zobsp.h"
174 #include "Zobsv.h"
175 
176  integer i
177 
178 ! for Job ='newflesh', we must call cfinTracking ourselves.
179  call cfintracking
180 ! end of 1 event; if you need to do some here is
181 ! the place
182 
183  call xenevent
184 
185 
186  end
187 
188 
189 ! ********************************* hook for end of a run
190 ! * all events have been created or time lacks
191 ! *
192  subroutine chookenrun
193  implicit none
194  call cprintstatus ! if don't like, comment out
195  end
196 ! ********************************* hook for trace
197 ! * This is called only when trace > 100
198 ! * User should manage the trace information here.
199 ! * If you use this, you may need some output for trace
200 ! * at the beginning of 1 event generatio and at the end of 1 event
201 ! * generation so that you can identfy each event.
202 ! *
203 ! *
204  subroutine chooktrace
205  implicit none
206 
207 #include "Ztrack.h"
208 #include "Ztrackv.h"
209 #include "Ztrackp.h"
210 #include "Zobs.h"
211 #include "Zobsv.h"
212 
213  real*4 h1, h2
214 !
215 ! Every time a particle is moved in the atmosphere, this routine is called,
216 ! if trace > 100
217 ! For a one track segment,
218 ! TrackBefMove has track information at the beginning of the segment.
219 ! MoveTrack has track information at the end of the segment.
220 !
221 ! You can know the information a track contains in the
222 ! chookObs routine. (Note however, no conversion of coordinate
223 ! has been done. The values are in the Earth xyz system.)
224 ! Besides quantities explained there, you can use, for a given 'track'
225 !
226 ! atrack.pos.xyz.x, atrack.pos.xyz.y, atrack.pos.xyz.z (x,y.z)
227 ! atrack.pos.radiallen (distance from the center of the earth)
228 ! atrack.pos.depth (vertical depth)
229 ! atrack.pos.height (vertical heigth from sea level)
230 !
231 
232  h1 = trackbefmove.pos.height- obssites(noofsites).pos.height
233  h2 = movedtrack.pos.height - obssites(noofsites).pos.height
234 
235  end
236 ! ********************* this is the hook called when
237 ! an electron made an interaction.
238 !
239  subroutine chookeint(never)
240  implicit none
241  integer never ! input & output
242  never = 1
243  end
244 
245 ! ********************* this is the hook called when
246 ! a gamma ray made an interaction.
247 !
248  subroutine chookgint(never)
249  implicit none
250  integer never ! input & output
251  never = 1
252  end
253 
254 ! ********************* this is the hook called when
255 ! non e-g particle made an interaction.
256 !
257  subroutine chooknepint(never)
258  implicit none
259  integer never ! input & output
260  never = 1
261  end
262 
263 
264  subroutine cgethes(from)
265  implicit none
266 #include "../../SkelFlesh/Zprivate.h"
267  integer from
268 
269  integer i
270 
271  read(from) np
272  do i = 1, np
273  read(from) o(i)
274  enddo
275  end
276 
277  subroutine cobshes
278  implicit none
279 #include "../../SkelFlesh/Zprivate.h"
280 #include "Ztrack.h"
281 !
282 ! memorized high energy showers at the skeleton making
283 ! time is put into the chookObs as if they are really observed
284  type(track):: aTrack
285 
286  integer i
287  logical HEobs ! if T, currently observing
288  common /zheobs/ heobs ! particles those obsrved at skeelton making time
289 
290  heobs = .true.
291  do i = 1, np
292  atrack.where = o(i).where
293  atrack.p.code = o(i).code
294  atrack.p.subcode = o(i).subcode
295  atrack.p.charge = o(i).charge
296  atrack.t = o(i).atime
297  atrack.p.fm.p(4) = o(i).erg
298  atrack.p.mass = o(i).mass
299  atrack.pos.xyz.r(1) = o(i).x
300  atrack.pos.xyz.r(2) = o(i).y
301  atrack.vec.w.r(1) = o(i).wx
302  atrack.vec.w.r(2) = o(i).wy
303  atrack.vec.w.r(3) = o(i).wz
304  atrack.vec.coszenith = o(i).zenith
305  call chookobs(atrack, 2)
306  enddo
307  heobs = .false.
308  end
309 
310 
311 ! push all low energy partilces in the skeleton in the stack
312 
313  subroutine cpushinci
314  implicit none
315 #include "../../SkelFlesh/Zprivate.h"
316 #include "Ztrack.h"
317 #include "Ztrackv.h"
318  integer i
319 
320  type(track)::aTrack
321 
322  call cinitstack ! empty the stack
323 
324  read(mdev) nooflowe
325  do i = 1, nooflowe
326  read(mdev) atrack
327 ! aTrack is already complete track so push it directly.
328  call cpush(atrack)
329  enddo
330 ! sort stack dscendent order
331  call csortstack
332 
333  end
subroutine cgetfname(fnin, fn)
Definition: copenf.f:275
subroutine cerrormsg(msg, needrtn)
Definition: cerrorMsg.f:4
subroutine cprintobs(io)
Definition: cprintObs.f:2
subroutine chookgint(never)
Definition: chook.f:191
subroutine cpushinci
Definition: chookFlesh.f:314
*Zfirst p fm *Zfirst p Zfirst p Zfirst p *Zfirst p *Zfirst pos xyz r
Definition: ZavoidUnionMap.h:1
latitude latitude this system is used *****************************************************************! type coord sequence union map real z z in m endmap xyz map real ! latitude in deg is to the north ! longitude in deg is to the east *h ! height in m endmap llh map real ! polar angle ! azimuthal angle *radius ! radial distance endmap sph endunion character *sys ! which system xyz
Definition: Zcoord.h:25
subroutine chookenrun
Definition: chook.f:147
Definition: Ztrack.h:44
subroutine chooknepint(never)
Definition: chook.f:219
subroutine cbegin1ev(nomore)
Definition: chookFlesh.f:67
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 cinitstack
Definition: cstack.f:76
subroutine cprintprim(out)
Definition: cprintPrim.f:3
subroutine cfintracking
Definition: cfinTracking.f:2
subroutine cmkincident(incident, fin)
Definition: cmkIncident.f:5
subroutine chookenevent
Definition: chook.f:116
subroutine rndc(u)
Definition: rnd.f:91
subroutine cwriteparam(io, force)
Definition: cwriteParam.f:4
********************block data cblkHeavy ********************integer j data *HeavyG2symbol p
Definition: cblkHeavy.h:7
subroutine copenfw2(io, fnin, form, icon)
Definition: copenf.f:205
*Zfirst p fm *Zfirst p Zfirst p code
Definition: ZavoidUnionMap.h:1
subroutine chooktrace
Definition: chook.f:275
subroutine cgethes(from)
Definition: chookFlesh.f:322
struct ob o[NpMax]
Definition: Zprivate.h:34
*Zfirst p fm *Zfirst p Zfirst p Zfirst p *Zfirst p *Zfirst pos xyz Zfirst pos xyz *Zfirst pos *Zfirst pos Zfirst pos height
Definition: ZavoidUnionMap.h:1
latitude latitude this system is used *****************************************************************! type coord sequence union map real y
Definition: Zcoord.h:25
subroutine cinitracking(incident)
Definition: ciniTracking.f:2
*Zfirst p fm *Zfirst p Zfirst p Zfirst p *Zfirst p *Zfirst pos xyz Zfirst pos xyz *Zfirst pos *Zfirst pos Zfirst pos *Zfirst pos *Zfirst Zfirst vec w *Zfirst vec w Zfirst vec coszenith
Definition: ZavoidUnionMap.h:1
float erg[maxp]
Definition: Zprivate.h:7
*Zfirst p fm *Zfirst p Zfirst p Zfirst p *Zfirst p charge
Definition: ZavoidUnionMap.h:1
subroutine chookobs(aTrack, id)
Definition: chook.f:59
subroutine chookbgevent
Definition: chook.f:39
subroutine cquhooki(i, iv)
Definition: cqUHookr.f:15
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
subroutine csortstack
Definition: cstack.f:102
nodes t
subroutine cquhookc(i, cv)
Definition: cqUHookr.f:28
Definition: Zcoord.h:43
*Zfirst p fm *Zfirst p mass
Definition: ZavoidUnionMap.h:1
subroutine cprintstatus
Definition: cendRun.f:35
subroutine chookbgrun
Definition: chook.f:15
subroutine cpush(a)
Definition: cstack.f:4
subroutine xbgrun
Definition: interface.f:10
*Zfirst p fm *Zfirst p Zfirst p Zfirst p *Zfirst p *Zfirst pos xyz Zfirst pos xyz *Zfirst pos *Zfirst pos Zfirst pos *Zfirst pos *Zfirst Zfirst vec w *Zfirst vec w Zfirst vec *Zfirst wgt
Definition: ZavoidUnionMap.h:1
subroutine chookeint(never)
Definition: chook.f:162
! 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
*Zfirst p fm *Zfirst p Zfirst p Zfirst p subcode
Definition: ZavoidUnionMap.h:1
subroutine cobshes
Definition: chookFlesh.f:335