]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVE/EveDet/blacklisted-classes/AliEveCascade.cxx
Temporary fix to avoid xrootd thrashing
[u/mrichter/AliRoot.git] / EVE / EveDet / blacklisted-classes / AliEveCascade.cxx
1 // $Id$
2 // Main authors: Matevz Tadel & Alja Mrak-Tadel: 2006, 2007
3
4 /**************************************************************************
5  * Copyright(c) 1998-2008, ALICE Experiment at CERN, all rights reserved. *
6  * See http://aliceinfo.cern.ch/Offline/AliRoot/License.html for          *
7  * full copyright notice.                                                 *
8  **************************************************************************/
9 /**************************************************************************
10  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
11  *                                                                        *
12  * Author: The ALICE Off-line Project.                                    *
13  * Contributors are mentioned in the code where appropriate.              *
14  *                                                                        *
15  * Permission to use, copy, modify and distribute this software and its   *
16  * documentation strictly for non-commercial purposes is hereby granted   *
17  * without fee, provided that the above copyright notice appears in all   *
18  * copies and that both the copyright notice and this permission notice   *
19  * appear in the supporting documentation. The authors make no claims     *
20  * about the suitability of this software for any purpose. It is          *
21  * provided "as is" without express or implied warranty.                  *
22  **************************************************************************/
23
24 /* $Id$ */
25
26
27 /***********************************************************************
28 *  This code defines the reconstructed cascades visualized with EVE
29 *
30 * Ludovic Gaudichet (gaudichet@to.infn.it)
31 ************************************************************************/
32
33 #include "AliEveCascade.h"
34
35 #include <TEveTrack.h>
36 #include <MCHelixLine.hi>
37
38 #include <TPolyLine3D.h>
39 #include <TPolyMarker3D.h>
40 #include <TColor.h>
41
42 // Updates
43 #include <TEveManager.h>
44 #include <TCanvas.h>
45 #include <TH1F.h>
46 #include <TH2F.h>
47
48 #include <vector>
49
50
51
52 /***********************************************************************
53 *
54 *  AliEveCascade class
55 *
56 ************************************************************************/
57
58 const Float_t AliEveCascade::fgkMassPion2 = 0.13956995*0.13956995;
59 const Float_t AliEveCascade::fgkMassKaon2 = 0.493677*0.493677;
60 const Float_t AliEveCascade::fgkMassProton2 = 0.93827231*0.93827231;
61 const Float_t AliEveCascade::fgkMassLambda2 = 1.115683*1.115683;
62
63 ClassImp(AliEveCascade)
64
65
66 AliEveCascade::AliEveCascade() :
67   TEveElement(),
68   TPolyMarker3D(1),
69   fV_neg(),
70   fP_neg(),
71   fV_pos(),
72   fP_pos(),
73   fV_bach(),
74   fP_bach(),
75   fV_decay(),
76   fV_birth(),
77   fPathMarksNeg(),
78   fPathMarksPos(),
79   fPathMarksBach(),
80   fRnrStyle(0),
81   fPolyLineNeg(),
82   fPolyLinePos(),
83   fPolyLineBach(),
84   fPolyLineV0(),
85   fPolyLineCas(),
86   fBeta_neg(0),
87   fBeta_pos(0),
88   fBeta_bach(0),
89   fESDIndex(-1),
90   fDCA_v0_Bach(999),
91   fCasCosPointingAngle(999),
92   fCasDecayLength(999)
93 {
94   fPolyLinePos.SetLineColor(2);  // red
95   fPolyLineNeg.SetLineColor(7);  // light blue
96   fPolyLineBach.SetLineColor(4);  //  blue
97
98   fMarkerStyle = 20;
99   fMarkerColor = 5;
100   fMarkerSize = 0.3;
101 }
102
103
104
105 AliEveCascade::AliEveCascade(TEveTrackPropagator* rs) :
106   TEveElement(),
107   TPolyMarker3D(1),
108   fV_neg(),
109   fP_neg(),
110   fV_pos(),
111   fP_pos(),
112   fV_bach(),
113   fP_bach(),
114   fV_decay(),
115   fV_birth(),
116   fPathMarksNeg(),
117   fPathMarksPos(),
118   fPathMarksBach(),
119   fRnrStyle(rs),
120   fPolyLineNeg(),
121   fPolyLinePos(),
122   fPolyLineBach(),
123   fPolyLineV0(),
124   fPolyLineCas(),
125   fBeta_neg(0),
126   fBeta_pos(0),
127   fBeta_bach(0),
128   fESDIndex(-1),
129   fDCA_v0_Bach(999),
130   fCasCosPointingAngle(999),
131   fCasDecayLength(999)
132 {
133   fMarkerColor = fRnrStyle->fFVAtt.GetMarkerColor();
134   fPolyLineV0.SetLineColor(fMarkerColor);
135   fPolyLinePos.SetLineColor(2);  // red
136   fPolyLineNeg.SetLineColor(7);  // light blue
137   fPolyLineBach.SetLineColor(4);  //  blue
138
139   fMainColorPtr = &fMarkerColor;
140   fMarkerStyle = 20;
141   fMarkerColor = 5;
142   fMarkerSize = 0.3;
143 }
144
145
146 AliEveCascade::~AliEveCascade()
147 {
148   for (vpPathMark_i i=fPathMarksNeg.begin(); i!=fPathMarksNeg.end(); ++i)
149     delete *i;
150   for (vpPathMark_i i=fPathMarksPos.begin(); i!=fPathMarksPos.end(); ++i)
151     delete *i;
152   for (vpPathMark_i i=fPathMarksBach.begin(); i!=fPathMarksBach.end(); ++i)
153     delete *i;
154 }
155 void AliEveCascade::Reset(TPolyLine3D* polyLine) {
156   //polyLine->SetPolyLine(n_points);
157   polyLine->SetPolyLine(0);
158 }
159
160
161 //______________________________________________________________________________
162 void AliEveCascade::SetDecayLength(Float_t primx, Float_t primy, Float_t primz) {
163
164
165   Float_t dx = fV_decay.x-primx;
166   Float_t dy = fV_decay.y-primy;
167   Float_t dz = fV_decay.z-primz;
168
169   fCasDecayLength = sqrt(dx*dx+dy*dy+dz*dz);
170   // This is probably wrong but I can only do this for now
171   Float_t distNorm = fCasDecayLength/GetMomentum();
172   fV_birth.x = fV_decay.x - distNorm*GetPx();
173   fV_birth.y = fV_decay.y - distNorm*GetPy();
174   fV_birth.z = fV_decay.z - distNorm*GetPz();
175
176   fV_bach.x = fV_decay.x;
177   fV_bach.y = fV_decay.y;
178   fV_bach.z = fV_decay.z;
179 }
180
181
182 //______________________________________________________________________________
183 void AliEveCascade::MakeTrack(vpPathMark_t& pathMark, TEveVector& vtx,  TEveVector& p,
184                    Int_t charge, Float_t beta, TPolyLine3D& polyLine) {
185
186   TEveTrackPropagator& RS((fRnrStyle != 0) ? *fRnrStyle : TEveTrackPropagator::fgDefStyle);
187
188   Float_t px = p.x, py = p.y, pz = p.z;
189
190   MCVertex  mc_v0;
191   mc_v0.x = vtx.x;
192   mc_v0.y = vtx.y;
193   mc_v0.z = vtx.z;
194   mc_v0.t = 0;
195
196   std::vector<MCVertex> track_points;
197   Bool_t decay = kFALSE;
198
199   if ((TMath::Abs(vtx.z) > RS.fMaxZ) || (vtx.x*vtx.x + vtx.y*vtx.y > RS.fMaxR*RS.fMaxR))
200     goto make_polyline;
201
202   if (TMath::Abs(RS.fMagField) > 1e-5) {
203
204     // Charged particle in magnetic field
205
206     Float_t a = RS.fgkB2C * RS.fMagField * charge;
207
208     MCHelix helix(fRnrStyle, &mc_v0, TMath::C()*beta, &track_points, a); //m->cm
209     helix.Init(TMath::Sqrt(px*px+py*py), pz);
210
211     if(!pathMark.empty()){
212       for(std::vector<TEvePathMark*>::iterator i=pathMark.begin();
213           i!=pathMark.end(); ++i) {
214
215         TEvePathMark* pm = *i;
216
217         if(RS.fFitDaughters &&  pm->type == TEvePathMark::Daughter){
218           if(TMath::Abs(pm->V.z) > RS.fMaxZ
219              || TMath::Sqrt(pm->V.x*pm->V.x + pm->V.y*pm->V.y) > RS.fMaxR )
220             goto helix_bounds;
221
222           //printf("%s fit daughter  \n", fName.Data());
223           helix.LoopToVertex(p.x, p.y, p.z, pm->V.x, pm->V.y, pm->V.z);
224           p.x -=  pm->P.x;
225           p.y -=  pm->P.y;
226           p.z -=  pm->P.z;
227         }
228         if(RS.fFitDecay &&  pm->type == TEvePathMark::Decay){
229
230           if(TMath::Abs(pm->V.z) > RS.fMaxZ
231              || TMath::Sqrt(pm->V.x*pm->V.x + pm->V.y*pm->V.y) > RS.fMaxR )
232             goto helix_bounds;
233           helix.LoopToVertex(p.x, p.y, p.z, pm->V.x, pm->V.y, pm->V.z);
234           decay = true;
235           break;
236         }
237       }
238     }
239   helix_bounds:
240     //go to bounds
241     if(!decay || RS.fFitDecay == kFALSE){
242       helix.LoopToBounds(px,py,pz);
243       // printf("%s loop to bounds  \n",fName.Data() );
244     }
245
246   } else {
247
248     // Neutral particle or no field
249
250     MCLine line(fRnrStyle, &mc_v0, TMath::C()*beta, &track_points);
251
252     if(!pathMark.empty()) {
253       for(std::vector<TEvePathMark*>::iterator i=pathMark.begin();
254           i!=pathMark.end(); ++i) {
255         TEvePathMark* pm = *i;
256
257         if(RS.fFitDaughters &&  pm->type == TEvePathMark::Daughter){
258           if(TMath::Abs(pm->V.z) > RS.fMaxZ
259              || TMath::Sqrt(pm->V.x*pm->V.x + pm->V.y*pm->V.y) > RS.fMaxR )
260             goto line_bounds;
261           line.GotoVertex(pm->V.x, pm->V.y, pm->V.z);
262           p.x -=  pm->P.x;
263           p.y -=  pm->P.y;
264           p.z -=  pm->P.z;
265         }
266
267         if(RS.fFitDecay &&  pm->type == TEvePathMark::Decay){
268           if(TMath::Abs(pm->V.z) > RS.fMaxZ
269              || TMath::Sqrt(pm->V.x*pm->V.x + pm->V.y*pm->V.y) > RS.fMaxR )
270             goto line_bounds;
271           line.GotoVertex(pm->V.x, pm->V.y, pm->V.z);
272           decay = true;
273           break;
274         }
275       }
276     }
277
278   line_bounds:
279     if(!decay || RS.fFitDecay == kFALSE)
280       line.GotoBounds(px,py,pz);
281
282   }
283 make_polyline:
284   Reset(&polyLine);
285   for(std::vector<MCVertex>::iterator i=track_points.begin();
286       i!=track_points.end(); ++i) {
287     polyLine.SetNextPoint(i->x,i->y, i->z);
288   }
289
290 }
291
292 //______________________________________________________________________________
293 void AliEveCascade::MakeV0path() {
294
295   MCVertex  mc_v0;
296   mc_v0.x = (fV_neg.x+fV_pos.x)/2;
297   mc_v0.y = (fV_neg.y+fV_pos.y)/2;
298   mc_v0.z = (fV_neg.z+fV_pos.z)/2;
299   mc_v0.t = 0;
300
301   std::vector<MCVertex> track_points;
302   MCLine line(fRnrStyle, &mc_v0, TMath::C()*0.99, &track_points);
303
304  line.GotoVertex(fV_decay.x,fV_decay.y,fV_decay.z);
305
306   Reset(&fPolyLineV0);
307   for(std::vector<MCVertex>::iterator i=track_points.begin();
308       i!=track_points.end(); ++i) {
309     fPolyLineV0.SetNextPoint(i->x,i->y, i->z);
310   }
311
312 }
313
314 //______________________________________________________________________________
315 void AliEveCascade::MakeCasPath() {
316
317   MCVertex  mc_v0;
318   mc_v0.x = fV_birth.x;
319   mc_v0.y = fV_birth.y;
320   mc_v0.z = fV_birth.z;
321   mc_v0.t = 0;
322
323   std::vector<MCVertex> track_points;
324   MCLine line(fRnrStyle, &mc_v0, TMath::C()*0.99, &track_points);
325
326  line.GotoVertex(fV_decay.x,fV_decay.y,fV_decay.z);
327
328   Reset(&fPolyLineCas);
329   for(std::vector<MCVertex>::iterator i=track_points.begin();
330       i!=track_points.end(); ++i) {
331     fPolyLineCas.SetNextPoint(i->x,i->y, i->z);
332   }
333 }
334
335
336 //______________________________________________________________________________
337 void AliEveCascade::MakeCascade()
338 {
339   SetNextPoint(fV_neg.x, fV_neg.y, fV_neg.z);
340   SetNextPoint(fV_decay.x, fV_decay.y, fV_decay.z);
341
342   MakeTrack(fPathMarksNeg, fV_neg, fP_neg, -1, fBeta_neg, fPolyLineNeg);
343   MakeTrack(fPathMarksPos, fV_pos, fP_pos,  1, fBeta_pos, fPolyLinePos);
344   if (fBeta_bach>0)
345     MakeTrack(fPathMarksBach, fV_bach, fP_bach,  1, fBeta_bach, fPolyLineBach);
346   else
347     MakeTrack(fPathMarksBach, fV_bach, fP_bach, -1, -fBeta_bach, fPolyLineBach);
348   MakeV0path();
349   MakeCasPath();
350 }
351
352
353
354
355 //______________________________________________________________________________
356 Float_t AliEveCascade::GetCasAlphaArmenteros() const
357 {
358   Float_t px = GetPx(), py = GetPy(), pz = GetPz();
359   Float_t posXcas, negXcas;
360
361   if (fBeta_bach>0) {
362     posXcas = fP_bach.x*px + fP_bach.y*py + fP_bach.z*pz;
363     negXcas = (fP_neg.x+fP_pos.x)*px + (fP_neg.y+fP_pos.y)*py + (fP_neg.z+fP_pos.z)*pz;
364   } else {
365     posXcas = (fP_neg.x+fP_pos.x)*px + (fP_neg.y+fP_pos.y)*py + (fP_neg.z+fP_pos.z)*pz;
366     negXcas = fP_bach.x*px + fP_bach.y*py + fP_bach.z*pz;
367   }
368
369   if (posXcas + negXcas > 1.e-39)
370     return (posXcas - negXcas)/(posXcas + negXcas);
371   else return -999;
372 }
373
374
375 //______________________________________________________________________________
376 Float_t AliEveCascade::GetCasPtArmenteros() const
377 {
378   Float_t px = GetPx(), py = GetPy(), pz = GetPz();
379   Float_t p2 = px*px + py*py + pz*pz;
380   if (p2 < 1.e-39) return  -999;
381
382   Float_t posXcas, posP2;
383
384   if (fBeta_bach>0) {
385     posXcas = fP_bach.x*px + fP_bach.y*py + fP_bach.z*pz;
386     posP2 = GetBachP2();
387   } else {
388     posXcas = (fP_neg.x+fP_pos.x)*px + (fP_neg.y+fP_pos.y)*py + (fP_neg.z+fP_pos.z)*pz;
389     posP2 = GetV0P2();
390   }
391   return sqrt( posP2 - posXcas*posXcas/p2 );
392 }
393
394
395
396 /***********************************************************************
397 *
398 *  CascadeList class
399 *
400 ************************************************************************/
401
402 ClassImp(CascadeList)
403
404
405 //______________________________________________________________________________
406 CascadeList::CascadeList(TEveTrackPropagator* rs) :
407   TEveElementList(),
408   fTitle(),
409   fRnrStyle(rs),
410   fRnrBach(kTRUE),
411   fRnrV0Daughters(kTRUE),
412   fRnrV0vtx(kTRUE),
413   fRnrV0path(kTRUE),
414   fRnrCasVtx(kTRUE),
415   fRnrCasPath(kTRUE),
416   fNegColor(0),
417   fPosColor(0),
418   fBachColor(0)
419 {
420   fChildClass = AliEveCascade::Class(); // override member from base TEveElementList
421
422   Init();
423 }
424
425
426 //______________________________________________________________________________
427 CascadeList::CascadeList(const Text_t* name, TEveTrackPropagator* rs) :
428   TEveElementList(),
429   fTitle(),
430   fRnrStyle(rs),
431   fRnrBach(kTRUE),
432   fRnrV0Daughters(kTRUE),
433   fRnrV0vtx(kTRUE),
434   fRnrV0path(kTRUE),
435   fRnrCasVtx(kTRUE),
436   fRnrCasPath(kTRUE),
437   fNegColor(0),
438   fPosColor(0),
439   fBachColor(0)
440 {
441   fChildClass = AliEveCascade::Class(); // override member from base TEveElementList
442
443   Init();
444   SetName(name);
445 }
446
447
448 //______________________________________________________________________________
449 void CascadeList::Init()
450 {
451   if (fRnrStyle== 0) fRnrStyle = new TEveTrackPropagator;
452
453   fMin[0]  =  0;     fMax[0]  = 5; // Xi mass
454   fMin[1]  =  0;     fMax[1]  = 5; // Omega mass
455   fMin[2]  =  0;     fMax[2]  = 1e5; // Index
456   fMin[3]  =  0.8;   fMax[3]  = 1; // cosPointingAngle
457   fMin[4]  =  0;     fMax[4]  = 5; // bachV0DCA
458   fMin[5]  =  0;     fMax[5]  = 100; // radius
459   fMin[6]  =  0;     fMax[6]  = 10; // Pt
460   fMin[7]  = -2;     fMax[7]  = 2; // PseudoRapidity
461   fMin[8]  =  0;     fMax[8]  = 10; // negPt
462   fMin[9]  = -2;     fMax[9]  = 2; // negEta
463   fMin[10] =  0;    fMax[10]  = 10; // posPt
464   fMin[11] = -2;    fMax[11]  = 2; // posEta
465   fMin[12] =  0;    fMax[12]  = 10; // bachPt
466   fMin[13] = -2;    fMax[13]  = 2; // backEta
467
468   char *ch = "XiMass";
469   fHist[0] = new TH1F(ch,ch, 100, fMin[0], fMax[0]);
470   ch = "OmegaMass";
471   fHist[1] = new TH1F(ch,ch, 100, fMin[1], fMax[1]);
472   ch = "Index";
473   fHist[2] = new TH1F(ch,ch, 100, fMin[2], fMax[2]);
474
475   ch = "cosPointingAngle";
476   fHist[3] = new TH1F(ch,ch, 100, fMin[3], fMax[3]);
477   ch = "bachV0DCA";
478   fHist[4] = new TH1F(ch,ch, 100, fMin[4], fMax[4]);
479   ch = "radius";
480   fHist[5] = new TH1F(ch,ch, 100, fMin[5], fMax[5]);
481   ch = "Pt";
482   fHist[6] = new TH1F(ch,ch, 100, fMin[6], fMax[6]);
483   ch = "PseudoRapidity";
484   fHist[7] = new TH1F(ch,ch, 100, fMin[7], fMax[7]);
485
486   ch = "negPt";
487   fHist[8] = new TH1F(ch,ch, 100, fMin[8], fMax[8]);
488   ch = "negEta";
489   fHist[9] = new TH1F(ch,ch, 100, fMin[9], fMax[9]);
490   ch = "posPt";
491   fHist[10] = new TH1F(ch,ch, 100, fMin[10], fMax[10]);
492   ch = "posEta";
493   fHist[11] = new TH1F(ch,ch, 100, fMin[11], fMax[11]);
494   ch = "bachPt";
495   fHist[12] = new TH1F(ch,ch, 100, fMin[12], fMax[12]);
496   ch = "backEta";
497   fHist[13] = new TH1F(ch,ch, 100, fMin[13], fMax[13]);
498
499   fMinX[0] = -1.2;
500   fMaxX[0] = 1.2;
501   fMinY[0] = 0;
502   fMaxY[0] = 0.4;
503   ch = "ArmenterosPodolansky";
504   fHist2D[0] = new TH2F(ch,ch, 70, fMinX[0], fMaxX[0], 70,
505                         fMinY[0], fMaxY[0]);
506
507   for (Int_t i=0; i<fgkNcutVar; i++) {
508     fHist[i]->GetXaxis()->SetLabelSize(0.07);
509     fHist[i]->GetYaxis()->SetLabelSize(0.07);
510     fHist[i]->SetStats(0);
511   }
512   for (Int_t i=0; i<fgkNcutVar2D; i++) {
513     fHist2D[i]->GetXaxis()->SetLabelSize(0.07);
514     fHist2D[i]->GetYaxis()->SetLabelSize(0.07);
515     fHist2D[i]->SetStats(0);
516   }
517 }
518
519
520 //______________________________________________________________________________
521 void CascadeList::Paint(Option_t* option) {
522   if(fRnrSelf) {
523
524     if(fRnrBach) {
525       for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
526         if((*i)->GetRnrSelf()) {
527           ((AliEveCascade*)(*i))->PaintBachelor(option);
528         }
529       }
530     }
531
532     if(fRnrV0Daughters) {
533       for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
534         if((*i)->GetRnrSelf()) {
535           ((AliEveCascade*)(*i))->PaintV0Daughters(option);
536         }
537       }
538     }
539
540     if(fRnrV0path) {
541       for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
542         if((*i)->GetRnrSelf()) {
543           ((AliEveCascade*)(*i))->PaintV0Path(option);
544         }
545       }
546     }
547
548     if(fRnrCasVtx) {
549       for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
550         if((*i)->GetRnrSelf()) {
551           ((AliEveCascade*)(*i))->Paint(option);
552         }
553       }
554     }
555
556     if(fRnrCasPath) {
557       for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
558         if((*i)->GetRnrSelf()) {
559           ((AliEveCascade*)(*i))->PaintCasPath(option);
560         }
561       }
562     }
563
564   } // end if(fRnrSelf)
565 }
566
567
568 //______________________________________________________________________________
569 void CascadeList::SetRnrV0vtx(Bool_t rnr)
570 {
571   fRnrV0vtx = rnr;
572   gEve->Redraw3D();
573 }
574
575 void CascadeList::SetRnrV0path(Bool_t rnr)
576 {
577   fRnrV0path = rnr;
578   gEve->Redraw3D();
579 }
580
581 void CascadeList::SetRnrV0Daughters(Bool_t rnr)
582 {
583   fRnrV0Daughters = rnr;
584   gEve->Redraw3D();
585 }
586
587
588 void CascadeList::SetRnrCasPath(Bool_t rnr)
589 {
590   fRnrCasPath = rnr;
591   gEve->Redraw3D();
592 }
593
594 void CascadeList::SetRnrCasVtx(Bool_t rnr)
595 {
596   fRnrCasVtx = rnr;
597   gEve->Redraw3D();
598 }
599
600 void CascadeList::SetRnrBachelor(Bool_t rnr)
601 {
602   fRnrBach = rnr;
603   gEve->Redraw3D();
604 }
605
606
607 //______________________________________________________________________________
608
609 void CascadeList::MakeCascades()
610 {
611   for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
612     ((AliEveCascade*)(*i))->MakeCascade();
613   }
614   gEve->Redraw3D();
615 }
616
617 //______________________________________________________________________________
618 void CascadeList::AdjustHist(Int_t iHist) {
619
620   if ((iHist<0)||(iHist>=fgkNcutVar)) return;
621   if (! fHist[iHist]) return;
622
623   TString name = fHist[iHist]->GetName();
624   Int_t nBin = fHist[iHist]->GetXaxis()->GetNbins();
625   delete fHist[iHist];
626   fHist[iHist] = new TH1F(name.Data(), name.Data(), nBin, GetMin(iHist),
627                           GetMax(iHist));
628   fHist[iHist]->GetXaxis()->SetLabelSize(0.07);
629   fHist[iHist]->GetYaxis()->SetLabelSize(0.07);
630   fHist[iHist]->SetStats(0);
631
632 }
633
634 //______________________________________________________________________________
635 void CascadeList::UnFill(AliEveCascade* cas) {
636
637
638     Int_t bin = fHist[0]->GetXaxis()->FindBin(cas->GetXiMass());
639     fHist[0]->SetBinContent( bin, fHist[0]->GetBinContent(bin)-1 );
640
641     bin = fHist[1]->GetXaxis()->FindBin( cas->GetOmegaMass() );
642     fHist[1]->SetBinContent( bin, fHist[1]->GetBinContent(bin)-1 );
643
644
645     bin = fHist[2]->GetXaxis()->FindBin( cas->GetESDIndex() );
646     fHist[2]->SetBinContent( bin, fHist[2]->GetBinContent(bin)-1 );
647
648     bin = fHist[3]->GetXaxis()->FindBin( cas->GetCasCosPointingAngle() );
649     fHist[3]->SetBinContent( bin, fHist[3]->GetBinContent(bin)-1 );
650
651     bin = fHist[4]->GetXaxis()->FindBin( cas->GetDCA_v0_Bach() );
652     fHist[4]->SetBinContent( bin, fHist[4]->GetBinContent(bin)-1 );
653
654     bin = fHist[5]->GetXaxis()->FindBin( cas->GetRadius() );
655     fHist[5]->SetBinContent( bin, fHist[5]->GetBinContent(bin)-1 );
656
657     bin = fHist[6]->GetXaxis()->FindBin( cas->GetPt() );
658     fHist[6]->SetBinContent( bin, fHist[6]->GetBinContent(bin)-1 );
659
660     bin = fHist[7]->GetXaxis()->FindBin( cas->GetPseudoRapidity() );
661     fHist[7]->SetBinContent( bin, fHist[7]->GetBinContent(bin)-1 );
662     //---
663
664     bin = fHist[8]->GetXaxis()->FindBin( cas->GetNegPt() );
665     fHist[8]->SetBinContent( bin, fHist[8]->GetBinContent(bin)-1 );
666
667     bin = fHist[9]->GetXaxis()->FindBin( cas->GetNegPseudoRapidity() );
668     fHist[9]->SetBinContent( bin, fHist[9]->GetBinContent(bin)-1 );
669
670     bin = fHist[10]->GetXaxis()->FindBin( cas->GetPosPt() );
671     fHist[10]->SetBinContent( bin, fHist[10]->GetBinContent(bin)-1 );
672
673     bin = fHist[11]->GetXaxis()->FindBin( cas->GetPosPseudoRapidity() );
674     fHist[11]->SetBinContent( bin, fHist[11]->GetBinContent(bin)-1 );
675
676     bin = fHist[12]->GetXaxis()->FindBin( cas->GetBachPt() );
677     fHist[12]->SetBinContent( bin, fHist[12]->GetBinContent(bin)-1 );
678
679     bin = fHist[13]->GetXaxis()->FindBin( cas->GetBachPseudoRapidity() );
680     fHist[13]->SetBinContent( bin, fHist[13]->GetBinContent(bin)-1 );
681
682     //---
683           bin  = fHist2D[0]->GetXaxis()->FindBin( cas->GetCasAlphaArmenteros() );
684     Int_t binY = fHist2D[0]->GetYaxis()->FindBin( cas->GetCasPtArmenteros() );
685     fHist2D[0]->SetBinContent( bin, binY, fHist2D[0]->GetBinContent(bin,binY)-1 );
686 }
687
688
689 //______________________________________________________________________________
690 void CascadeList::Filter(AliEveCascade* cas) {
691
692   Float_t xiMass = cas->GetXiMass();
693   if ((xiMass<fMin[0])||(xiMass>fMax[0])) return;
694
695   Float_t omegaMass = cas->GetOmegaMass();
696   if ( (omegaMass<fMin[1])||(omegaMass>fMax[1]) ) return;
697
698
699   Float_t index = cas->GetESDIndex();
700   if ( (index<fMin[2])||(index>fMax[2]) ) return;
701
702   Float_t cosPointingAngle = cas->GetCasCosPointingAngle();
703   if ( (cosPointingAngle<fMin[3])||(cosPointingAngle>fMax[3]) ) return;
704
705   Float_t bachV0DCA = cas->GetDCA_v0_Bach();
706   if ( (bachV0DCA<fMin[4])||(bachV0DCA>fMax[4]) ) return;
707
708   Float_t radius = cas->GetRadius();
709   if ( (radius<fMin[5])||(radius>fMax[5]) ) return;
710
711   Float_t pt = cas->GetPt();
712   if ( (pt<fMin[6])||(pt>fMax[6]) ) return;
713
714   Float_t pseudoRapidity = cas->GetPseudoRapidity();
715   if ( (pseudoRapidity<fMin[7])||(pseudoRapidity>fMax[7]) ) return;
716
717   Float_t negPt = cas->GetNegPt();
718   if ( (negPt<fMin[8])||(negPt>fMax[8]) ) return;
719
720   Float_t negEta = cas->GetNegPseudoRapidity();
721   if ( (negEta<fMin[9])||(negEta>fMax[9]) ) return;
722
723   Float_t posPt = cas->GetPosPt();
724   if ( (posPt<fMin[10])||(posPt>fMax[10]) ) return;
725
726   Float_t posEta = cas->GetPosPseudoRapidity();
727   if ( (posEta<fMin[11])||(posEta>fMax[11]) ) return;
728
729   Float_t bachPt = cas->GetBachPt();
730   if ( (bachPt<fMin[12])||(bachPt>fMax[12]) ) return;
731
732   Float_t bachEta = cas->GetBachPseudoRapidity();
733   if ( (bachEta<fMin[13])||(bachEta>fMax[13]) ) return;
734
735   cas->SetRnrSelf(kTRUE);
736   fHist[0]->Fill(xiMass);
737   fHist[1]->Fill(omegaMass);
738   fHist[2]->Fill(index);
739   fHist[3]->Fill(cosPointingAngle);
740   fHist[4]->Fill(bachV0DCA);
741   fHist[5]->Fill(radius);
742   fHist[6]->Fill(pt);
743   fHist[7]->Fill(pseudoRapidity);
744   fHist[8]->Fill(negPt);
745   fHist[9]->Fill(negEta);
746   fHist[10]->Fill(posPt);
747   fHist[11]->Fill(posEta);
748   fHist[12]->Fill(bachPt);
749   fHist[13]->Fill(bachEta);
750
751   fHist2D[0]->Fill(cas->GetCasAlphaArmenteros(), cas->GetCasPtArmenteros() );
752 }
753
754 //______________________________________________________________________________
755 void CascadeList::FilterAll() {
756
757   for (Int_t i=0; i<fgkNcutVar; i++)
758     fHist[i]->Reset();
759
760   for (Int_t i=0; i<fgkNcutVar2D; i++)
761     fHist2D[i]->Reset();
762
763   AliEveCascade* myCas;
764   for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
765
766     myCas = (AliEveCascade*)(*i);
767     Filter(myCas);
768   }
769 }
770
771
772 //______________________________________________________________________________
773 void CascadeList::GetCasIndexRange(Int_t &imin, Int_t &imax) {
774
775   Int_t index;
776   AliEveCascade* myCas;
777   List_i i = fChildren.begin();
778   myCas = (AliEveCascade*)(*i);
779   index = myCas->GetESDIndex();
780   imin = index;
781   imax = index;
782
783   for(; i!=fChildren.end(); ++i) {
784
785     myCas = (AliEveCascade*)(*i);
786     index = myCas->GetESDIndex();
787     if (index<imin) imin = index;
788     if (index>imax) imax = index;
789   }
790 }
791
792 //______________________________________________________________________________
793 void CascadeList::XiMassFilter(Float_t min, Float_t max) {
794
795   fMin[0] = min;
796   fMax[0] = max;
797
798   Float_t val;
799   Bool_t wasSelected;
800   Bool_t isSelected;
801   AliEveCascade* myCas;
802
803   for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
804
805     myCas = (AliEveCascade*)(*i);
806     val = myCas->GetXiMass();
807     wasSelected = myCas->GetRnrSelf();
808     isSelected = ( (val>=min) && (val<=max) );
809
810     if (wasSelected) {
811       if (! isSelected) {
812         UnFill(myCas);
813         myCas->SetRnrSelf(isSelected);
814       }
815     } else {
816       if (isSelected) Filter(myCas);
817     }
818   }
819 }
820
821 //______________________________________________________________________________
822 void CascadeList::OmegaMassFilter(Float_t min, Float_t max) {
823
824   fMin[1] = min;
825   fMax[1] = max;
826
827   Float_t val;
828   Bool_t wasSelected;
829   Bool_t isSelected;
830   AliEveCascade* myCas;
831
832   for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
833
834     myCas = (AliEveCascade*)(*i);
835     val = myCas->GetOmegaMass();
836     wasSelected = myCas->GetRnrSelf();
837     isSelected = ( (val>=min) && (val<=max) );
838
839     if (wasSelected) {
840       if (! isSelected) {
841         UnFill(myCas);
842         myCas->SetRnrSelf(isSelected);
843       }
844     } else {
845       if (isSelected) Filter(myCas);
846     }
847   }
848 }
849
850 //______________________________________________________________________________
851 void CascadeList::IndexFilter(Float_t min, Float_t max) {
852
853   fMin[2] = min;
854   fMax[2] = max;
855
856   Float_t val;
857   Bool_t wasSelected;
858   Bool_t isSelected;
859   AliEveCascade* myCas;
860
861   for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
862
863     myCas = (AliEveCascade*)(*i);
864     val = myCas->GetESDIndex();
865     wasSelected = myCas->GetRnrSelf();
866     isSelected = ( (val>=min) && (val<=max) );
867
868     if (wasSelected) {
869       if (! isSelected) {
870         UnFill(myCas);
871         myCas->SetRnrSelf(isSelected);
872       }
873     } else {
874       if (isSelected) Filter(myCas);
875     }
876   }
877 }
878
879 //______________________________________________________________________________
880 void CascadeList::CosPointingFilter(Float_t min, Float_t max) {
881
882   fMin[3] = min;
883   fMax[3] = max;
884
885   Float_t val;
886   Bool_t wasSelected;
887   Bool_t isSelected;
888   AliEveCascade* myCas;
889
890   for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
891
892     myCas = (AliEveCascade*)(*i);
893     val = myCas->GetCasCosPointingAngle();
894     wasSelected = myCas->GetRnrSelf();
895     isSelected = ( (val>=min) && (val<=max) );
896
897     if (wasSelected) {
898       if (! isSelected) {
899         UnFill(myCas);
900         myCas->SetRnrSelf(isSelected);
901       }
902     } else {
903       if (isSelected) Filter(myCas);
904     }
905   }
906 }
907
908
909 //______________________________________________________________________________
910 void CascadeList::BachV0DCAFilter(Float_t min, Float_t max) {
911
912   fMin[4] = min;
913   fMax[4] = max;
914
915   Float_t val;
916   Bool_t wasSelected;
917   Bool_t isSelected;
918   AliEveCascade* myCas;
919
920   for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
921
922     myCas = (AliEveCascade*)(*i);
923     val = myCas->GetDCA_v0_Bach();
924     wasSelected = myCas->GetRnrSelf();
925     isSelected = ( (val>=min) && (val<=max) );
926
927     if (wasSelected) {
928       if (! isSelected) {
929         UnFill(myCas);
930         myCas->SetRnrSelf(isSelected);
931       }
932     } else {
933       if (isSelected) Filter(myCas);
934     }
935   }
936 }
937
938
939 //______________________________________________________________________________
940 void CascadeList::RadiusFilter(Float_t min, Float_t max) {
941
942   fMin[5] = min;
943   fMax[5] = max;
944
945   Float_t val;
946   Bool_t wasSelected;
947   Bool_t isSelected;
948   AliEveCascade* myCas;
949
950   for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
951
952     myCas = (AliEveCascade*)(*i);
953     val = myCas->GetRadius();
954     wasSelected = myCas->GetRnrSelf();
955     isSelected = ( (val>=min) && (val<=max) );
956
957     if (wasSelected) {
958       if (! isSelected) {
959         UnFill(myCas);
960         myCas->SetRnrSelf(isSelected);
961       }
962     } else {
963       if (isSelected) Filter(myCas);
964     }
965   }
966 }
967
968
969 //______________________________________________________________________________
970 void CascadeList::PtFilter(Float_t min, Float_t max) {
971
972   fMin[6] = min;
973   fMax[6] = max;
974
975   Float_t val;
976   Bool_t wasSelected;
977   Bool_t isSelected;
978   AliEveCascade* myCas;
979
980   for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
981
982     myCas = (AliEveCascade*)(*i);
983     val = myCas->GetPt();
984     wasSelected = myCas->GetRnrSelf();
985     isSelected = ( (val>=min) && (val<=max) );
986
987     if (wasSelected) {
988       if (! isSelected) {
989         UnFill(myCas);
990         myCas->SetRnrSelf(isSelected);
991       }
992     } else {
993       if (isSelected) Filter(myCas);
994     }
995   }
996 }
997
998
999 //______________________________________________________________________________
1000 void CascadeList::PseudoRapFilter(Float_t min, Float_t max) {
1001
1002   fMin[7] = min;
1003   fMax[7] = max;
1004
1005   Float_t val;
1006   Bool_t wasSelected;
1007   Bool_t isSelected;
1008   AliEveCascade* myCas;
1009
1010   for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
1011
1012     myCas = (AliEveCascade*)(*i);
1013     val = myCas->GetPseudoRapidity();
1014     wasSelected = myCas->GetRnrSelf();
1015     isSelected = ( (val>=min) && (val<=max) );
1016
1017     if (wasSelected) {
1018       if (! isSelected) {
1019         UnFill(myCas);
1020         myCas->SetRnrSelf(isSelected);
1021       }
1022     } else {
1023       if (isSelected) Filter(myCas);
1024     }
1025   }
1026 }
1027
1028
1029 //______________________________________________________________________________
1030 void CascadeList::NegPtFilter(Float_t min, Float_t max) {
1031
1032   fMin[8] = min;
1033   fMax[8] = max;
1034
1035   Float_t val;
1036   Bool_t wasSelected;
1037   Bool_t isSelected;
1038   AliEveCascade* myCas;
1039
1040   for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
1041
1042     myCas = (AliEveCascade*)(*i);
1043     val = myCas->GetNegPt();
1044     wasSelected = myCas->GetRnrSelf();
1045     isSelected = ( (val>=min) && (val<=max) );
1046
1047     if (wasSelected) {
1048       if (! isSelected) {
1049         UnFill(myCas);
1050         myCas->SetRnrSelf(isSelected);
1051       }
1052     } else {
1053       if (isSelected) Filter(myCas);
1054     }
1055   }
1056 }
1057
1058
1059 //______________________________________________________________________________
1060 void CascadeList::NegEtaFilter(Float_t min, Float_t max) {
1061
1062   fMin[9] = min;
1063   fMax[9] = max;
1064
1065   Float_t val;
1066   Bool_t wasSelected;
1067   Bool_t isSelected;
1068   AliEveCascade* myCas;
1069
1070   for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
1071
1072     myCas = (AliEveCascade*)(*i);
1073     val = myCas->GetNegPseudoRapidity();
1074     wasSelected = myCas->GetRnrSelf();
1075     isSelected = ( (val>=min) && (val<=max) );
1076
1077     if (wasSelected) {
1078       if (! isSelected) {
1079         UnFill(myCas);
1080         myCas->SetRnrSelf(isSelected);
1081       }
1082     } else {
1083       if (isSelected) Filter(myCas);
1084     }
1085   }
1086 }
1087
1088
1089 //______________________________________________________________________________
1090 void CascadeList::PosPtFilter(Float_t min, Float_t max) {
1091
1092   fMin[10] = min;
1093   fMax[10] = max;
1094
1095   Float_t val;
1096   Bool_t wasSelected;
1097   Bool_t isSelected;
1098   AliEveCascade* myCas;
1099
1100   for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
1101
1102     myCas = (AliEveCascade*)(*i);
1103     val = myCas->GetPosPt();
1104     wasSelected = myCas->GetRnrSelf();
1105     isSelected = ( (val>=min) && (val<=max) );
1106
1107     if (wasSelected) {
1108       if (! isSelected) {
1109         UnFill(myCas);
1110         myCas->SetRnrSelf(isSelected);
1111       }
1112     } else {
1113       if (isSelected) Filter(myCas);
1114     }
1115   }
1116 }
1117
1118 //______________________________________________________________________________
1119 void CascadeList::PosEtaFilter(Float_t min, Float_t max) {
1120
1121   fMin[11] = min;
1122   fMax[11] = max;
1123
1124   Float_t val;
1125   Bool_t wasSelected;
1126   Bool_t isSelected;
1127   AliEveCascade* myCas;
1128
1129   for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
1130
1131     myCas = (AliEveCascade*)(*i);
1132     val = myCas->GetPosPseudoRapidity();
1133     wasSelected = myCas->GetRnrSelf();
1134     isSelected = ( (val>=min) && (val<=max) );
1135
1136     if (wasSelected) {
1137       if (! isSelected) {
1138         UnFill(myCas);
1139         myCas->SetRnrSelf(isSelected);
1140       }
1141     } else {
1142       if (isSelected) Filter(myCas);
1143     }
1144   }
1145 }
1146
1147
1148 //______________________________________________________________________________
1149 void CascadeList::BachPtFilter(Float_t min, Float_t max) {
1150
1151   fMin[12] = min;
1152   fMax[12] = max;
1153
1154   Float_t val;
1155   Bool_t wasSelected;
1156   Bool_t isSelected;
1157   AliEveCascade* myCas;
1158
1159   for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
1160
1161     myCas = (AliEveCascade*)(*i);
1162     val = myCas->GetBachPt();
1163     wasSelected = myCas->GetRnrSelf();
1164     isSelected = ( (val>=min) && (val<=max) );
1165
1166     if (wasSelected) {
1167       if (! isSelected) {
1168         UnFill(myCas);
1169         myCas->SetRnrSelf(isSelected);
1170       }
1171     } else {
1172       if (isSelected) Filter(myCas);
1173     }
1174   }
1175 }
1176
1177 //______________________________________________________________________________
1178 void CascadeList::BachEtaFilter(Float_t min, Float_t max) {
1179
1180   fMin[13] = min;
1181   fMax[13] = max;
1182
1183   Float_t val;
1184   Bool_t wasSelected;
1185   Bool_t isSelected;
1186   AliEveCascade* myCas;
1187
1188   for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
1189
1190     myCas = (AliEveCascade*)(*i);
1191     val = myCas->GetBachPseudoRapidity();
1192     wasSelected = myCas->GetRnrSelf();
1193     isSelected = ( (val>=min) && (val<=max) );
1194
1195     if (wasSelected) {
1196       if (! isSelected) {
1197         UnFill(myCas);
1198         myCas->SetRnrSelf(isSelected);
1199       }
1200     } else {
1201       if (isSelected) Filter(myCas);
1202     }
1203   }
1204 }
1205