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