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