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 "Zprivate.h"
14 
15  real*8 temp
16  character*100 msg
17  integer klena
18  integer i
19 ! ==================================================
20 
21  integer seed(2)
22 ! ==================================================
23 
24  eventno = 0
25  realbegin = .true.
26  topofnode = .true.
27 
28 ! namelist output
29  call cwriteparam(errorout, 0)
30 ! primary information
31  call cprintprim(errorout)
32 ! observation level information
33  call cprintobs(errorout)
34 
35  call cquhooki(1, mdev) ! get skeleton memo dev #
36  call cquhooki(6, howflesh) ! how to flesh. non zero-->orthodox
37  call cquhookc(1, msg) ! get file name for sekelton memo
38  call cgetfname(msg, mskel) ! add host name etc if needed
39 
40 ! do i = 1, NoOfSites
41 ! write(*, *) ' depth ',
42 ! * sngl(ObsSites(i).pos.depth/10.)
43 ! enddo
44 
45  open(mdev, file=mskel(1:klena(mskel)), form='unformatted',
46  * status='old' )
47 
48  end
49 
50 ! *********************************** hook for Beginning of 1 event
51 ! * All system-level initialization for 1 event generation has been
52 ! * eneded at this moment.
53 ! * After this is executed, event generation starts.
54 ! *
55  subroutine chookbgevent
56  implicit none
57 #include "Zprivate.h"
58 
59 
60  integer nomore
61  if( realbegin ) then
62  call cbegin1ev( nomore )
63  if( nomore .eq. 1) then
64  call cerrormsg('all events are fleshed', 1)
65  stop !!!!!!!!!!!!
66  endif
67  topofnode = .true.
68  endif
69  call c1by1
70 
71 
72  end
73  subroutine cbegin1ev(nomore)
74  implicit none
75 #include "Zprivate.h"
76 #include "Ztrack.h"
77 #include "Ztrackv.h"
78 #include "Ztrackp.h"
79 #include "Zobs.h"
80 #include "Zobsp.h"
81 #include "Zobsv.h"
82 #include "Zcode.h"
83 #include "Zmanager.h"
84 #include "Zmanagerp.h"
85 
86  integer nomore ! output. 0 still there are showers
87  ! 1 no more skeleton showers to be fleshed
88 ! event number, primary
89 
90  type(track):: incident, zsave
91  type(coord):: angle
92 
93  integer i
94  integer seed(2)
95  integer cumnum, num, jeof, fin
96  read( mdev, end=1000 ) cumnum, num, seedsave, zfirst
97 
98  eventsintherun = eventsintherun + 1
99  eventno = eventno + 1
100 ! get random seed at skelton making; this can work
101 ! if seed file is supplied
102 ! call creadSeed(SeedSave, EventNo, jeof)
103 ! if( jeof .ne. 0 ) goto 1000
104 
105 ! reset the seed.
106  call rnd1r(seedsave)
107 
108 ! next incident; confirmed to be the same one as preserved one
109  call cmkincident(incident, fin)
110 
111  if(fin .ne. 0 ) goto 1000
112  zsave = zfirst ! save; this is reset in next
113  call cinitracking( incident )
114 ! set first interaction pos
115  zfirst = zsave
116  call cresettimer(zfirst)
117 
118  realbegin = .false.
119 
120 ! do your own init for a one event here
121 ! ==========================================================
122  call cqincident( incident, angle)
123  write(*,999)
124  * eventno,
125  * incident.p.code,
126  * incident.p.subcode,
127  * incident.p.charge,
128  * incident.p.fm.e,
129  * -angle.r(1),
130  * -angle.r(2),
131  * -angle.r(3)
132  999 format(i9,3i4,e15.5,3(1x,f12.8))
133 
134 
135 ! ==========================================================
136 !
137 
138  call cgethes(mdev) ! get high energy ptlcs
139  call cobshes ! imitate their observation
140  nomore = 0
141  return
142 
143  1000 continue
144  nomore = 1
145  end
146 
147 ! ************************************ hook for observation
148 ! * One particel information is brought here by the system.
149 ! * All information of the particle is in aTrack
150 ! *
151  subroutine chookobs(aTrack, id)
152 !
153 ! Note that every real variable is in double precision so
154 ! that you may output it in sigle precision to save the memory.
155 ! In some cases it is essential to put it in sigle (say,
156 ! for gnuplot).
157 !
158  implicit none
159 #include "Zcode.h"
160 #include "Ztrack.h"
161 #include "Zprivate.h"
162  integer id ! input. 2 ==> reached at an observation level
163 ! 1 ==> aTrack is going out from
164 ! outer boundery.
165 ! 2 ==> reached at an observation level
166 ! 3 ==> reached at inner boundery.
167  type(track):: aTrack
168 !
169 ! For id =2, you need not output the z value, because it is always
170 ! 0 (within the computational accuracy).
171 !
172  if(id .eq. 2 .and. atrack.p.code .ne. kneumu .and.
173  * atrack.p.code .ne. kneue) then
174 
175 ! ===================================================
176 
177  if( atrack.p.code .le. 6 .and. atrack.p.code .ne. 3 ) then
178 ! write(*, 959)
179 ! * aTrack.where,
180 ! * aTrack.p.code,
181 ! * aTrack.p.charge,
182 ! * sngl( aTrack.p.fm.p(4) - aTrack.p.mass ),
183 ! * sngl( aTrack.pos.xyz.r(1) ),
184 ! * sngl( aTrack.pos.xyz.r(2) ) ,
185 ! * sngl( aTrack.vec.w.r(1) ),
186 ! * sngl( aTrack.vec.w.r(2) ),
187 ! * sngl( aTrack.vec.w.r(3) ),
188 ! * sngl( aTrack.vec.coszenith )
189 ! 959 format(3i3,f12.3,2f16.6,4(1x,f12.8))
190  endif
191 
192 ! ===================================================
193 
194 ! write(*,'(4i5, g15.4,g15.3)')
195 ! * aTrack.where, aTrack.p.code, aTrack.p.subcode,
196 ! * aTrack.p.charge, sngl( aTrack.t ),
197 ! * sngl( aTrack.p.fm.p(4) - aTrack.p.mass)
198 ! * sngl( aTrack.pos.xyz.r(1) ), sngl( aTrack.pos.xyz.r(2) ),
199 ! * sngl( aTrack.vec.w.r(1) ), sngl(aTrack.vec.w.r(2) ),
200 ! * sngl(aTrack.vec.w.r(3) ),
201 ! * sngl(aTrack.vec.coszenith)
202 
203  endif
204  end
205 
206 ! *********************************** hook for end of 1 event
207 ! * At this moment, 1 event generation has been ended.
208 ! *
209  subroutine chookenevent
211  implicit none
212 #include "Zprivate.h"
213 #include "Ztrack.h"
214 #include "Ztrackv.h"
215 #include "Zobs.h"
216 #include "Zobsp.h"
217 #include "Zobsv.h"
218 
219  integer i
220 
221  if(realend) then
222  call cfintracking
223 ! real end of 1 event; if you need to do some here is
224 ! the place
225 ! ========================================================
226 
227  if(observeas) then
228 ! electron size in B approx.
229  do i = 1, noofassites
230  write(*, *) asobssites(i).age, asobssites(i).esize
231  enddo
232  endif
233 
234 
235 
236 ! ========================================================
237 !
238  else
239 ! there is still low energy skeleton ptcls
240 ! nothing to do here
241  endif
242 
243  end
244 
245 
246 ! ********************************* hook for end of a run
247 ! * all events have been created or time lacks
248 ! *
249  subroutine chookenrun
250  implicit none
251 #include "Zprivate.h"
252 ! =========================================================
253 
254 ! =========================================================
255  call cprintstatus ! if don't like, comment out
256  end
257 ! ********************************* hook for trace
258 ! * This is called only when trace > 100
259 ! * User should manage the trace information here.
260 ! * If you use this, you may need some output for trace
261 ! * at the beginning of 1 event generatio and at the end of 1 event
262 ! * generation so that you can identfy each event.
263 ! *
264 ! *
265  subroutine chooktrace
266  implicit none
267 
268 #include "Ztrack.h"
269 #include "Ztrackv.h"
270 #include "Ztrackp.h"
271 #include "Zobs.h"
272 #include "Zobsv.h"
273 
274  real*4 h1, h2
275 !
276 ! Every time a particle is moved in the atmosphere, this routine is called,
277 ! if trace > 100
278 ! For a one track segment,
279 ! TrackBefMove has track information at the beginning of the segment.
280 ! MoveTrack has track information at the end of the segment.
281 !
282 ! You can know the information a track contains in the
283 ! chookObs routine. (Note however, no conversion of coordinate
284 ! has been done. The values are in the Earth xyz system.)
285 ! Besides quantities explained there, you can use, for a given 'track'
286 !
287 ! atrack.pos.xyz.x, atrack.pos.xyz.y, atrack.pos.xyz.z (x,y.z)
288 ! atrack.pos.radiallen (distance from the center of the earth)
289 ! atrack.pos.depth (vertical depth)
290 ! atrack.pos.height (vertical heigth from sea level)
291 !
292 
293  h1 = trackbefmove.pos.height- obssites(noofsites).pos.height
294  h2 = movedtrack.pos.height - obssites(noofsites).pos.height
295 
296  end
297 ! ********************* this is the hook called when
298 ! an electron made an interaction.
299 !
300  subroutine chookeint(never)
301  implicit none
302  integer never ! input & output
303  never = 1
304  end
305 
306 ! ********************* this is the hook called when
307 ! a gamma ray made an interaction.
308 !
309  subroutine chookgint(never)
310  implicit none
311  integer never ! input & output
312  never = 1
313  end
314 
315 ! ********************* this is the hook called when
316 ! non e-g particle made an interaction.
317 !
318  subroutine chooknepint(never)
319  implicit none
320  integer never ! input & output
321  never = 1
322  end
323 
324 
325  subroutine cgethes(from)
326  implicit none
327 #include "Zprivate.h"
328  integer from
329 
330  integer i
331 
332  read(from) np
333  do i = 1, np
334  read(from) o(i)
335  enddo
336  end
337 
338  subroutine cobshes
339  implicit none
340 #include "Zprivate.h"
341 #include "Ztrack.h"
342 !
343 ! memorized high energy showers at the skeleton making
344 ! time is put into the chookObs as if they are really observed
345  type(track):: aTrack
346 
347  integer i
348 
349  do i = 1, np
350  atrack.where = o(i).where
351  atrack.p.code = o(i).code
352  atrack.p.subcode = o(i).subcode
353  atrack.p.charge = o(i).charge
354  atrack.t = o(i).atime
355  atrack.p.fm.p(4) = o(i).erg
356  atrack.p.mass = o(i).mass
357  atrack.pos.xyz.r(1) = o(i).x
358  atrack.pos.xyz.r(2) = o(i).y
359  atrack.vec.w.r(1) = o(i).wx
360  atrack.vec.w.r(2) = o(i).wy
361  atrack.vec.w.r(3) = o(i).wz
362  atrack.vec.coszenith = o(i).zenith
363  call chookobs(atrack, 2)
364  enddo
365  end
366 
367 
368 ! process low energy partilces in the skeleton 1 by 1
369 
370  subroutine c1by1
371  implicit none
372 #include "Zprivate.h"
373 #include "Ztrack.h"
374 #include "Ztrackv.h"
375 
376  character*100 msg
377 
378  call cinitstack ! empty the stack
379 
380  if( topofnode ) then
381  read(mdev) nooflowe, p
382  if( p.asflag .eq. -1 .and. observeas
383  * .and. howflesh .ne. 0 ) then
384  call embedas
385  endif
386  nlowcounter = 0
387  if( nooflowe .eq. -1 ) then
388  realend = .true.
389  realbegin = .true.
390  return ! ************
391  endif
392  endif
393 
394  realbegin = .false.
395  realend = .false.
396 
397 
398  if( nlowcounter .eq. nooflowe ) then
399  topofnode =.true.
400  return ! **************
401  endif
402 
403  topofnode = .false.
404 ! still not the end of 1 event
405 
406  read(mdev) c
407 
408  nlowcounter = nlowcounter + 1
409  call cmove_c_stack ! move c into stack
410 
411  end
412 !
413  subroutine embedas
414  implicit none
415 #include "Zprivate.h"
416 #include "Ztrack.h"
417 #include "Zearth.h"
418 
419  type(track):: el
420 
421  el.pos.depth = p.depth
422  el.vec.coszenith = p.coszenith
423  el.pos.radiallen = p.height + eradius
424  el.pos.height = p.height
425  el.p.fm.p(4) = p.erg
426  el.wgt = 1.0
427  call cobas(el)
428  end
429 
430 
431 
432 
433  subroutine cmove_c_stack
434  implicit none
435 
436 #include "Zprivate.h"
437 #include "Ztrack.h"
438 #include "Zearth.h"
439 
440  type(track):: aTrack
441 !
442 ! a child of the current parent is moved to stack
443 ! as a track info.
444 !
445  atrack.pos.xyz.r(1) = p.posx
446  atrack.pos.xyz.r(2) = p.posy
447  atrack.pos.xyz.r(3) = p.posz
448  atrack.pos.depth = p.depth
449  atrack.pos.height = p.height
450  atrack.pos.colheight = p.colheight
451  atrack.t = p.atime
452 
453  atrack.where = p.where
454 
455  atrack.p.code = c.code
456  atrack.p.subcode = c.subcode
457  atrack.p.charge = c.charge
458  atrack.p.fm.p(1) = c.fm(1)
459  atrack.p.fm.p(2) = c.fm(2)
460  atrack.p.fm.p(3) = c.fm(3)
461  atrack.p.fm.p(4) = c.fm(4)
462  atrack.p.mass = c.mass
463 
464 ! --------------- next must be compute here
465 
466  atrack.pos.radiallen = eradius +atrack.pos.height
467  atrack.pos.xyz.sys = 'xyz'
468  atrack.vec.w.sys = 'xyz'
469  atrack.wgt = 1.0
470  if(p.asflag .ne. 0 .and. howflesh .ne. 0 ) then
471 ! asflag= -1 or -2 may happen
472  atrack.asflag = 1
473  else
474 ! parent is made to be unconceivable, so
475 ! don't put flag for AS generation
476  atrack.asflag = 0
477  endif
478 
479 
480 
481 
482  call cresetdirec( atrack ) ! set vec.w and coszenith
483 
484  call cpush(atrack)
485  end
subroutine cgetfname(fnin, fn)
Definition: copenf.f:275
*Zfirst p fm *Zfirst p Zfirst p Zfirst p *Zfirst p *Zfirst pos xyz Zfirst pos xyz *Zfirst pos *Zfirst pos depth
Definition: ZavoidUnionMap.h:1
subroutine cerrormsg(msg, needrtn)
Definition: cerrorMsg.f:4
subroutine cprintobs(io)
Definition: cprintObs.f:2
subroutine cqincident(incident, AngleAtObs)
Definition: cmkIncident.f:486
subroutine chookgint(never)
Definition: chook.f:191
dE dx *! Nuc Int sampling table e
Definition: cblkMuInt.h:130
*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 colheight
Definition: ZavoidUnionMap.h:1
*Zfirst p fm *Zfirst p Zfirst p Zfirst p *Zfirst p *Zfirst pos xyz r
Definition: ZavoidUnionMap.h:1
subroutine embedas
Definition: chookFlesh.f:409
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
*Zfirst p fm *Zfirst p Zfirst p Zfirst p *Zfirst p *Zfirst pos xyz Zfirst pos xyz sys
Definition: ZavoidUnionMap.h:1
subroutine cfintracking
Definition: cfinTracking.f:2
max ptcl codes in the kneue
Definition: Zcode.h:2
subroutine cmkincident(incident, fin)
Definition: cmkIncident.f:5
subroutine chookenevent
Definition: chook.f:116
subroutine cwriteparam(io, force)
Definition: cwriteParam.f:4
subroutine cresetdirec(aTrack)
Definition: cresetDirec.f:5
********************block data cblkHeavy ********************integer j data *HeavyG2symbol p
Definition: cblkHeavy.h:7
*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 Zfirst Zfirst asflag
Definition: ZavoidUnionMap.h:1
*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
max ptcl codes in the kneumu
Definition: Zcode.h:2
*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 cresettimer(aTrack)
Definition: ctracking.f:302
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
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
*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 c1by1
Definition: chookFlesh.f:367
subroutine cmove_c_stack
Definition: chookFlesh.f:429
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 cobas(el)
Definition: cobAS.f:9
*Zfirst p fm *Zfirst p Zfirst p Zfirst p *Zfirst p *Zfirst pos xyz Zfirst pos xyz *Zfirst pos radiallen
Definition: ZavoidUnionMap.h:1
subroutine cobshes
Definition: chookFlesh.f:335
dE dx *! Nuc Int sampling table c
Definition: cblkMuInt.h:130