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