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