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.11 2000/03/22 13:42:26 fca
19 SetGenerator does not replace an existing generator, ResetGenerator does
21 Revision 1.10 2000/02/23 16:25:22 fca
22 AliVMC and AliGeant3 classes introduced
23 ReadEuclid moved from AliRun to AliModule
25 Revision 1.9 1999/12/03 10:54:01 fca
28 Revision 1.8 1999/10/01 09:54:33 fca
29 Correct logics for Lego StepManager
31 Revision 1.7 1999/09/29 09:24:29 fca
32 Introduction of the Copyright and cvs Log
36 //////////////////////////////////////////////////////////////
37 //////////////////////////////////////////////////////////////
39 // Utility class to evaluate the material budget from
40 // a given radius to the surface of an arbitrary cylinder
41 // along radial directions from the centre:
44 // - Interaction length
47 // Geantinos are shot in the bins in the fNtheta bins in theta
48 // and fNphi bins in phi with specified rectangular limits.
49 // The statistics are accumulated per
50 // fRadMin < r < fRadMax and <0 < z < fZMax
52 // To activate this option, you can do:
53 // Root > gAlice.RunLego();
54 // or Root > .x menu.C then select menu item "RunLego"
55 // Note that when calling gAlice->RunLego, an optional list
56 // of arguments may be specified.
60 <img src="picts/alilego.gif">
64 //////////////////////////////////////////////////////////////
75 //___________________________________________
84 //___________________________________________
85 AliLego::AliLego(const char *title, Int_t ntheta, Float_t themin, Float_t themax,
86 Int_t nphi, Float_t phimin, Float_t phimax,
87 Float_t rmin, Float_t rmax, Float_t zmax)
88 : TNamed("Lego Generator",title)
90 // specify the angular limits and the size of the rectangular box
92 fGener = new AliLegoGenerator(ntheta, themin, themax,
93 nphi, phimin, phimax, rmin, rmax, zmax);
95 gAlice->ResetGenerator(fGener);
97 Float_t etamin = -TMath::Log(TMath::Tan(TMath::Min((Double_t)themax*kDegrad/2,TMath::Pi()/2-1.e-10)));
98 Float_t etamax = -TMath::Log(TMath::Tan(TMath::Max((Double_t)themin*kDegrad/2, 1.e-10)));
100 fHistRadl = new TH2F("hradl","Radiation length map",
101 nphi,phimin,phimax,ntheta,themin,themax);
102 fHistAbso = new TH2F("habso","Interaction length map",
103 nphi,phimin,phimax,ntheta,themin,themax);
104 fHistGcm2 = new TH2F("hgcm2","g/cm2 length map",
105 nphi,phimin,phimax,ntheta,themin,themax);
106 fHistReta = new TH2F("hetar","Radiation length vs. eta",
107 nphi,phimin,phimax,ntheta,etamin,etamax);
111 //___________________________________________
121 //___________________________________________
124 // loop on phi,theta bins
126 Float_t thed, phid, eta;
127 for (Int_t i=0; i<=fGener->Nphi()*fGener->Ntheta(); ++i) {
128 // --- Set to 0 radiation length, absorption length and g/cm2 ---
135 thed = fGener->CurTheta()*kRaddeg;
136 phid = fGener->CurPhi()*kRaddeg;
137 eta = -TMath::Log(TMath::Tan(TMath::Max(
138 TMath::Min((Double_t)(fGener->CurTheta())/2,
139 TMath::Pi()/2-1.e-10),1.e-10)));
141 fHistRadl->Fill(phid,thed,fTotRadl);
142 fHistAbso->Fill(phid,thed,fTotAbso);
143 fHistGcm2->Fill(phid,thed,fTotGcm2);
144 fHistReta->Fill(phid,eta,fTotRadl);
145 gAlice->FinishEvent();
147 // store histograms in current Root file
154 //___________________________________________
155 void AliLego::StepManager()
157 // called from AliRun::Stepmanager from gustep.
158 // Accumulate the 3 parameters step by step
161 Float_t a,z,dens,radl,absl;
164 Float_t step = gMC->TrackStep();
166 Float_t vect[3], dir[3];
167 TLorentzVector pos, mom;
169 gMC->TrackPosition(pos);
170 gMC->TrackMomentum(mom);
171 gMC->CurrentMaterial(a,z,dens,radl,absl);
175 // --- See if we have to stop now
176 if (TMath::Abs(pos[2]) > fGener->ZMax() ||
177 pos[0]*pos[0] +pos[1]*pos[1] > fGener->RadMax()*fGener->RadMax()) {
178 if (gMC->TrackLength()) {
179 // Not the first step, add past contribution
188 // --- See how long we have to go
194 t = fGener->PropagateCylinder(vect,dir,fGener->RadMax(),fGener->ZMax());
197 fTotAbso += step/absl;
198 fTotRadl += step/radl;
199 fTotGcm2 += step*dens;
203 ClassImp(AliLegoGenerator)
205 //___________________________________________
206 AliLegoGenerator::AliLegoGenerator(Int_t ntheta, Float_t themin,
207 Float_t themax, Int_t nphi,
208 Float_t phimin, Float_t phimax,
209 Float_t rmin, Float_t rmax, Float_t zmax) :
210 AliGenerator(0), fRadMin(rmin), fRadMax(rmax), fZMax(zmax), fNtheta(ntheta),
211 fNphi(nphi), fThetaBin(ntheta), fPhiBin(-1), fCurTheta(0), fCurPhi(0)
214 SetPhiRange(phimin,phimax);
215 SetThetaRange(themin,themax);
220 //___________________________________________
221 void AliLegoGenerator::Generate()
223 // Create a geantino with kinematics corresponding to the current
224 // bins in theta and phi.
228 const Int_t mpart = 0;
229 Float_t orig[3], pmom[3];
230 Float_t t, cost, sint, cosp, sinp;
232 // Prepare for next step
233 if(fThetaBin>=fNtheta-1)
234 if(fPhiBin>=fNphi-1) {
235 Warning("Generate","End of Lego Generation");
239 printf("Generating rays in phi bin:%d\n",fPhiBin);
243 fCurTheta = (fThetaMin+(fThetaBin+0.5)*(fThetaMax-fThetaMin)/fNtheta);
244 fCurPhi = (fPhiMin+(fPhiBin+0.5)*(fPhiMax-fPhiMin)/fNphi);
245 cost = TMath::Cos(fCurTheta);
246 sint = TMath::Sin(fCurTheta);
247 cosp = TMath::Cos(fCurPhi);
248 sinp = TMath::Sin(fCurPhi);
254 // --- Where to start
255 orig[0] = orig[1] = orig[2] = 0;
256 Float_t dalicz = 3000;
258 t = PropagateCylinder(orig,pmom,fRadMin,dalicz);
262 if (TMath::Abs(orig[2]) > fZMax) return;
265 Float_t polar[3]={0.,0.,0.};
267 gAlice->SetTrack(1, 0, mpart, pmom, orig, polar, 0, "LEGO ray", ntr);
271 //___________________________________________
272 Float_t AliLegoGenerator::PropagateCylinder(Float_t *x, Float_t *v, Float_t r, Float_t z)
274 // Propagate to cylinder from inside
276 Double_t hnorm, sz, t, t1, t2, t3, sr;
278 const Float_t kSmall = 1e-8;
279 const Float_t kSmall2 = kSmall*kSmall;
281 // ---> Find intesection with Z planes
285 hnorm = TMath::Sqrt(1/(d[0]*d[0]+d[1]*d[1]+d[2]*d[2]));
289 if (d[2] > kSmall) sz = (z-x[2])/d[2];
290 else if (d[2] < -kSmall) sz = -(z+x[2])/d[2];
291 else sz = 1.e10; // ---> Direction parallel to X-Y, no intersection
293 // ---> Intersection with cylinders
294 // Intersection point (x,y,z)
295 // (x,y,z) is on track : x=X(1)+t*D(1)
298 // (x,y,z) is on cylinder : x**2 + y**2 = R**2
300 // (D(1)**2+D(2)**2)*t**2
301 // +2.*(X(1)*D(1)+X(2)*D(2))*t
302 // +X(1)**2+X(2)**2-R**2=0
303 // ---> Solve second degree equation
304 t1 = d[0]*d[0] + d[1]*d[1];
306 t = sz; // ---> Track parallel to the z-axis, take distance to planes
308 t2 = x[0]*d[0] + x[1]*d[1];
309 t3 = x[0]*x[0] + x[1]*x[1];
310 // ---> It should be positive, but there may be numerical problems
311 sr = (t2 +TMath::Sqrt(TMath::Max(t2*t2-(t3-r*r)*t1,0.)))/t1;
312 // ---> Find minimum distance between planes and cylinder
313 t = TMath::Min(sz,sr);