]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TGeant3/AliGeant3.cxx
Introducing new Rndm and QA classes
[u/mrichter/AliRoot.git] / TGeant3 / AliGeant3.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /*
17 $Log$
18 Revision 1.12  2000/11/06 11:35:46  morsch
19 Call BuildGeometry() after Init() to be able to share common detector parameters.
20
21 Revision 1.11  2000/10/04 16:30:22  fca
22 Add include for exit()
23
24 Revision 1.10  2000/10/02 21:28:16  fca
25 Removal of useless dependecies via forward declarations
26
27 Revision 1.9  2000/09/12 14:36:17  morsch
28 In gudcay(): call ForceDecay() before Decay()
29
30 Revision 1.8  2000/09/06 14:56:34  morsch
31 gudcay() method implemented.
32 Decays are performed using  the AliDecayer interface. The pointer to the instance of AliDecayer
33 is obtained from geant3 (will be gMC once it has been implemented for Geant4).
34
35 Revision 1.7  2000/07/13 16:19:10  fca
36 Mainly coding conventions + some small bug fixes
37
38 Revision 1.6  2000/07/11 18:24:59  fca
39 Coding convention corrections + few minor bug fixes
40
41 Revision 1.5  2000/05/20 14:49:48  fca
42 Call gdebug at the end of gustep
43
44 Revision 1.4  2000/04/26 10:17:32  fca
45 Changes in Lego for G4 compatibility
46
47 Revision 1.3  2000/04/07 11:12:35  fca
48 G4 compatibility changes
49
50 Revision 1.2  2000/02/29 19:11:17  fca
51 Move gucode into AliGeant3.cxx
52
53 Revision 1.1  2000/02/23 16:25:25  fca
54 AliVMC and AliGeant3 classes introduced
55 ReadEuclid moved from AliRun to AliModule
56
57 */
58
59 #include <stdlib.h>
60
61 #include <TParticle.h>
62
63 #include "AliDecayer.h"
64 #include "AliGeant3.h"
65 #include "AliRun.h"
66 #include "TGeant3.h"
67 #include "AliCallf77.h" 
68 #include "AliModule.h"
69 #include "AliMagF.h"
70 #include "AliGenerator.h"
71
72 #ifndef WIN32
73
74 # define rxgtrak rxgtrak_
75 # define rxouth  rxouth_
76 #else
77
78 # define rxgtrak RXGTRAK 
79 # define rxouth  RXOUTH
80 #endif
81
82 ClassImp(AliGeant3)
83
84 AliGeant3::AliGeant3(const char *title) : 
85   TGeant3(title) {}
86
87 //____________________________________________________________________________
88 void AliGeant3::FinishGeometry()
89 {
90   //
91   // Finalise geometry construction
92   //
93   TGeant3::FinishGeometry();
94   //Create the color table
95   SetColors();
96 }
97
98 //____________________________________________________________________________
99 void AliGeant3::Init()
100 {
101   //
102   //=================Create Materials and geometry
103   //
104   TObjArray *modules = gAlice->Modules();
105   TIter next(modules);
106   AliModule *detector;
107   while((detector = (AliModule*)next())) {
108     // Initialise detector materials and geometry
109     detector->CreateMaterials();
110     detector->CreateGeometry();
111     detector->Init();
112     detector->BuildGeometry();
113   }
114
115   //Terminate building of geometry
116   FinishGeometry();
117 }
118
119 //____________________________________________________________________________
120 void AliGeant3::ProcessRun(Int_t nevent)
121 {
122   //
123   // Process the run
124   //
125   Int_t todo = TMath::Abs(nevent);
126   for (Int_t i=0; i<todo; i++) {
127   // Process one run (one run = one event)
128      gAlice->BeginEvent();
129      ProcessEvent();
130      gAlice->FinishEvent();
131   }
132 }
133
134 //_____________________________________________________________________________
135 void AliGeant3::ProcessEvent()
136 {
137   //
138   // Process one event
139   //
140   Gtrigi();
141   Gtrigc();
142   Gtrig();
143 }
144
145 //_____________________________________________________________________________
146 void AliGeant3::SetColors()
147 {
148   //
149   // Set the colors for all the volumes
150   // this is done sequentially for all volumes
151   // based on the number of their medium
152   //
153   Int_t kv, icol;
154   Int_t jvolum=fGclink->jvolum;
155   //Int_t jtmed=fGclink->jtmed;
156   //Int_t jmate=fGclink->jmate;
157   Int_t nvolum=fGcnum->nvolum;
158   char name[5];
159   //
160   //    Now for all the volumes
161   for(kv=1;kv<=nvolum;kv++) {
162     //     Get the tracking medium
163     Int_t itm=Int_t (fZq[fZlq[jvolum-kv]+4]);
164     //     Get the material
165     //Int_t ima=Int_t (fZq[fZlq[jtmed-itm]+6]);
166     //     Get z
167     //Float_t z=fZq[fZlq[jmate-ima]+7];
168     //     Find color number
169     //icol = Int_t(z)%6+2;
170     //icol = 17+Int_t(z*150./92.);
171     //icol = kv%6+2;
172     icol = itm%6+2;
173     strncpy(name,(char*)&fZiq[jvolum+kv],4);
174     name[4]='\0';
175     Gsatt(name,"COLO",icol);
176   }
177 }
178
179 //_____________________________________________________________________________
180 //
181 //                 Interfaces to Fortran
182 //
183 //_____________________________________________________________________________
184
185 extern "C" void type_of_call  rxgtrak (Int_t &mtrack, Int_t &ipart, Float_t *pmom, 
186                                        Float_t &e, Float_t *vpos, Float_t *polar,
187                                        Float_t &tof)
188 {
189   //
190   //     Fetches next track from the ROOT stack for transport. Called by the
191   //     modified version of GTREVE.
192   //
193   //              Track number in the ROOT stack. If MTRACK=0 no
194   //      mtrack  more tracks are left in the stack to be
195   //              transported.
196   //      ipart   Particle code in the GEANT conventions.
197   //      pmom[3] Particle momentum in GeV/c
198   //      e       Particle energy in GeV
199   //      vpos[3] Particle position
200   //      tof     Particle time of flight in seconds
201   //
202   Int_t pdg;
203   gAlice->GetNextTrack(mtrack, pdg, pmom, e, vpos, polar, tof);
204   ipart = gMC->IdFromPDG(pdg);
205   mtrack++;
206 }
207
208 //_____________________________________________________________________________
209 extern "C" void type_of_call  rxouth ()
210 {
211   //
212   // Called by Gtreve at the end of each primary track
213   //
214   gAlice->FinishPrimary();
215 }
216
217 #ifndef WIN32
218 #  define gudigi gudigi_
219 #  define guhadr guhadr_
220 #  define guout  guout_
221 #  define guphad guphad_
222 #  define gudcay gudcay_
223 #  define guiget guiget_
224 #  define guinme guinme_
225 #  define guinti guinti_
226 #  define gunear gunear_
227 #  define guskip guskip_
228 #  define guview guview_
229 #  define gupara gupara_
230 #  define gudtim gudtim_
231 #  define guplsh guplsh_
232 #  define gutrev gutrev_
233 #  define gutrak gutrak_
234 #  define guswim guswim_
235 #  define gufld  gufld_
236 #  define gustep gustep_
237 #  define gukine gukine_
238 #  define uglast uglast_
239
240 #  define gheish gheish_
241 #  define flufin flufin_
242 #  define gfmfin gfmfin_
243 #  define gpghei gpghei_
244 #  define fldist fldist_
245 #  define gfmdis gfmdis_
246 #  define ghelx3 ghelx3_
247 #  define ghelix ghelix_
248 #  define grkuta grkuta_
249 #  define gtrack gtrack_
250 #  define gtreveroot gtreveroot_
251 #  define glast  glast_
252
253 #else
254 #  define gudigi GUDIGI
255 #  define guhadr GUHADR
256 #  define guout  GUOUT
257 #  define guphad GUPHAD
258 #  define gudcay GUDCAY
259 #  define guiget GUIGET
260 #  define guinme GUINME
261 #  define guinti GUINTI
262 #  define gunear GUNEAR
263 #  define guskip GUSKIP
264 #  define guview GUVIEW
265 #  define gupara GUPARA
266 #  define gudtim GUDTIM
267 #  define guplsh GUPLSH
268 #  define gutrev GUTREV
269 #  define gutrak GUTRAK
270 #  define guswim GUSWIM
271 #  define gufld  GUFLD
272 #  define gustep GUSTEP
273 #  define gukine GUKINE
274 #  define uglast UGLAST
275
276 #  define gheish GHEISH
277 #  define flufin FLUFIN
278 #  define gfmfin GFMFIN
279 #  define gpghei GPGHEI
280 #  define fldist FLDIST
281 #  define gfmdis GFMDIS
282 #  define ghelx3 GHELX3
283 #  define ghelix GHELIX
284 #  define grkuta GRKUTA
285 #  define gtrack GTRACK
286 #  define gtreveroot GTREVEROOT
287 #  define glast  GLAST
288
289 #endif
290
291 extern "C" type_of_call void gheish();
292 extern "C" type_of_call void flufin();
293 extern "C" type_of_call void gfmfin();
294 extern "C" type_of_call void gpghei();
295 extern "C" type_of_call void fldist();
296 extern "C" type_of_call void gfmdis();
297 extern "C" type_of_call void ghelx3(Float_t&, Float_t&, Float_t*, Float_t*);
298 extern "C" type_of_call void ghelix(Float_t&, Float_t&, Float_t*, Float_t*);
299 extern "C" type_of_call void grkuta(Float_t&, Float_t&, Float_t*, Float_t*);
300 extern "C" type_of_call void gtrack();
301 extern "C" type_of_call void gtreveroot();
302 extern "C" type_of_call void glast();
303
304 extern "C" type_of_call {
305
306 //______________________________________________________________________
307 void gudigi() 
308 {
309 //
310 //    ******************************************************************
311 //    *                                                                *
312 //    *       User routine to digitize one event                       *
313 //    *                                                                *
314 //    *    ==>Called by : GTRIG                                        *
315 //    *                                                                *
316 //    ******************************************************************
317
318 }
319
320
321 //______________________________________________________________________
322 void guhadr()
323 {
324 //
325 //    ******************************************************************
326 //    *                                                                *
327 //    *       User routine to generate one hadronic interaction        *
328 //    *                                                                *
329 //    *    ==>Called by : GTHADR,GTNEUT                                *
330 //    *                                                                *
331 //    ******************************************************************
332 //
333 //
334 //    ------------------------------------------------------------------
335 //
336       TGeant3* geant3 = (TGeant3*) gMC;
337       Int_t ihadr=geant3->Gcphys()->ihadr;
338       if (ihadr<4)       gheish();
339       else if (ihadr==4) flufin();
340       else               gfmfin();
341 }
342
343 //______________________________________________________________________
344 void guout()
345 {
346 //
347 //    ******************************************************************
348 //    *                                                                *
349 //    *       User routine called at the end of each event             *
350 //    *                                                                *
351 //    *    ==>Called by : GTRIG                                        *
352 //    *                                                                *
353 //    ******************************************************************
354 //
355 //
356 //    ------------------------------------------------------------------
357 //
358 }
359
360 //______________________________________________________________________
361 void guphad()
362 {
363 //
364 //    ******************************************************************
365 //    *                                                                *
366 //    *       User routine to compute Hadron. inter. probabilities     *
367 //    *                                                                *
368 //    *    ==>Called by : GTHADR,GTNEUT                                *
369 //    *                                                                *
370 //    ******************************************************************
371 //
372 //
373 //    ------------------------------------------------------------------
374 //
375       TGeant3* geant3 = (TGeant3*) gMC;
376       Int_t ihadr=geant3->Gcphys()->ihadr;
377       if (ihadr<4)       gpghei();
378       else if (ihadr==4) fldist();
379       else               gfmdis();
380 }
381
382 //______________________________________________________________________
383 void gudcay()
384 {
385 //
386 //    ******************************************************************
387 //    *                                                                *
388 //    *       User routine to decay particles                          *
389 //    *                                                                *
390 //    *    ==>Called by : GDECAY                                       *
391 //    *                                                                *
392 //    ******************************************************************
393 //
394 //
395 //    ------------------------------------------------------------------
396 //
397     
398     TGeant3* geant3=(TGeant3*) gMC;
399   // set decay table
400   gMC->Decayer()->ForceDecay();
401
402 // Initialize 4-momentum vector    
403     Int_t ipart = geant3->Gckine()->ipart;
404     TLorentzVector p;
405     
406     p[0]=geant3->Gctrak()->vect[3];
407     p[1]=geant3->Gctrak()->vect[4];
408     p[2]=geant3->Gctrak()->vect[5];
409     p[3]=geant3->Gctrak()->vect[6];    
410     
411 // Convert from geant to lund particle code
412     Int_t iplund=gMC->PDGFromId(ipart);
413 // Particle list
414     static TClonesArray *particles;
415     if(!particles) particles=new TClonesArray("TParticle",1000);
416 // Decay
417     gMC->Decayer()->Decay(iplund, &p);
418     
419 // Fetch Particles
420     Int_t np = geant3->Decayer()->ImportParticles(particles);
421     if (np <=1) return;
422     
423     for (Int_t i=0; i< np; i++) 
424     {
425         TParticle *  iparticle = (TParticle *) particles->At(i);
426         Int_t ks = iparticle->GetStatusCode();
427 // Final state particles only
428         if (ks < 1 || ks >  10) continue;
429 // Skip neutrinos
430         Int_t  kf=iparticle->GetPdgCode();
431         if (kf==12 || kf ==-12) continue;
432         if (kf==14 || kf ==-14) continue;
433         if (kf==16 || kf ==-16) continue;
434         
435         Int_t index=geant3->Gcking()->ngkine;
436 // Put particle on geant stack
437 // momentum vector
438
439         (geant3->Gcking()->gkin[index][0]) = iparticle->Px();
440         (geant3->Gcking()->gkin[index][1]) = iparticle->Py();
441         (geant3->Gcking()->gkin[index][2]) = iparticle->Pz();
442         (geant3->Gcking()->gkin[index][3]) = iparticle->Energy();
443         Int_t ilu = gMC->IdFromPDG(kf);
444
445 // particle type        
446         (geant3->Gcking()->gkin[index][4]) = Float_t(ilu);
447 // position
448         (geant3->Gckin3()->gpos[index][0]) = geant3->Gctrak()->vect[0];
449         (geant3->Gckin3()->gpos[index][1]) = geant3->Gctrak()->vect[1];
450         (geant3->Gckin3()->gpos[index][2]) = geant3->Gctrak()->vect[2];
451 // time of flight offset (mm)
452         (geant3->Gcking()->tofd[index])    = 0.;
453 //      (geant3->Gcking()->tofd[index])    = iparticle->T()/(10*kSpeedOfLight);
454 // increase stack counter
455         (geant3->Gcking()->ngkine)=index+1;
456     }
457 }
458
459 //______________________________________________________________________
460 void guiget(Int_t&, Int_t&, Int_t&)
461 {
462 //
463 //    ******************************************************************
464 //    *                                                                *
465 //    *       User routine for interactive control of GEANT            *
466 //    *                                                                *
467 //    *    ==>Called by : <GXINT>, GINCOM                              *
468 //    *                                                                *
469 //    ******************************************************************
470 //
471 //
472 //    ------------------------------------------------------------------
473 //
474 }
475
476 //______________________________________________________________________
477 void guinme(Float_t*, Int_t&, Float_t*, Int_t& IYES)
478 {
479 //
480 //    **********************************************
481 //    *                                            *
482 //    *    USER ROUTINE TO PROVIDE GINME FUNCTION  *
483 //    *    FOR ALL USER SHAPES IDENTIFIED BY THE   *
484 //    *    SHAPE NUMBER SH. POINT IS GIVEN IN X    *
485 //    *    THE PARAMETERS ARE GIVEN IN P. IYES IS  *
486 //    *    RETURNED 1 IF POINT IS IN, 0 IF POINT   *
487 //    *    IS OUT AND LESS THAN ZERO IF SHAPE      *
488 //    *    NUMBER IS NOT SUPPORTED.                *
489 //    *                                            *
490 //    *    ==>Called by : GINME                    *
491 //    *                                            *
492 //    **********************************************
493 //
494       IYES=-1;
495 }
496
497 //______________________________________________________________________
498 void guinti()
499 {
500 //
501 //    ******************************************************************
502 //    *                                                                *
503 //    *       User routine for interactive version                     *
504 //    *                                                                *
505 //    *    ==>Called by : <GXINT>,  GINTRI                             *
506 //    *                                                                *
507 //    ******************************************************************
508 //
509 //
510 //    ------------------------------------------------------------------
511 //
512 }
513
514 //______________________________________________________________________
515 void gunear(Int_t&, Int_t&, Float_t*, Int_t&)
516 {
517 //
518 //    ******************************************************************
519 //    *                                                                *
520 //    *    User search                                                 *
521 //    *       ISEARC to identify the given volume                      *
522 //    *       ICALL  to identify the calling routine                   *
523 //    *              1 GMEDIA like                                     *
524 //    *              2 GNEXT like                                      *
525 //    *       X      coordinates (+direction for ICALL=2)              *
526 //    *       JNEAR  address of default list of neighbours             *
527 //    *              (list to be overwriten by user)                   *
528 //    *                                                                *
529 //    *    Called by : GFTRAC, GINVOL, GTMEDI, GTNEXT, GNEXT, GMEDIA   *
530 //    *                                                                *
531 //    ******************************************************************
532 //
533 //
534 //    ------------------------------------------------------------------
535 //
536 }
537
538 //______________________________________________________________________
539 void guskip(Int_t& ISKIP)
540 {
541 //
542 //    ******************************************************************
543 //    *                                                                *
544 //    *   User routine to skip unwanted tracks                         *
545 //    *                                                                *
546 //    *   Called by : GSSTAK                                           *
547 //    *   Author    : F.Bruyant                                        *
548 //    *                                                                *
549 //    ******************************************************************
550 //
551 //
552 //    ------------------------------------------------------------------
553 //
554       ISKIP = 0;
555 }
556
557 //______________________________________________________________________
558 void guswim(Float_t& CHARGE, Float_t& STEP, Float_t* VECT, Float_t* VOUT)
559 {
560 //
561 //    ******************************************************************
562 //    *                                                                *
563 //    *       User routine to control tracking of one track            *
564 //    *       in a magnetic field                                      *
565 //    *                                                                *
566 //    *    ==>Called by : GTELEC,GTHADR,GTMUON                         *
567 //    *                                                                *
568 //    ******************************************************************
569 //
570 //
571 //    ------------------------------------------------------------------
572 //
573   TGeant3* geant3 = (TGeant3*) gMC;
574   Int_t ifield=geant3->Gctmed()->ifield;
575   Float_t fieldm=geant3->Gctmed()->fieldm;
576
577   if (ifield==3) {
578     Float_t fldcharge = fieldm*CHARGE;
579     ghelx3(fldcharge,STEP,VECT,VOUT);
580   }
581   else if (ifield==2) ghelix(CHARGE,STEP,VECT,VOUT);
582   else                grkuta(CHARGE,STEP,VECT,VOUT);
583 }
584
585 //______________________________________________________________________
586 void guview(Int_t&, Int_t&, DEFCHARD, Int_t& DEFCHARL)
587 {
588 //
589 //    ******************************************************************
590 //    *                                                                *
591 //    *       User routine for interactive version                     *
592 //    *                                                                *
593 //    *    ==>Called by : <GXINT>, GINC1                               *
594 //    *                                                                *
595 //    ******************************************************************
596 //
597 //
598 //    ------------------------------------------------------------------
599 //
600 }
601
602 //______________________________________________________________________
603 void gupara()
604 {
605 //
606 //    ******************************************************************
607 //    *                                                                *
608 //    *       User routine called every time a particle falls below    *
609 //    *       parametrization threshold. This routine should create    *
610 //    *       the parametrization stack, and, when this is full,       *
611 //    *       parametrize the shower and track the geantinos.          *
612 //    *                                                                *
613 //    *    ==>Called by : GTRACK                                       *
614 //    *                                                                *
615 //    ******************************************************************
616 //
617 //
618 //    ------------------------------------------------------------------
619 //
620 }
621
622 //______________________________________________________________________
623 Float_t gudtim(Float_t&, Float_t&, Int_t&, Int_t&)
624 {
625 //
626 //    ******************************************************************
627 //    *                                                                *
628 //    *       User function called by GCDRIF to return drift time      *
629 //    *                                                                *
630 //    *    ==>Called by : GCDRIF                                       *
631 //    *                                                                *
632 //    ******************************************************************
633 //
634 //
635 //    ------------------------------------------------------------------
636 //
637       return 0;
638 }
639
640
641 //______________________________________________________________________
642 Float_t guplsh(Int_t&, Int_t&)
643 {
644 //
645 //    ******************************************************************
646 //    *                                                                *
647 //    *                                                                *
648 //    *    ==>Called by : GLISUR                                       *
649 //    *                                                                *
650 //    ******************************************************************
651 //
652 //
653 //    ------------------------------------------------------------------
654 //
655 //
656 //*** By default this defines perfect smoothness
657       return 1;
658 }
659
660 //______________________________________________________________________
661 void gutrak()
662 {
663 //
664 //    ******************************************************************
665 //    *                                                                *
666 //    *       User routine to control tracking of one track            *
667 //    *                                                                *
668 //    *    ==>Called by : GTREVE                                       *
669 //    *                                                                *
670 //    ******************************************************************
671 //
672 //
673 //    ------------------------------------------------------------------
674 //
675      gAlice->PreTrack();
676
677      gtrack();
678
679      gAlice->PostTrack();
680 }
681
682 //______________________________________________________________________
683 void gutrev()
684 {
685 //
686 //    ******************************************************************
687 //    *                                                                *
688 //    *       User routine to control tracking of one event            *
689 //    *                                                                *
690 //    *    ==>Called by : GTRIG                                        *
691 //    *                                                                *
692 //    ******************************************************************
693 //
694 //
695 //    ------------------------------------------------------------------
696 //
697   gtreveroot();
698 }
699
700
701 //______________________________________________________________________
702 void gufld(Float_t *x, Float_t *b)
703 {
704       if(gAlice->Field()) {
705          gAlice->Field()->Field(x,b);
706       } else {
707          printf("No mag field defined!\n");
708          b[0]=b[1]=b[2]=0;
709       }
710 }
711
712 //______________________________________________________________________
713 void gustep()
714 {
715 //
716 //    ******************************************************************
717 //    *                                                                *
718 //    *       User routine called at the end of each tracking step     *
719 //    *       INWVOL is different from 0 when the track has reached    *
720 //    *              a volume boundary                                 *
721 //    *       ISTOP is different from 0 if the track has stopped       *
722 //    *                                                                *
723 //    *    ==>Called by : GTRACK                                       *
724 //    *                                                                *
725 //    ******************************************************************
726 //
727
728
729   TLorentzVector x;
730   Float_t r;
731   Int_t ipp, jk, id, nt;
732   Float_t polar[3]={0,0,0};
733   Float_t mom[3];
734   AliMCProcess pProc;
735
736   
737   TGeant3* geant3 = (TGeant3*) gMC;
738   //     Stop particle if outside user defined tracking region 
739   gMC->TrackPosition(x);
740   r=TMath::Sqrt(x[0]*x[0]+x[1]*x[1]);
741   if (r > gAlice->TrackingRmax() || TMath::Abs(x[2]) > gAlice->TrackingZmax()) {
742         gMC->StopTrack();
743   }
744
745   // --- Add new created particles 
746   if (gMC->NSecondaries() > 0) {
747     pProc=gMC->ProdProcess();
748     for (jk = 0; jk < geant3->Gcking()->ngkine; ++jk) {
749       ipp = Int_t (geant3->Gcking()->gkin[jk][4]+0.5);
750       // --- Skip neutrinos! 
751       if (ipp != 4) {
752         gAlice->SetTrack(1,gAlice->CurrentTrack(),gMC->PDGFromId(ipp), geant3->Gcking()->gkin[jk], 
753                          geant3->Gckin3()->gpos[jk], polar,geant3->Gctrak()->tofg, pProc, nt);
754       }
755     }
756   }
757   // Cherenkov photons here
758   if ( geant3->Gckin2()->ngphot ) {
759     for (jk = 0; jk < geant3->Gckin2()->ngphot; ++jk) {
760       mom[0]=geant3->Gckin2()->xphot[jk][3]*geant3->Gckin2()->xphot[jk][6];
761       mom[1]=geant3->Gckin2()->xphot[jk][4]*geant3->Gckin2()->xphot[jk][6];
762       mom[2]=geant3->Gckin2()->xphot[jk][5]*geant3->Gckin2()->xphot[jk][6];
763       gAlice->SetTrack(1, gAlice->CurrentTrack(), gMC->PDGFromId(50),
764                        mom,                             //momentum
765                        geant3->Gckin2()->xphot[jk],     //position
766                        &geant3->Gckin2()->xphot[jk][7], //polarisation
767                        geant3->Gckin2()->xphot[jk][10], //time of flight
768                        kPCerenkov, nt);
769       }
770   }
771   // --- Particle leaving the setup ?
772   if (!gMC->IsTrackOut()) 
773     if ((id=gAlice->DetFromMate(geant3->Gctmed()->numed)) >= 0) gAlice->StepManager(id);
774
775   // --- Standard GEANT debug routine 
776   if(geant3->Gcflag()->idebug) geant3->Gdebug();
777 }
778
779 //______________________________________________________________________
780 void gukine ()
781 {
782 //
783 //    ******************************************************************
784 //    *                                                                *
785 //    *       Read or Generates Kinematics for primary tracks          *
786 //    *                                                                *
787 //    *    ==>Called by : GTRIG                                        *
788 //    *                                                                *
789 //    ******************************************************************
790 //
791 //
792 //    ------------------------------------------------------------------
793 //
794   gAlice->Generator()->Generate();
795 }
796
797
798 //______________________________________________________________________
799 void uglast()
800 {
801 //
802 //    ******************************************************************
803 //    *                                                                *
804 //    *       User routine called at the end of the run                *
805 //    *                                                                *
806 //    *    ==>Called by : GRUN                                         *
807 //    *                                                                *
808 //    ******************************************************************
809 //
810 //
811 }
812 }
813