COSMOS v7.655  COSMOSv7655
(AirShowerMC)
chookSkel.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  character*10 append
16  character*100 msg
17  logical ex, opn
18  integer klena
19 
20 ! namelist output
21  call cwriteparam(errorout, 0)
22 ! primary information
23  call cprintprim(errorout)
24 ! observation level information
25  call cprintobs(errorout)
26 
27 
28  call cquhooki(1, mdev) ! get skeleton memo dev #
29  call cquhooki(2, wdev) ! get working disk dev #
30 
31 
32  call cquhooki(3, ngmin) ! get Nh min
33  call cquhooki(4, nhmin) ! get Ng min
34  call cquhooki(5, where) ! where to check
35 
36  call cquhookr(1, sumegmin) ! sum E min
37  call cquhookr(2, sumehmin)
38 
39 
40 
41  call cquhookc(1, msg) ! get file name for sekelton memo
42  call cgetfname(msg, mskel) ! add host name etc if needed
43  call cquhookc(2, msg) ! get file name for working
44  call cgetfname(msg, wskel) ! add host name etc if neededn
45  call cquhookc(3, append) ! append data, if Mskel already exists
46 
47  write(msg, *) 'Skeleton is judged at obs.pos=', Where
48  call cerrormsg(msg, 1)
49  write(msg, *) ' Ngmin=',ngmin,
50  * ' SumEGmin=',sumegmin/1000.,'TeV'
51  call cerrormsg(msg, 1)
52  write(msg, *) ' Nhmin=',nhmin,
53  * ' SumEHmin=',sumehmin/1000.,'TeV'
54  call cerrormsg(msg, 1)
55 
56 !
57  inquire(file=mskel(1:klena(mskel)), opened=opn, exist=ex)
58  if(opn) then
59  call cerrormsg(mskel, 1)
60  call cerrormsg(' already opened: starange', 0)
61  elseif(ex .and.
62  * append(1:klena(append)) .eq. 'append' ) then
63  open(mdev, file=mskel, form='unformatted',status='old')
64  call cerrormsg('skeleton node info. will be appended', 1)
65 ! skip to the end of file
66  do while( .true. )
67  read(mdev, end=100)
68  enddo
69  else
70  if(ex .and.
71  * append(1:klena(append)) .ne. 'append' ) then
72  call cerrormsg(
73  * 'Old skeleton node info. file exists', 1)
74  call cerrormsg(
75  * 'but node info. will NOT be appended', 1)
76  endif
77  open(mdev, file=mskel(1:klena(mskel)), form='unformatted',
78  * status='unknown')
79  endif
80 
81  100 continue
82 
83  open(wdev, file=wskel(1:klena(wskel)), form='unformatted',
84  * status='unknown' )
85 
86  accepted = 0 ! counter; accepted as skeleton
87 !
88  end
89 
90 
91 ! *********************************** hook for Beginning of 1 event
92 ! * All system-level initialization for 1 event generation has been
93 ! * eneded at this moment.
94 ! * After this is executed, event generation starts.
95 ! *
96  subroutine chookbgevent
97  implicit none
98 #include "Zprivate.h"
99 #include "Ztrack.h"
100 #include "Ztrackv.h"
101 #include "Ztrackp.h"
102 #include "Zobs.h"
103 #include "Zobsp.h"
104 #include "Zobsv.h"
105 #include "Zcode.h"
106 
107 
108 
109  type(track)::incident
110  type(coord)::angle
111 
112 
113  integer EventNo
114  integer i, j
115  integer seed(2)
116  real*8 svEasWait, svEthin, kepn
117 
118 
119  np = 0
120  call cqincident( incident, angle)
121  kepn = incident.p.fm.p(4)
122  if( incident.p.code .eq. kgnuc ) then
123  kepn = kepn/ incident.p.subcode
124  endif
125  ethresh = kepn * waitratio
126 
127  sveaswait = easwait ! for safety save
128  svethin = ethin ! //
129  call csetemin(generate2, keminobs2, cutneg, cuteg)
130  easwait = sveaswait ! restore
131  ethin = svethin
132 
133  rewind wdev
134 
135 ! ===================================================
136 
137  eventno = eventno + 1
138  do i = 1, noofsites
139  write(*, 1000)
140  * sngl(obssites(i).pos.depth),
141  * eventno,
142  * incident.p.code,
143  * incident.p.subcode,
144  * incident.p.charge,
145  * incident.p.fm.e,
146  * -angle.r(1),
147  * -angle.r(2),
148  * -angle.r(3)
149  1000 format(f10.3,i9,3i4,e15.5,3(1x,f12.8))
150  enddo
151 
152 
153 ! ===================================================
154  end
155 
156 
157 ! ************************************ hook for observation
158 ! * One particel information is brought here by the system.
159 ! * All information of the particle is in aTrack
160 ! *
161  subroutine chookobs(aTrack, id)
162 !
163 ! Note that every real variable is in double precision so
164 ! that you may output it in sigle precision to save the memory.
165 ! In some cases it is essential to put it in sigle (say,
166 ! for gnuplot).
167 !
168  implicit none
169 #include "Zcode.h"
170 #include "Ztrack.h"
171 #include "Zprivate.h"
172  integer id ! input. 2 ==> reached at an observation level
173 ! 1 ==> aTrack is going out from
174 ! outer boundery.
175 ! 2 ==> reached at an observation level
176 ! 3 ==> reached at inner boundery.
177  type(track)::aTrack
178 !
179 ! For id =2, you need not output the z value, because it is always
180 ! 0 (within the computational accuracy).
181 !
182 ! ===================================================
183  common /testcos/sumg, ng(20), nth, eventno
184  real*8 sumg
185  integer ng, nth, EventNo
186 
187 ! ===================================================
188  if( id .eq. 2 .and. atrack.p.code .ne. kneumu .and.
189  * atrack.p.code .ne. kneue) then
190  np = np + 1
191  if( np .gt. npmax) then
192  call cerrormsg(
193  * '# of particles >NpMax in observation', 0)
194  endif
195  o(np).where = atrack.where
196  o(np).code = atrack.p.code
197  o(np).subcode = atrack.p.subcode
198  o(np).charge = atrack.p.charge
199  o(np).atime = atrack.t
200  o(np).erg = atrack.p.fm.p(4)
201  o(np).mass = atrack.p.mass
202  o(np).x = atrack.pos.xyz.r(1)
203  o(np).y = atrack.pos.xyz.r(2)
204  o(np).wx =atrack.vec.w.r(1)
205  o(np).wy =atrack.vec.w.r(2)
206  o(np).wz =atrack.vec.w.r(3)
207  o(np).zenith = atrack.vec.coszenith
208 ! ===================================================
209  if( o(np).code .le. 6 .and. o(np).code .ne. 3 ) then
210  write(*, 959)
211  * o(np).where,
212  * o(np).code,
213  * o(np).charge,
214  * ( o(np).erg - o(np).mass ),
215  * ( o(np).x ), ( o(np).y ) ,
216  * ( o(np).wx ),
217  * ( o(np).wy ),
218  * ( o(np).wz ),
219  * ( o(np).zenith )
220  959 format(3i3,f12.3,2f16.6,4(1x,f12.8))
221  endif
222 
223 ! ===================================================
224  endif
225 ! ===================================================
226 
227 ! ===================================================
228  end
229 
230 ! *********************************** hook for end of 1 vent
231 ! * At this moment, 1 event generation has been ended.
232 ! *
233  subroutine chookenevent
235  implicit none
236 #include "Ztrack.h"
237 #include "Zcode.h"
238 #include "Ztrackv.h"
239 #include "Zobs.h"
240 #include "Zobsp.h"
241 #include "Zobsv.h"
242 #include "Zmanagerp.h"
243 #include "Zprivate.h"
244 
245  integer i
246  integer seed(2)
247  real sumeg, sumeh
248  logical memorize
249 
250  integer ng, nh
251 
252  ng = 0
253  nh = 0
254 
255  sumeg = 0
256  sumeh =0
257 ! count sum Eg etc
258  do i = 1, np
259  if(o(i).where .eq. where) then
260  if(o(i).code .le. kelec) then
261  ng = ng + 1
262  sumeg = sumeg + o(i).erg
263  elseif( ( o(i).code .ge. kpion .and.
264  * o(i).code .le. knuc )
265  * .or.
266  * o(i).code .eq. kgnuc) then
267  nh = nh + 1
268  sumeh = sumeh + o(i).erg
269  endif
270  endif
271  enddo
272 
273 ! ===================================================
274 
275 
276  memorize =(ng .ge. ngmin .and. sumeg .ge. sumegmin) .or.
277  * ( nh .ge. nhmin .and. sumeh .ge. sumehmin)
278 
279 
280 ! ===================================================
281  if(memorize) then
282  accepted = accepted + 1
283  call cwriteseed ! SeedFile
284 ! flag for end of 1 event on working disk
285  write(wdev) -1, p
286  call cmemorize(wdev, mdev) ! reocord this event
287  endif
288  end
289 
290 
291 ! ********************************* hook for end of a run
292 ! * all events have been created or time lacks
293 ! *
294  subroutine chookenrun
295  implicit none
296 #include "Zprivate.h"
297  character*100 msg
298 ! =========================================
299 
300 
301 ! =========================================
302  call cerrormsg('++++++++++++', 1)
303  write(msg, '(i8, a)') accepted,
304  * ' events are memorized as skeleton'
305  call cerrormsg(msg, 1)
306  call cerrormsg('their seeds are also memorized', 1)
307  call cerrormsg('-----------',1)
308  call cprintstatus ! if don't like, comment out
309  end
310 ! ********************************* hook for trace
311 ! * This is called only when trace > 100
312 ! * User should manage the trace information here.
313 ! * If you use this, you may need some output for trace
314 ! * at the beginning of 1 event generatio and at the end of 1 event
315 ! * generation so that you can identfy each event.
316 ! *
317 ! *
318  subroutine chooktrace
319  implicit none
320 
321 #include "Ztrack.h"
322 #include "Ztrackv.h"
323 #include "Ztrackp.h"
324 #include "Zobs.h"
325 #include "Zobsv.h"
326 
327  real*4 h1, h2
328 !
329 ! Every time a particle is moved in the atmosphere, this routine is called,
330 ! if trace > 100
331 ! For a one track segment,
332 ! TrackBefMove has track information at the beginning of the segment.
333 ! MoveTrack has track information at the end of the segment.
334 !
335 ! You can know the information a track contains in the
336 ! chookObs routine. (Note however, no conversion of coordinate
337 ! has been done. The values are in the Earth xyz system.)
338 ! Besides quantities explained there, you can use, for a given 'track'
339 !
340 ! atrack.pos.xyz.x, atrack.pos.xyz.y, atrack.pos.xyz.z (x,y.z)
341 ! atrack.pos.radiallen (distance from the center of the earth)
342 ! atrack.pos.depth (vertical depth)
343 ! atrack.pos.height (vertical heigth from sea level)
344 !
345 
346  h1 = trackbefmove.pos.height- obssites(noofsites).pos.height
347  h2 = movedtrack.pos.height - obssites(noofsites).pos.height
348 
349  end
350 ! ********************* this is the hook called when
351 ! an electron made an interaction.
352 !
353  subroutine chookeint(never)
354  implicit none
355 #include "Zprivate.h"
356 #include "Ztrack.h"
357 #include "Ztrackv.h"
358 #include "Zcode.h"
359 ! #include "Ztrackp.h"
360  integer never ! input & output
361 !
362 ! If Job = 'newskel', input "never" may be < 0, and MovedTrack
363 ! may not be an electron.
364 ! never = -1 : MovedTrack K.E becomes < KEminObs due to
365 ! energy loss during traveling.
366 ! never = -2 : The same as above, but the particle crosses an
367 ! observation depth, so the calling point to this
368 ! routine is different from the never = -1 case.
369 ! never = -3 : K.E >= KEminiObs. The ptcl is observed at the
370 ! current deepest Obs. level. but at the flesh time
371 ! the deepest level will be more deep so that
372 ! this must be memorized.
373 !
374 ! For above cases, the product is only MovedTrack and
375 ! no particle is in PWork.
376 ! Otherwise,
377 ! MovedTrack is the electron that made interaction
378 ! Pwork contains produced particles (normally gamma, but mayb e).
379 ! Nproduced has the number of particles in Pwork (normally 1)
380 ! IntInfArray(ProcessNo) contains the type of interaction
381 ! IntInfArray(ProcessNo).process will have one of
382 ! 'brems', 'mscat', 'bscat' 'anihi' or 'mbrem'
383 !
384  if(never .lt. 0) then
385  nproduced = 1
386 ! IBM user must modify next
387  pwork(1) = movedtrack.p
388  endif
389 
390 ! high energy parent at node might be used
391 ! for hybrid AS generation if it is in some
392 ! energy region.
393  if( movedtrack.p.code .eq. kelec ) then
394  if( movedtrack.asflag .eq. 0 ) then
395  if( movedtrack.p.fm.p(4) .lt. ethresh ) then
396  movedtrack.asflag = -1
397  endif
398  endif
399  endif
400 
401  call cmemonode(wdev, never)
402 
403  if(movedtrack.asflag .eq. -1) then
404 ! af flesh time, decendent should not be used to
405 ! generae A.S
406  movedtrack.asflag = -2
407  endif
408 
409  10 continue
410  never = 0
411  end
412 
413 ! ********************* this is the hook called when
414 ! a gamma ray made an interaction.
415 !
416  subroutine chookgint(never)
417  implicit none
418 #include "Zprivate.h"
419 #include "Ztrack.h"
420 #include "Ztrackv.h"
421 ! #include "Ztrackp.h"
422 
423  integer never ! input & output
424 
425 ! don't make never = 1, if you want to get
426 ! information after a gamma ray made interaction
427 ! if this is made non zero, this routine will never be called.
428 !
429 ! MovedTrack is the gamma that made interaction
430 ! Pwork contains produced particles.
431 ! Nproduced has the number of particles in Pwork
432 ! IntInfArray(ProcessNo) contains the type of interaction
433 ! IntInfArray(ProcessNo).process will have one of
434 ! 'pair', 'comp', 'photoe' 'photop' 'mpair'
435  call cmemonode(wdev, 1)
436  never = 0
437  end
438 
439 ! ********************* this is the hook called when
440 ! non e-g particle made an interaction.
441 !
442  subroutine chooknepint(never)
443  implicit none
444 #include "Zprivate.h"
445 #include "Ztrack.h"
446 #include "Ztrackv.h"
447 #include "Ztrackp.h"
448 
449  integer never ! input & output
450 
451 
452 ! don't make never = 1, if you want to get
453 ! information after a non-e-g particle made interaction
454 ! if this is made non zero, this routine will never be called.
455 !
456 ! MovedTrack is the particle that made interaction
457 ! Pwork contains produced particles.
458 ! Nproduced has the number of particles in Pwork
459 ! IntInfArray(ProcessNo) contains the type of interaction
460 !
461 
462  call cmemonode(wdev, 1)
463  never = 0
464  end
465 
466  subroutine cmemonode(dev, flag)
467  implicit none
468 #include "Zcode.h"
469 #include "Ztrack.h"
470 #include "Ztrackv.h"
471 #include "Ztrackp.h"
472 #include "Zprivate.h"
473 !
474  integer dev ! input. fortran logical dev. number for data write
475  integer flag ! input. If this is -3, high energy particles must
476  ! also be memorized.
477 !
478 !
479 ! memorize nodal information at interaction.
480 !
481 ! \
482 ! \ projectile = MovedTrack
483 ! /|\ produced particles.: Pwork(i), i=1~Nproduced
484 ! / | \
485 !
486 ! From projectile, the track information is extracted
487 ! From produced particles, only those of K.E< KEminObs is
488 ! extracted and information needed for further tracking is
489 ! obtaned and memorized. The position information is common
490 ! to all the children.
491 !
492 ! Track information; pos =
493 ! xyz = r(1~3), sys
494 ! radiallen, depth, height, colheight
495 ! t
496 ! vec =
497 ! w = r(1~3), sys
498 ! coszenith
499 ! wgt
500 ! where
501 ! asflag
502 !
503 ! Among these, we don't memorize
504 ! sys which is always 'xyz'
505 ! radiallen: can be computed from height
506 ! vec; children knows their direction
507 ! wgt: normally not needed, it should be 1 for skeleton making
508 ! So thinsamping must not be used when making skeleton.
509 ! asflab: should be always F, (assume for skeleton making, hybrid
510 ! air shower is not requested)
511 !
512 ! write track info. of projectile
513 !
514  integer i, nlow
515  real*8 ke
516 
517 !
518 !
519  nlow = 0
520  do i = 1, nproduced
521  ke = pwork(i).fm.p(4) - pwork(i).mass
522  if( ( pwork(i).code .le. kelec .and.
523  * ke .ge. cuteg ) .or.
524  * ( pwork(i).code .gt. kelec .and.
525  * ke .ge. cutneg ) ) then
526 ! count low energy ptcls
527  if(flag .ne. -3) then
528  if( ke .lt. keminobs) then
529  nlow = nlow + 1
530  endif
531  else
532 ! all ptcl must be memorized
533  nlow = nlow + 1
534  endif
535  endif
536  enddo
537 
538  if(nlow .eq. 0 ) return ! *************
539 
540  p.posx = movedtrack.pos.xyz.r(1)
541  p.posy = movedtrack.pos.xyz.r(2)
542  p.posz = movedtrack.pos.xyz.r(3)
543  p.depth = movedtrack.pos.depth
544  p.height = movedtrack.pos.height
545  if( movedtrack.pos.colheight .gt. 1.e36) then
546  p.colheight = 1.e36 ! no col. yet.
547  else
548  p.colheight = movedtrack.pos.colheight
549  endif
550  p.atime = movedtrack.t
551  p.where = movedtrack.where
552  p.coszenith = movedtrack.vec.coszenith
553  p.code = movedtrack.p.code
554  p.erg = movedtrack.p.fm.p(4)
555  p.asflag = movedtrack.asflag
556 
557 
558 !c
559 !c * p.posx, p.posy, p.posz, p.depth, p.height,
560 !c * p.colHeight, p.atime, p.where
561 !c
562 ! write particle info
563 ! p(1~4)
564 ! mass
565 ! code
566 ! subcode
567 ! charge
568 !
569 
570 
571  write(dev) nlow, p
572 
573 
574  do i = 1, nproduced
575  ke = pwork(i).fm.p(4) - pwork(i).mass
576  if( ( pwork(i).code .le. kelec .and.
577  * ke .ge. cuteg ) .or.
578  * ( pwork(i).code .gt. kelec .and.
579  * ke .ge. cutneg ) ) then
580 
581  if(flag .eq. -3 .or. ke .lt. keminobs) then
582  c.code = pwork(i).code
583  c.subcode = pwork(i).subcode
584  c.charge = pwork(i).charge
585  c.fm(1) = pwork(i).fm.p(1)
586  c.fm(2) = pwork(i).fm.p(2)
587  c.fm(3) = pwork(i).fm.p(3)
588  c.fm(4) = pwork(i).fm.p(4)
589  c.mass = pwork(i).mass
590 
591  write(dev) c
592  endif
593  endif
594  enddo
595  end
596 
597 !
598  subroutine cmemorize(from, to)
599  implicit none
600 #include "Ztrack.h"
601 #include "Ztrackv.h"
602 
603  integer from ! working disk file
604  integer to ! permanent disk file where skeleton is sotered
605 
606 
607 
608 
609  integer num, cumnum, irevent(2)
610 ! record /track/ incident
611 ! record /coord/ angle
612 
613 
614  rewind from
615 ! need not memorize, can be generated at flesh
616 ! call cqIncident(incident, angle)
617 
618  call cqeventno(num, cumnum)
619 
620  call cqinirn(irevent) ! seed of the event
621 
622  write(to) cumnum, num, irevent, zfirst
623  call cputhes(to) ! put high energy showers.
624 
625  call cputnodinfo(from, to) ! put nordal info.
626 
627  end
628 
629  subroutine cputhes(to)
630  implicit none
631 #include "Zprivate.h"
632  integer to
633 !
634 ! write high energy sekeleton data into 'to'
635 !
636  integer i
637 
638  write(to) np
639 
640  do i = 1, np
641  write(to) o(i)
642  enddo
643 
644  end
645 
646  subroutine cputnodinfo(from, to)
647  implicit none
648 #include "Zprivate.h"
649 
650  integer from, to
651 
652  integer i, nlow
653 
654  nlow = 1
655  do while ( nlow .ge. 0 )
656  read(from) nlow, p
657  write(to) nlow, p
658  do i = 1, nlow
659  read(from) c
660  write(to) c
661  enddo
662  enddo
663 
664  end
665 
666 
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
subroutine cputhes(to)
Definition: chookSkel.f:630
subroutine cputnodinfo(from, to)
Definition: chookSkel.f:647
dE dx *! Nuc Int sampling table e
Definition: cblkMuInt.h:130
others if is ng
Definition: cblkManager.h:9
*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
max ptcl codes in the kgnuc
Definition: Zcode.h:2
subroutine cmemorize(from, to)
Definition: chookSkel.f:599
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
max ptcl codes in the kelec
Definition: Zcode.h:2
subroutine chooknepint(never)
Definition: chook.f:219
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 cprintprim(out)
Definition: cprintPrim.f:3
max ptcl codes in the kneue
Definition: Zcode.h:2
!onst int nth
Definition: Zprivate.h:2
subroutine chookenevent
Definition: chook.f:116
subroutine cwriteparam(io, force)
Definition: cwriteParam.f:4
subroutine csetemin(gen, eminob, emin, emCas)
Definition: ciniTracking.f:242
********************block data cblkHeavy ********************integer j data *HeavyG2symbol p
Definition: cblkHeavy.h:7
subroutine cwriteseed
Definition: cwriteSeed.f:15
subroutine cquhookr(i, rv)
Definition: cqUHookr.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 *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 cmemonode(dev, flag)
Definition: chookSkel.f:467
********************block data cblkHeavy ********************integer j data *HeavyG2symbol *data *HeavyG2code knuc
Definition: cblkHeavy.h:7
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 cqeventno(num, cumnum)
Definition: cqEventNo.f:3
*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
subroutine cqinirn(ir)
Definition: cwriteSeed.f:4
nodes t
*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 where
Definition: ZavoidUnionMap.h:1
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
max ptcl codes in the kpion
Definition: Zcode.h:2
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
dE dx *! Nuc Int sampling table c
Definition: cblkMuInt.h:130