]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliLego.cxx
AliGUIMaterial and AliGUIMedium added
[u/mrichter/AliRoot.git] / STEER / AliLego.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.8  1999/10/01 09:54:33  fca
19 Correct logics for Lego StepManager
20
21 Revision 1.7  1999/09/29 09:24:29  fca
22 Introduction of the Copyright and cvs Log
23
24 */
25
26 //////////////////////////////////////////////////////////////
27 //////////////////////////////////////////////////////////////
28 //                            
29 //  Utility class to evaluate the material budget from
30 //  a given radius to the surface of an arbitrary cylinder
31 //  along radial directions from the centre:
32 // 
33 //   - radiation length
34 //   - Interaction length
35 //   - g/cm2
36 // 
37 //  Geantinos are shot in the bins in the fNtheta bins in theta
38 //  and fNphi bins in phi with specified rectangular limits.
39 //  The statistics are accumulated per
40 //    fRadMin < r < fRadMax    and  <0 < z < fZMax
41 //
42 //  To activate this option, you can do:
43 //      Root > gAlice.RunLego();
44 //  or  Root > .x menu.C  then select menu item "RunLego"
45 //  Note that when calling gAlice->RunLego, an optional list
46 //  of arguments may be specified.
47 //
48 //Begin_Html
49 /*
50 <img src="picts/alilego.gif">
51 */
52 //End_Html
53 //
54 //////////////////////////////////////////////////////////////
55
56 #include "TMath.h"
57 #include "AliLego.h"
58 #include "AliRun.h"
59 #include "AliConst.h"
60
61 ClassImp(AliLego)
62
63
64 //___________________________________________
65 AliLego::AliLego()
66 {
67    fHistRadl = 0;
68    fHistAbso = 0;
69    fHistGcm2 = 0;
70    fHistReta = 0;
71 }
72
73 //___________________________________________
74 AliLego::AliLego(const char *name, const char *title) 
75         : TNamed(name,title)
76 {
77    fHistRadl = 0;
78    fHistAbso = 0;
79    fHistGcm2 = 0;
80    fHistReta = 0;
81 }
82
83 //___________________________________________
84 AliLego::~AliLego()
85 {
86    delete fHistRadl;
87    delete fHistAbso;
88    delete fHistGcm2;
89    delete fHistReta;
90 }
91
92 //___________________________________________
93 void AliLego::GenerateKinematics()
94 {
95 // Create a geantino with kinematics corresponding to the current
96 // bins in theta and phi.
97    
98   //
99   // Rootinos are 0
100    const Int_t mpart = 0;
101    Float_t orig[3], pmom[3];
102    Float_t t, cost, sint, cosp, sinp;
103    
104 // --- Set to 0 radiation length, absorption length and g/cm2 ---
105    fTotRadl = 0;
106    fTotAbso = 0;
107    fTotGcm2 = 0;
108
109    fCurTheta = (fThetaMin+(fThetaBin-0.5)*(fThetaMax-fThetaMin)/fNtheta)*kDegrad;
110    fCurPhi   = (fPhiMin+(fPhiBin-0.5)*(fPhiMax-fPhiMin)/fNphi)*kDegrad;
111    cost      = TMath::Cos(fCurTheta);
112    sint      = TMath::Sin(fCurTheta);
113    cosp      = TMath::Cos(fCurPhi);
114    sinp      = TMath::Sin(fCurPhi);
115
116    pmom[0] = cosp*sint;
117    pmom[1] = sinp*sint;
118    pmom[2] = cost;
119
120 // --- Where to start
121    orig[0] = orig[1] = orig[2] = 0;
122    Float_t dalicz = 3000;
123    if (fRadMin > 0) {
124       t = PropagateCylinder(orig,pmom,fRadMin,dalicz);
125       orig[0] = pmom[0]*t;
126       orig[1] = pmom[1]*t;
127       orig[2] = pmom[2]*t;
128       if (TMath::Abs(orig[2]) > fZMax) return;
129    }
130
131 // --- We do start here
132    Float_t polar[3]={0.,0.,0.};
133    Int_t ntr;
134    gAlice->SetTrack(1, 0, mpart, pmom, orig, polar, 0, "LEGO ray", ntr);
135 }
136
137 //___________________________________________
138 void AliLego::Init(Int_t ntheta,Float_t themin,Float_t themax,
139           Int_t nphi,Float_t phimin,Float_t phimax,Float_t rmin,Float_t rmax,
140                   Float_t zmax)
141 {
142 // specify the angular limits and the size of the rectangular box
143    fNtheta   = ntheta;
144    fThetaMin = themin;
145    fThetaMax = themax;
146    fNphi     = nphi;
147    fPhiMin   = phimin;
148    fPhiMax   = phimax;
149    fRadMin   = rmin;
150    fRadMax   = rmax;
151    fZMax     = zmax;
152    Float_t etamin = -TMath::Log(TMath::Tan(TMath::Min((Double_t)fThetaMax*kDegrad/2,TMath::Pi()/2-1.e-10)));
153    Float_t etamax = -TMath::Log(TMath::Tan(TMath::Max((Double_t)fThetaMin*kDegrad/2,              1.e-10)));
154
155    fHistRadl = new TH2F("hradl","Radiation length map",    nphi,phimin,phimax,ntheta,themin,themax);
156    fHistAbso = new TH2F("habso","Interaction length map",  nphi,phimin,phimax,ntheta,themin,themax);
157    fHistGcm2 = new TH2F("hgcm2","g/cm2 length map",        nphi,phimin,phimax,ntheta,themin,themax);
158    fHistReta = new TH2F("hetar","Radiation length vs. eta",nphi,phimin,phimax,ntheta,etamin,etamax);
159
160 }
161
162 //___________________________________________
163 Float_t AliLego::PropagateCylinder(Float_t *x, Float_t *v, Float_t r, Float_t z)
164 {
165 // Propagate to cylinder from inside
166
167    Double_t hnorm, sz, t, t1, t2, t3, sr;
168    Double_t d[3];
169    const Float_t kSmall  = 1e-8;
170    const Float_t kSmall2 = kSmall*kSmall;
171
172 // ---> Find intesection with Z planes
173    d[0]  = v[0];
174    d[1]  = v[1];
175    d[2]  = v[2];
176    hnorm = TMath::Sqrt(1/(d[0]*d[0]+d[1]*d[1]+d[2]*d[2]));
177    d[0] *= hnorm;
178    d[1] *= hnorm;
179    d[2] *= hnorm;
180    if (d[2] > kSmall)       sz = (z-x[2])/d[2];
181    else if (d[2] < -kSmall) sz = -(z+x[2])/d[2];
182    else                     sz = 1.e10;  // ---> Direction parallel to X-Y, no intersection
183
184 // ---> Intersection with cylinders
185 //      Intersection point (x,y,z)
186 //      (x,y,z) is on track :    x=X(1)+t*D(1)
187 //                               y=X(2)+t*D(2)
188 //                               z=X(3)+t*D(3)
189 //      (x,y,z) is on cylinder : x**2 + y**2 = R**2
190 //
191 //      (D(1)**2+D(2)**2)*t**2
192 //      +2.*(X(1)*D(1)+X(2)*D(2))*t
193 //      +X(1)**2+X(2)**2-R**2=0
194 // ---> Solve second degree equation
195    t1 = d[0]*d[0] + d[1]*d[1];
196    if (t1 <= kSmall2) {
197       t = sz;  // ---> Track parallel to the z-axis, take distance to planes
198    } else {
199       t2 = x[0]*d[0] + x[1]*d[1];
200       t3 = x[0]*x[0] + x[1]*x[1];
201       // ---> It should be positive, but there may be numerical problems
202       sr = (t2 +TMath::Sqrt(TMath::Max(t2*t2-(t3-r*r)*t1,0.)))/t1;
203       // ---> Find minimum distance between planes and cylinder
204       t  = TMath::Min(sz,sr);
205    }
206    return t;
207 }
208
209 //___________________________________________
210 void AliLego::Run()
211 {
212    // loop on phi,theta bins
213    gMC->InitLego();
214    Float_t thed, phid, eta;
215    for (fPhiBin=1; fPhiBin<=fNphi; fPhiBin++) {
216       printf("AliLego Generating rays in phi bin:%d\n",fPhiBin);
217       for (fThetaBin=1; fThetaBin<=fNtheta; fThetaBin++) {
218          gMC->Gtrigi();
219          gMC->Gtrigc();
220          GenerateKinematics();
221          gMC->Gtreve_root();
222
223          thed = fCurTheta*kRaddeg;
224          phid = fCurPhi*kRaddeg;
225          eta  = -TMath::Log(TMath::Tan(TMath::Max(
226                      TMath::Min((Double_t)fCurTheta/2,TMath::Pi()/2-1.e-10),1.e-10)));
227          fHistRadl->Fill(phid,thed,fTotRadl);
228          fHistAbso->Fill(phid,thed,fTotAbso);
229          fHistGcm2->Fill(phid,thed,fTotGcm2);
230          fHistReta->Fill(phid,eta,fTotRadl);
231          gAlice->FinishEvent();
232       }
233    }
234    // store histograms in current Root file
235    fHistRadl->Write();
236    fHistAbso->Write();
237    fHistGcm2->Write();
238    fHistReta->Write();
239 }
240
241 //___________________________________________
242 void AliLego::StepManager()
243 {
244 // called from AliRun::Stepmanager from gustep.
245 // Accumulate the 3 parameters step by step
246   
247    static Float_t t;
248    Float_t a,z,dens,radl,absl;
249    Int_t i;
250    
251    Float_t step  = gMC->TrackStep();
252        
253    Float_t vect[3], dir[3];
254    TLorentzVector pos, mom;
255
256    gMC->TrackPosition(pos);  
257    gMC->TrackMomentum(mom);
258    gMC->CurrentMaterial(a,z,dens,radl,absl);
259    
260    if (z < 1) return;
261    
262 // --- See if we have to stop now
263    if (TMath::Abs(pos[2]) > fZMax  || 
264        pos[0]*pos[0] +pos[1]*pos[1] > fRadMax*fRadMax) {
265      if (gMC->TrackLength()) {
266        // Not the first step, add past contribution
267        fTotAbso += t/absl;
268        fTotRadl += t/radl;
269        fTotGcm2 += t*dens;
270      }
271      gMC->StopTrack();
272      return;
273    }
274
275 // --- See how long we have to go
276    for(i=0;i<3;++i) {
277      vect[i]=pos[i];
278      dir[i]=mom[i];
279    }
280
281    t  = PropagateCylinder(vect,dir,fRadMax,fZMax);
282
283    fTotAbso += step/absl;
284    fTotRadl += step/radl;
285    fTotGcm2 += step*dens;
286 }
287