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