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.10 2000/02/23 16:25:22 fca
19 AliVMC and AliGeant3 classes introduced
20 ReadEuclid moved from AliRun to AliModule
22 Revision 1.9 1999/12/03 10:54:01 fca
25 Revision 1.8 1999/10/01 09:54:33 fca
26 Correct logics for Lego StepManager
28 Revision 1.7 1999/09/29 09:24:29 fca
29 Introduction of the Copyright and cvs Log
33 //////////////////////////////////////////////////////////////
34 //////////////////////////////////////////////////////////////
36 // Utility class to evaluate the material budget from
37 // a given radius to the surface of an arbitrary cylinder
38 // along radial directions from the centre:
41 // - Interaction length
44 // Geantinos are shot in the bins in the fNtheta bins in theta
45 // and fNphi bins in phi with specified rectangular limits.
46 // The statistics are accumulated per
47 // fRadMin < r < fRadMax and <0 < z < fZMax
49 // To activate this option, you can do:
50 // Root > gAlice.RunLego();
51 // or Root > .x menu.C then select menu item "RunLego"
52 // Note that when calling gAlice->RunLego, an optional list
53 // of arguments may be specified.
57 <img src="picts/alilego.gif">
61 //////////////////////////////////////////////////////////////
72 //___________________________________________
81 //___________________________________________
82 AliLego::AliLego(const char *title, Int_t ntheta, Float_t themin, Float_t themax,
83 Int_t nphi, Float_t phimin, Float_t phimax,
84 Float_t rmin, Float_t rmax, Float_t zmax)
85 : TNamed("Lego Generator",title)
87 // specify the angular limits and the size of the rectangular box
89 fGener = new AliLegoGenerator(ntheta, themin, themax,
90 nphi, phimin, phimax, rmin, rmax, zmax);
92 gAlice->ResetGenerator(fGener);
94 Float_t etamin = -TMath::Log(TMath::Tan(TMath::Min((Double_t)themax*kDegrad/2,TMath::Pi()/2-1.e-10)));
95 Float_t etamax = -TMath::Log(TMath::Tan(TMath::Max((Double_t)themin*kDegrad/2, 1.e-10)));
97 fHistRadl = new TH2F("hradl","Radiation length map",
98 nphi,phimin,phimax,ntheta,themin,themax);
99 fHistAbso = new TH2F("habso","Interaction length map",
100 nphi,phimin,phimax,ntheta,themin,themax);
101 fHistGcm2 = new TH2F("hgcm2","g/cm2 length map",
102 nphi,phimin,phimax,ntheta,themin,themax);
103 fHistReta = new TH2F("hetar","Radiation length vs. eta",
104 nphi,phimin,phimax,ntheta,etamin,etamax);
108 //___________________________________________
118 //___________________________________________
121 // loop on phi,theta bins
123 Float_t thed, phid, eta;
124 for (Int_t i=0; i<=fGener->Nphi()*fGener->Ntheta(); ++i) {
125 // --- Set to 0 radiation length, absorption length and g/cm2 ---
130 gVMC->ProcessEvent();
132 thed = fGener->CurTheta()*kRaddeg;
133 phid = fGener->CurPhi()*kRaddeg;
134 eta = -TMath::Log(TMath::Tan(TMath::Max(
135 TMath::Min((Double_t)(fGener->CurTheta())/2,
136 TMath::Pi()/2-1.e-10),1.e-10)));
138 fHistRadl->Fill(phid,thed,fTotRadl);
139 fHistAbso->Fill(phid,thed,fTotAbso);
140 fHistGcm2->Fill(phid,thed,fTotGcm2);
141 fHistReta->Fill(phid,eta,fTotRadl);
142 gAlice->FinishEvent();
144 // store histograms in current Root file
151 //___________________________________________
152 void AliLego::StepManager()
154 // called from AliRun::Stepmanager from gustep.
155 // Accumulate the 3 parameters step by step
158 Float_t a,z,dens,radl,absl;
161 Float_t step = gMC->TrackStep();
163 Float_t vect[3], dir[3];
164 TLorentzVector pos, mom;
166 gMC->TrackPosition(pos);
167 gMC->TrackMomentum(mom);
168 gMC->CurrentMaterial(a,z,dens,radl,absl);
172 // --- See if we have to stop now
173 if (TMath::Abs(pos[2]) > fGener->ZMax() ||
174 pos[0]*pos[0] +pos[1]*pos[1] > fGener->RadMax()*fGener->RadMax()) {
175 if (gMC->TrackLength()) {
176 // Not the first step, add past contribution
185 // --- See how long we have to go
191 t = fGener->PropagateCylinder(vect,dir,fGener->RadMax(),fGener->ZMax());
194 fTotAbso += step/absl;
195 fTotRadl += step/radl;
196 fTotGcm2 += step*dens;
200 ClassImp(AliLegoGenerator)
202 //___________________________________________
203 AliLegoGenerator::AliLegoGenerator(Int_t ntheta, Float_t themin,
204 Float_t themax, Int_t nphi,
205 Float_t phimin, Float_t phimax,
206 Float_t rmin, Float_t rmax, Float_t zmax) :
207 AliGenerator(0), fRadMin(rmin), fRadMax(rmax), fZMax(zmax), fNtheta(ntheta),
208 fNphi(nphi), fThetaBin(ntheta), fPhiBin(-1), fCurTheta(0), fCurPhi(0)
211 SetPhiRange(phimin,phimax);
212 SetThetaRange(themin,themax);
217 //___________________________________________
218 void AliLegoGenerator::Generate()
220 // Create a geantino with kinematics corresponding to the current
221 // bins in theta and phi.
225 const Int_t mpart = 0;
226 Float_t orig[3], pmom[3];
227 Float_t t, cost, sint, cosp, sinp;
229 // Prepare for next step
230 if(fThetaBin>=fNtheta-1)
231 if(fPhiBin>=fNphi-1) {
232 Warning("Generate","End of Lego Generation");
236 printf("Generating rays in phi bin:%d\n",fPhiBin);
240 fCurTheta = (fThetaMin+(fThetaBin+0.5)*(fThetaMax-fThetaMin)/fNtheta);
241 fCurPhi = (fPhiMin+(fPhiBin+0.5)*(fPhiMax-fPhiMin)/fNphi);
242 cost = TMath::Cos(fCurTheta);
243 sint = TMath::Sin(fCurTheta);
244 cosp = TMath::Cos(fCurPhi);
245 sinp = TMath::Sin(fCurPhi);
251 // --- Where to start
252 orig[0] = orig[1] = orig[2] = 0;
253 Float_t dalicz = 3000;
255 t = PropagateCylinder(orig,pmom,fRadMin,dalicz);
259 if (TMath::Abs(orig[2]) > fZMax) return;
262 Float_t polar[3]={0.,0.,0.};
264 gAlice->SetTrack(1, 0, mpart, pmom, orig, polar, 0, "LEGO ray", ntr);
268 //___________________________________________
269 Float_t AliLegoGenerator::PropagateCylinder(Float_t *x, Float_t *v, Float_t r, Float_t z)
271 // Propagate to cylinder from inside
273 Double_t hnorm, sz, t, t1, t2, t3, sr;
275 const Float_t kSmall = 1e-8;
276 const Float_t kSmall2 = kSmall*kSmall;
278 // ---> Find intesection with Z planes
282 hnorm = TMath::Sqrt(1/(d[0]*d[0]+d[1]*d[1]+d[2]*d[2]));
286 if (d[2] > kSmall) sz = (z-x[2])/d[2];
287 else if (d[2] < -kSmall) sz = -(z+x[2])/d[2];
288 else sz = 1.e10; // ---> Direction parallel to X-Y, no intersection
290 // ---> Intersection with cylinders
291 // Intersection point (x,y,z)
292 // (x,y,z) is on track : x=X(1)+t*D(1)
295 // (x,y,z) is on cylinder : x**2 + y**2 = R**2
297 // (D(1)**2+D(2)**2)*t**2
298 // +2.*(X(1)*D(1)+X(2)*D(2))*t
299 // +X(1)**2+X(2)**2-R**2=0
300 // ---> Solve second degree equation
301 t1 = d[0]*d[0] + d[1]*d[1];
303 t = sz; // ---> Track parallel to the z-axis, take distance to planes
305 t2 = x[0]*d[0] + x[1]*d[1];
306 t3 = x[0]*x[0] + x[1]*x[1];
307 // ---> It should be positive, but there may be numerical problems
308 sr = (t2 +TMath::Sqrt(TMath::Max(t2*t2-(t3-r*r)*t1,0.)))/t1;
309 // ---> Find minimum distance between planes and cylinder
310 t = TMath::Min(sz,sr);