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