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