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