1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 Revision 1.15 2000/05/16 13:10:40 fca
19 New method IsNewTrack and fix for a problem in Father-Daughter relations
21 Revision 1.14 2000/04/27 10:38:21 fca
22 Correct termination of Lego Run and introduce Lego getter in AliRun
24 Revision 1.13 2000/04/26 10:17:31 fca
25 Changes in Lego for G4 compatibility
27 Revision 1.12 2000/04/07 11:12:33 fca
28 G4 compatibility changes
30 Revision 1.11 2000/03/22 13:42:26 fca
31 SetGenerator does not replace an existing generator, ResetGenerator does
33 Revision 1.10 2000/02/23 16:25:22 fca
34 AliVMC and AliGeant3 classes introduced
35 ReadEuclid moved from AliRun to AliModule
37 Revision 1.9 1999/12/03 10:54:01 fca
40 Revision 1.8 1999/10/01 09:54:33 fca
41 Correct logics for Lego StepManager
43 Revision 1.7 1999/09/29 09:24:29 fca
44 Introduction of the Copyright and cvs Log
48 //////////////////////////////////////////////////////////////
49 //////////////////////////////////////////////////////////////
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:
56 // - Interaction length
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
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.
72 <img src="picts/alilego.gif">
76 //////////////////////////////////////////////////////////////
87 //___________________________________________
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)
102 // specify the angular limits and the size of the rectangular box
104 fGener = new AliLegoGenerator(ntheta, themin, themax,
105 nphi, phimin, phimax, rmin, rmax, zmax);
107 gAlice->ResetGenerator(fGener);
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)));
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);
123 //___________________________________________
130 gAlice->ResetGenerator(0);
134 //___________________________________________
135 void AliLego::BeginEvent()
137 // --- Set to 0 radiation length, absorption length and g/cm2 ---
143 //___________________________________________
144 void AliLego::FinishEvent()
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)));
153 fHistRadl->Fill(phid,thed,fTotRadl);
154 fHistAbso->Fill(phid,thed,fTotAbso);
155 fHistGcm2->Fill(phid,thed,fTotGcm2);
156 fHistReta->Fill(phid,eta,fTotRadl);
159 //___________________________________________
160 void AliLego::FinishRun()
162 // Store histograms in current Root file
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;
177 //___________________________________________
178 void AliLego::StepManager()
180 // called from AliRun::Stepmanager from gustep.
181 // Accumulate the 3 parameters step by step
184 Float_t a,z,dens,radl,absl;
187 Float_t step = gMC->TrackStep();
189 Float_t vect[3], dir[3];
190 TLorentzVector pos, mom;
192 gMC->CurrentMaterial(a,z,dens,radl,absl);
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
211 // --- See how long we have to go
217 t = fGener->PropagateCylinder(vect,dir,fGener->RadMax(),fGener->ZMax());
220 fTotAbso += step/absl;
221 fTotRadl += step/radl;
222 fTotGcm2 += step*dens;
226 ClassImp(AliLegoGenerator)
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)
237 SetPhiRange(phimin,phimax);
238 SetThetaRange(themin,themax);
243 //___________________________________________
244 void AliLegoGenerator::Generate()
246 // Create a geantino with kinematics corresponding to the current
247 // bins in theta and phi.
251 const Int_t mpart = 0;
252 Float_t orig[3], pmom[3];
253 Float_t t, cost, sint, cosp, sinp;
255 // Prepare for next step
256 if(fThetaBin>=fNtheta-1)
257 if(fPhiBin>=fNphi-1) {
258 Warning("Generate","End of Lego Generation");
262 printf("Generating rays in phi bin:%d\n",fPhiBin);
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);
277 // --- Where to start
278 orig[0] = orig[1] = orig[2] = 0;
279 Float_t dalicz = 3000;
281 t = PropagateCylinder(orig,pmom,fRadMin,dalicz);
285 if (TMath::Abs(orig[2]) > fZMax) return;
288 Float_t polar[3]={0.,0.,0.};
290 gAlice->SetTrack(1, 0, mpart, pmom, orig, polar, 0, "LEGO ray", ntr);
294 //___________________________________________
295 Float_t AliLegoGenerator::PropagateCylinder(Float_t *x, Float_t *v, Float_t r, Float_t z)
297 // Propagate to cylinder from inside
299 Double_t hnorm, sz, t, t1, t2, t3, sr;
301 const Float_t kSmall = 1e-8;
302 const Float_t kSmall2 = kSmall*kSmall;
304 // ---> Find intesection with Z planes
308 hnorm = TMath::Sqrt(1/(d[0]*d[0]+d[1]*d[1]+d[2]*d[2]));
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
316 // ---> Intersection with cylinders
317 // Intersection point (x,y,z)
318 // (x,y,z) is on track : x=X(1)+t*D(1)
321 // (x,y,z) is on cylinder : x**2 + y**2 = R**2
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];
329 t = sz; // ---> Track parallel to the z-axis, take distance to planes
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);