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