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