COSMOS v7.655  COSMOSv7655
(AirShowerMC)
chookFlesh.f File Reference

Go to the source code of this file.

Functions/Subroutines

subroutine chookbgrun
 
subroutine chookbgevent
 
subroutine cbegin1ev (nomore)
 
subroutine chookobs (aTrack, id)
 
subroutine chookenevent
 
subroutine chookenrun
 
subroutine chooktrace
 
subroutine chookeint (never)
 
subroutine chookgint (never)
 
subroutine chooknepint (never)
 
subroutine cgethes (from)
 
subroutine cobshes
 
subroutine cpushinci
 

Function/Subroutine Documentation

◆ cbegin1ev()

subroutine cbegin1ev ( integer  nomore)

Definition at line 82 of file chookFlesh.f.

References cgethes(), cinitracking(), cmkincident(), and cobshes().

82  implicit none
83 #include "../../SkelFlesh/Zprivate.h"
84 #include "Ztrack.h"
85 #include "Ztrackv.h"
86 #include "Ztrackp.h"
87 #include "Zobs.h"
88 #include "Zobsp.h"
89 #include "Zobsv.h"
90 #include "Zcode.h"
91 #include "Zmanager.h"
92 #include "Zmanagerp.h"
93  integer nomore ! output. 0 still there are showers
94  ! 1 no more skeleton showers to be fleshed
95 ! event number, primary
96 
97  type(track):: incident, zsave
98  type(coord):: angle
99 
100  integer i
101  integer seed(2)
102  integer cumnum, num, jeof, fin
103  read( mdev, end=1000, err=999 ) cumnum, num, seedsave,
104 #if defined (KEKB) || (KEKA)
105 #include "ZavoidUnionMap.h"
106 !////////////
107 ! write(0,*) ' rank=',mpirank,' ir=',SeedSave,' Zfirst=',
108 ! * Zfirst.pos.depth*0.1, 'g/cm2'
109 !//////////////
110 #else
111  * zfirst
112 #endif
113  eventsintherun = eventsintherun + 1
114  eventno = eventno + 1
115 ! reset the seed.
116  call rnd1r(seedsave)
117 ! next incident; confirmed to be the same one as preserved one
118  call cmkincident(incident, fin)
119  if(fin .ne. 0 ) goto 1000
120  zsave = zfirst ! save; this is reset in next
121  call cinitracking( incident )
122 ! set first interaction pos
123  zfirst = zsave
124 ! call cresetTimer(Zfirst)
125 
126 
127 
128 ! do your own init for a one event here
129 ! ==========================================================
130 
131 
132 ! ==========================================================
133 !
134 
135  call cgethes(mdev) ! get high energy ptlcs
136  call cobshes ! imitate their observation
137  nomore = 0
138  return
139 
140  1000 continue
141  nomore = 1
142  return
143  999 continue
144  write(0,*) ' Mdev read err'
145  stop 1111
nodes i
Definition: Ztrack.h:44
real(8), dimension(:,:,:,:), allocatable, save num
Definition: chook.f:40
subroutine cmkincident(incident, fin)
Definition: cmkIncident.f:5
subroutine cgethes(from)
Definition: chookFlesh.f:322
subroutine cinitracking(incident)
Definition: ciniTracking.f:2
Definition: Zcoord.h:43
subroutine cobshes
Definition: chookFlesh.f:335
Here is the call graph for this function:

◆ cgethes()

subroutine cgethes ( integer  from)

Definition at line 297 of file chookFlesh.f.

References o.

297  implicit none
298 #include "../../SkelFlesh/Zprivate.h"
299 #include "Zmaxdef.h"
300 #include "Ztrack.h"
301  integer from
302 
303  integer i
304 
305  read(from) np
306  do i = 1, np
307  read(from) o(i)
308  enddo
nodes i
struct ob o[NpMax]
Definition: Zprivate.h:34
integer maxbin nregion c minsave drx2 ! drawing region real maxdep integer maxpos integer from
Definition: Zfit.h:15

◆ chookbgevent()

subroutine chookbgevent ( )

Definition at line 65 of file chookFlesh.f.

References cbegin1ev(), cerrormsg(), and cpushinci().

65  implicit none
66 #include "../../SkelFlesh/Zprivate.h"
67 
68 
69  integer nomore
70  call cbegin1ev( nomore )
71  if( nomore .eq. 1) then
72  call cerrormsg('all events have been fleshed', 1)
73  stop !!!!!!!!!!!!
74  endif
75  call cpushinci
76 
77  call xbgevent
78  call xihist ! instanciate histogram
79  call xclearnrfai ! clear Nrfai region
subroutine cerrormsg(msg, needrtn)
Definition: cerrorMsg.f:4
subroutine cpushinci
Definition: chookFlesh.f:314
subroutine cbegin1ev(nomore)
Definition: chookFlesh.f:67
Here is the call graph for this function:

◆ chookbgrun()

subroutine chookbgrun ( )

Definition at line 19 of file chookFlesh.f.

References cerrormsg(), cgetfname(), copenfw2(), cprintobs(), cprintprim(), cquhookc(), cquhooki(), cwriteparam(), parameter(), and xbgrun().

19  implicit none
20 #include "Zmaxdef.h"
21 #include "Zmanagerp.h"
22 #include "../../SkelFlesh/Zprivate.h"
23 #include "Ztrack.h"
24 #if defined (DOMPI)
25 #include "mpif.h"
26 #include "Zmpi.h"
27 #else
28  integer mpirank
30 #endif
31  real*8 temp
32  character*100 msg
33  integer icon
34  integer i
35  eventno = 0
36 
37  if(mpirank .eq. 0) then
38 ! namelist output
39  call cwriteparam(errorout, 0)
40 ! primary information
41  call cprintprim(errorout)
42 ! observation level information
43  call cprintobs(errorout)
44  endif
45 
46  call cquhooki(1, mdev) ! get skeleton memo dev #
47  call cquhookc(1, msg) ! get file name for sekelton data
48  call cgetfname(msg, mskel) ! add host name etc if needed
49  call copenfw2(mdev, mskel, 2, icon)
50  if(icon .ne. 1) then
51  call cerrormsg(mskel,1)
52  call cerrormsg(' could not be opened',0)
53  endif
54 
55  call xbgrun
56 
integer npitbl real *nx parameter(n=101, npitbl=46, nx=n-1) real *8 uconst
subroutine cgetfname(fnin, fn)
Definition: copenf.f:275
subroutine cerrormsg(msg, needrtn)
Definition: cerrorMsg.f:4
subroutine cprintobs(io)
Definition: cprintObs.f:2
nodes i
subroutine cprintprim(out)
Definition: cprintPrim.f:3
subroutine cwriteparam(io, force)
Definition: cwriteParam.f:4
integer mpirank
Definition: Zmpibasic.h:1
subroutine copenfw2(io, fnin, form, icon)
Definition: copenf.f:205
real(4), dimension(:), allocatable, save temp
Definition: cNRLAtmos.f:29
subroutine cquhooki(i, iv)
Definition: cqUHookr.f:15
subroutine cquhookc(i, cv)
Definition: cqUHookr.f:28
subroutine xbgrun
Definition: interface.f:10
Here is the call graph for this function:

◆ chookeint()

subroutine chookeint ( integer  never)

Definition at line 272 of file chookFlesh.f.

272  implicit none
273  integer never ! input & output
274  never = 1

◆ chookenevent()

subroutine chookenevent ( )

Definition at line 200 of file chookFlesh.f.

References cfintracking().

200 
201  implicit none
202 #include "Ztrack.h"
203 #include "Ztrackv.h"
204 #include "Zobs.h"
205 #include "Zobsp.h"
206 #include "Zobsv.h"
207 
208  integer i
209 
210 ! for Job ='newflesh', we must call cfinTracking ourselves.
211  call cfintracking
212 ! end of 1 event; if you need to do some here is
213 ! the place
214 
215  call xenevent
216 
217 
nodes i
subroutine cfintracking
Definition: cfinTracking.f:2
Here is the call graph for this function:

◆ chookenrun()

subroutine chookenrun ( )

Definition at line 225 of file chookFlesh.f.

References cprintstatus().

225  implicit none
226  call cprintstatus ! if don't like, comment out
subroutine cprintstatus
Definition: cendRun.f:35
Here is the call graph for this function:

◆ chookgint()

subroutine chookgint ( integer  never)

Definition at line 281 of file chookFlesh.f.

281  implicit none
282  integer never ! input & output
283  never = 1

◆ chooknepint()

subroutine chooknepint ( integer  never)

Definition at line 290 of file chookFlesh.f.

290  implicit none
291  integer never ! input & output
292  never = 1

◆ chookobs()

subroutine chookobs ( type(track aTrack,
integer  id 
)

Definition at line 153 of file chookFlesh.f.

References rndc(), and wgt.

153 !
154 ! Note that every real variable is in double precision so
155 ! that you may output it in sigle precision to save the memory.
156 ! In some cases it is essential to put it in sigle (say,
157 ! for gnuplot).
158 !
159  implicit none
160 #include "Zcode.h"
161 #include "Ztrack.h"
162 #include "Zprivate4.h"
163 
164  integer id ! input. 2 ==> reached at an observation level
165 ! 1 ==> aTrack is going out from
166 ! outer boundery.
167 ! 2 ==> reached at an observation level
168 ! 3 ==> reached at inner boundery.
169  type(track):: atrack
170 !
171 !
172  integer n, i
173  real*8 eps, u
174 
175  if( keepweight ) then
176  n = 1
177  else
178  if(atrack.wgt .gt. 1.) then
179  n=atrack.wgt
180  eps = atrack.wgt - n
181  call rndc(u)
182  if(u .lt. eps) then
183  n = n + 1
184  endif
185  else
186  n = 1
187  endif
188  endif
189 
190  do i = 1, n
191  call xobs(atrack, id)
192  enddo
193 
nodes i
Definition: Ztrack.h:44
subroutine rndc(u)
Definition: rnd.f:91
integer n
Definition: Zcinippxc.h:1
*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
Here is the call graph for this function:

◆ chooktrace()

subroutine chooktrace ( )

Definition at line 237 of file chookFlesh.f.

References height.

237  implicit none
238 
239 #include "Ztrack.h"
240 #include "Ztrackv.h"
241 #include "Ztrackp.h"
242 #include "Zobs.h"
243 #include "Zobsv.h"
244 
245  real*4 h1, h2
246 !
247 ! Every time a particle is moved in the atmosphere, this routine is called,
248 ! if trace > 100
249 ! For a one track segment,
250 ! TrackBefMove has track information at the beginning of the segment.
251 ! MoveTrack has track information at the end of the segment.
252 !
253 ! You can know the information a track contains in the
254 ! chookObs routine. (Note however, no conversion of coordinate
255 ! has been done. The values are in the Earth xyz system.)
256 ! Besides quantities explained there, you can use, for a given 'track'
257 !
258 ! atrack.pos.xyz.x, atrack.pos.xyz.y, atrack.pos.xyz.z (x,y.z)
259 ! atrack.pos.radiallen (distance from the center of the earth)
260 ! atrack.pos.depth (vertical depth)
261 ! atrack.pos.height (vertical heigth from sea level)
262 !
263 
264  h1 = trackbefmove.pos.height- obssites(noofsites).pos.height
265  h2 = movedtrack.pos.height - obssites(noofsites).pos.height
266 
*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

◆ cobshes()

subroutine cobshes ( )

Definition at line 312 of file chookFlesh.f.

References charge, chookobs(), code, coszenith, erg, false, mass, o, p, r, subcode, t, true, wgt, x, xyz, and y.

312  implicit none
313 #include "../../SkelFlesh/Zprivate.h"
314 #include "Ztrack.h"
315 !
316 ! memorized high energy showers at the skeleton making
317 ! time is put into the chookObs as if they are really observed
318  type(track):: atrack
319 
320  integer i
321  logical heobs ! if T, currently observing
322  common /zheobs/ heobs ! particles those obsrved at skeelton making time
323 
324  heobs = .true.
325  do i = 1, np
326  atrack.where = o(i).where
327  atrack.p.code = o(i).code
328  atrack.p.subcode = o(i).subcode
329  atrack.p.charge = o(i).charge
330  atrack.t = o(i).atime
331  atrack.p.fm.p(4) = o(i).erg
332  atrack.p.mass = o(i).mass
333  atrack.pos.xyz.r(1) = o(i).x
334  atrack.pos.xyz.r(2) = o(i).y
335  atrack.vec.w.r(1) = o(i).wx
336  atrack.vec.w.r(2) = o(i).wy
337  atrack.vec.w.r(3) = o(i).wz
338  atrack.vec.coszenith = o(i).zenith
339  atrack.wgt = 1.
340  call chookobs(atrack, 2)
341  enddo
342  heobs = .false.
*Zfirst p fm *Zfirst p Zfirst p Zfirst p *Zfirst p *Zfirst pos xyz r
Definition: ZavoidUnionMap.h:1
nodes i
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
Definition: Ztrack.h:44
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
********************block data cblkHeavy ********************integer j data *HeavyG2symbol p
Definition: cblkHeavy.h:7
*Zfirst p fm *Zfirst p Zfirst p code
Definition: ZavoidUnionMap.h:1
struct ob o[NpMax]
Definition: Zprivate.h:34
latitude latitude this system is used *****************************************************************! type coord sequence union map real y
Definition: Zcoord.h:25
*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
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
*Zfirst p fm *Zfirst p mass
Definition: ZavoidUnionMap.h:1
*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
! 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
Here is the call graph for this function:

◆ cpushinci()

subroutine cpushinci ( )

Definition at line 349 of file chookFlesh.f.

References cinitstack(), cpush(), and csortstack().

349  implicit none
350 #include "../../SkelFlesh/Zprivate.h"
351 #include "Ztrack.h"
352 #include "Ztrackv.h"
353  integer i
354 
355  type(track)::atrack
356  type(track)::zsave
357  call cinitstack ! empty the stack
358 
359  read(mdev) nooflowe
360 !/////////
361 ! write(0,*) ' rank=',mpirank,' NoOfLowE=',NoOfLowE
362 !//////////
363  do i = 1, nooflowe
364 #if defined (KEKA) || (KEKB)
365  zsave = zfirst ! save; next destroys Zfirst
366  read(mdev)
367 #include "ZavoidUnionMap.h"
368  atrack = zfirst
369  zfirst = zsave ! restore
370 #else
371  read(mdev) atrack
372 #endif
373 ! aTrack is already complete track so push it directly.
374  call cpush(atrack)
375  enddo
376 ! sort stack dscendent order
377  call csortstack
378 
nodes i
Definition: Ztrack.h:44
subroutine cinitstack
Definition: cstack.f:76
subroutine csortstack
Definition: cstack.f:102
subroutine cpush(a)
Definition: cstack.f:4
Here is the call graph for this function: