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