]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ALIFAST/AliFParticle.cxx
User stepping methods added (E. Futo)
[u/mrichter/AliRoot.git] / ALIFAST / AliFParticle.cxx
1
2 //////////////////////////////////////////////////////////////////////////
3 //                                                                      //
4 // AliFParticle                                                         //
5 //                                                                      //
6 //  Graphics interface to event generators particle                     //
7 //////////////////////////////////////////////////////////////////////////
8
9 #include <TROOT.h>
10 #include "TMCParticle.h"
11 #include "TClonesArray.h"
12 #include "TPolyLine3D.h"
13 #include "TCanvas.h"
14 #include <TList.h>
15 #include <TMath.h>
16
17 #include "AliFParticle.h"
18 #include "AliFMCMaker.h"
19 #include "AliFast.h"
20 #include "AliFDisplay.h"
21
22 const Int_t kRECONS = BIT(16);
23
24 ClassImp(AliFParticle)
25
26
27 //_____________________________________________________________________________
28 AliFParticle::AliFParticle(const char * name) :TNamed(name,name)
29 {
30    // Create list to support list of particles
31    fParticles = new TList();
32    fDisplay   = (AliFDisplay*)gAliFast->Display();
33 }
34
35 //_____________________________________________________________________________
36 AliFParticle::~AliFParticle()
37 {
38    if (fParticles) fParticles->Delete();
39    delete fParticles;
40 }
41
42 //_____________________________________________________________________________
43 void AliFParticle::Clear(Option_t *)
44 {
45 //    Delete graphics temporary objects
46
47    fParticles->Delete();
48
49 }
50
51 //_____________________________________________________________________________
52 void AliFParticle::Delete(Option_t *)
53 {
54 //    Dummy
55
56 }
57
58 //_____________________________________________________________________________
59 Int_t AliFParticle::DistancetoPrimitive(Int_t px, Int_t py)
60 {
61     // scan list of particles
62    TClonesArray *particles = gAliFast->MCMaker()->Fruits();
63    Int_t uid, dist;
64    TIter next(fParticles);
65    TPolyLine3D *line;
66    while((line=(TPolyLine3D*)next())) {
67       dist = line->DistancetoPrimitive(px, py);
68       if (dist < 2) {
69          uid = line->GetUniqueID();
70          fMCParticle = (TMCParticle*)particles->UncheckedAt(uid);
71          if (!fMCParticle) continue;
72          fLine = line;
73          SetName(fMCParticle->GetName());
74          gPad->SetSelected(this);
75          gPad->SetCursor(kCross);
76          return 0;
77       }
78    }
79    return 999;
80 }
81
82 //______________________________________________________________________________
83 void AliFParticle::ExecuteEvent(Int_t event, Int_t , Int_t )
84 {
85    switch (event) {
86
87    case kButton1Down:
88       gGXW->SetLineColor(-1);
89       gPad->AbsCoordinates(kTRUE);
90       fLine->SetLineColor(6);
91       fLine->SetLineWidth(6);
92       fLine->Paint(); 
93       break;
94
95    case kMouseMotion:
96       break;
97
98    case kButton1Motion:
99       break;
100
101    case kButton1Up:
102       gGXW->SetLineColor(-1);
103       fLine->SetLineColor(kYellow);
104       fLine->SetLineWidth(1);
105       fLine->Paint(); 
106       gPad->AbsCoordinates(kFALSE);
107    }
108 }
109
110 //______________________________________________________________________________
111 char *AliFParticle::GetObjectInfo(Int_t , Int_t )
112 {
113    static char info[100];
114    sprintf(info,"px=%f, py=%f, pz=%f, E=%f",
115             fMCParticle->GetPx(),
116             fMCParticle->GetPy(),
117             fMCParticle->GetPz(),
118             fMCParticle->GetEnergy());
119
120
121    return info;
122 }
123
124
125 //______________________________________________________________________________
126 TPolyLine3D *AliFParticle::HelixCurve(Float_t field, Float_t pmom, Float_t *vin)
127 {
128 //    Estimate step size in function of field.
129 //    Create a 3-D polyline with points computed with this step size
130
131    Float_t step = 10;
132    const Int_t kMAXSTEP = 2000;
133    Float_t sx[kMAXSTEP], sy[kMAXSTEP], sz[kMAXSTEP];
134    if (pmom > 0)   step = 10*pmom;
135    if (step > 100) step = 100;
136    if (step <  .5) step = .5;
137
138    Float_t vout[6];
139    Int_t i,j;
140    sx[0] = vin[0];
141    sy[0] = vin[1];
142    sz[0] = vin[2];
143    Int_t nsteps = 1;
144    Float_t rin  = Display()->Rin();
145    Float_t zout = Display()->Zout();
146    for (i=1;i<kMAXSTEP;i++) {
147       HelixStep(field,step,pmom,vin,vout);
148       if (TMath::Abs(vout[2]) > zout) break;
149       if (vout[0]*vout[0] + vout[1]*vout[1] > 1.1*rin*rin) break;
150       sx[nsteps] = vout[0];
151       sy[nsteps] = vout[1];
152       sz[nsteps] = vout[2];
153       nsteps++;
154       for (j=0;j<6;j++) vin[j] = vout[j];
155    }
156    if (nsteps < 2) return 0;
157    TPolyLine3D *line = new TPolyLine3D(nsteps,sx, sy,sz);
158    line->SetBit(kCanDelete);
159    return line;
160 }
161
162 //______________________________________________________________________________
163 void AliFParticle::HelixStep(Float_t field, Float_t step, Float_t pmom, Float_t *vin, Float_t *vout)
164 {
165 //     extrapolate track with parameters in vector vin in a constant field
166 //     oriented along Z axis (in tesla/meters).
167 //     Output in vector vout
168 //     vin[0-->6] = x,y,z,px,py,pz
169 //     translated to C++ from GEANT3 routine GHELX3
170
171    Float_t sint, sintt, tsint, cos1t, sin2;
172 //      units are tesla,centimeters,gev/c
173    const Float_t ec = 2.9979251e-3;
174    Float_t h4  = field*ec;
175    Float_t hp  = vin[2];
176    Float_t tet = -h4*step/pmom;
177    if (TMath::Abs(tet) > 0.15) {
178       sint  = TMath::Sin(tet);
179       sintt = sint/tet;
180       tsint = (tet-sint)/tet;
181       sin2  = TMath::Sin(0.5*tet);
182       cos1t = 2*sin2*sin2/tet;
183    } else {
184       tsint = tet*tet/6;
185       sintt = 1 - tsint;
186       sint  = tet*sintt;
187       cos1t = 0.5*tet;
188    }
189    Float_t f1 = step*sintt;
190    Float_t f2 = step*cos1t;
191    Float_t f3 = step*tsint*hp;
192    Float_t f4 = -tet*cos1t;
193    Float_t f5 = sint;
194    Float_t f6 = tet*cos1t*hp;
195
196    vout[0] = vin[0] + (f1*vin[3] - f2*vin[4]);
197    vout[1] = vin[1] + (f1*vin[4] + f2*vin[3]);
198    vout[2] = vin[2] + (f1*vin[5] + f3);
199
200    vout[3] = vin[3] + (f4*vin[3] - f5*vin[4]);
201    vout[4] = vin[4] + (f4*vin[4] + f5*vin[3]);
202    vout[5] = vin[5] + (f4*vin[5] + f6);
203 }
204
205 //_____________________________________________________________________________
206 void AliFParticle::Paint(Option_t *option)
207 {
208 //    Paint particles generated by AliFMCMaker
209 //    Only particles above fPTcut are drawn
210 //    Particle trajectory is computed along an helix in a constant field
211
212
213    // clean list of particles
214    fParticles->Delete();
215
216    TClonesArray *particles = gAliFast->MCMaker()->Fruits();
217    Int_t nparticles = particles->GetEntriesFast();
218    TMCParticle *part;
219    Float_t pmom, vx, vy, vz, pt;
220    Float_t vin[6];
221    Float_t field = 2; // 2 tesla
222    Int_t KF, aKF, charge;
223    for (Int_t i=0;i<nparticles;i++) {
224       part = (TMCParticle*)particles->UncheckedAt(i);     
225       if (part->GetKS() != 1) continue;
226       KF = part->GetKF();
227       aKF = TMath::Abs(KF);
228       charge = gAliFast->MCMaker()->Charge(KF);
229       pt = TMath::Sqrt(part->GetPx()*part->GetPx() + part->GetPy()*part->GetPy());
230       if (pt < gAliFast->Display()->PTcut()) continue;
231       if (charge == 0 && pt < gAliFast->Display()->PTcutEGMUNU()) continue;
232       pmom = TMath::Sqrt(part->GetPx()*part->GetPx() + part->GetPy()*part->GetPy() + part->GetPz()*part->GetPz());
233       vx = part->GetVx();
234       vy = part->GetVy();
235       vz = part->GetVz();
236       vin[0] = vx;
237       vin[1] = vy;
238       vin[2] = vz;
239       vin[3] = part->GetPx()/pmom;
240       vin[4] = part->GetPy()/pmom;
241       vin[5] = part->GetPz()/pmom;
242       TPolyLine3D *line = HelixCurve(charge*field/3, pmom, vin);
243       if (line == 0) continue;
244       fParticles->Add(line);
245       if (part->GetLineColor() == 1) {
246          if (charge == 0) part->SetLineStyle(2);
247          part->SetLineColor(kYellow);
248          if (aKF == 11) part->SetLineColor(kRed);
249          if (aKF == 13) part->SetLineColor(kMagenta);
250          if (aKF == 22) part->SetLineColor(kGreen);
251       }
252       line->SetUniqueID(i);
253       if (!part->TestBit(kRECONS)) {
254          part->SetLineColor(7);
255          part->SetLineWidth(1);
256       }
257       line->SetLineColor(part->GetLineColor());
258       line->SetLineStyle(part->GetLineStyle());
259       line->SetLineWidth(part->GetLineWidth());
260       line->Paint(option);
261    }
262 }
263
264 //______________________________________________________________________________
265 void AliFParticle::SetLineAttributes()
266 {
267 //*-*-*-*-*-*-*-*-*Invoke the DialogCanvas Line attributes*-*-*-*-*-*-*
268 //*-*              =======================================
269
270    gROOT->SetSelectedPrimitive(fMCParticle);
271    fMCParticle->SetLineAttributes();
272 }
273
274
275 //______________________________________________________________________________
276 void AliFParticle::SizeParticles() const
277 {
278    TIter next(fParticles);
279    TPolyLine3D *line;
280    while((line=(TPolyLine3D*)next())) {
281       line->Sizeof3D();
282    }
283 }
284
285
286
287
288
289
290
291