]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVE/Alieve/blacklisted-classes/V0.cxx
First big commit of the mchview program and its accompanying library,
[u/mrichter/AliRoot.git] / EVE / Alieve / blacklisted-classes / V0.cxx
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 /***********************************************************************
18 *  This code defines the reconstructed v0 visualized with EVE
19 *
20 * Ludovic Gaudichet (gaudichet@to.infn.it)
21 ************************************************************************/
22
23 #include "V0.h"
24
25 #include <Reve/Track.h>
26 #include <Reve/MCHelixLine.hi>
27
28 #include <TPolyLine3D.h>
29 #include <TPolyMarker3D.h>
30 #include <TColor.h>
31
32 #include <Reve/ReveManager.h>
33 #include <TCanvas.h>
34 #include <TH1F.h>
35 #include <TH2F.h>
36
37 #include <vector>
38
39 using namespace Reve;
40 using namespace Alieve;
41
42
43 /***********************************************************************
44 *
45 *  V0 class
46 *
47 ************************************************************************/
48
49 const Float_t V0::fgkMassPion2 = 0.13956995*0.13956995;
50 const Float_t V0::fgkMassProton2 = 0.93827231*0.93827231;
51
52 ClassImp(Alieve::V0)
53
54 V0::V0() :
55   RenderElement(),
56   TPolyMarker3D(1),
57   fV_neg(),
58   fP_neg(),
59   fV_pos(),
60   fP_pos(),
61   fV_v0(),
62   fV0_birth(),
63   fBeta_neg(0),
64   fBeta_pos(0),
65   fLabel_neg(0),
66   fLabel_pos(0),
67   fPathMarksNeg(),
68   fPathMarksPos(),
69   fRnrStyle(0),
70   fPolyLineNeg(),
71   fPolyLinePos(),
72   fPolyLineV0(),
73   fESDIndex(-1),
74   fDaughterDCA(999),
75   fCosPointingAngle(999),
76   fDecayLength(0)
77 {}
78
79
80 V0::V0(Reve::RecTrack* tNeg, Reve::RecTrack* tPos,
81        Reve::RecV0* v0, TrackRnrStyle* rs) :
82   RenderElement(),
83   TPolyMarker3D(1),
84   fV_neg(v0->V_neg),
85   fP_neg(v0->P_neg ),
86   fV_pos(v0->V_pos ),
87   fP_pos(v0->P_pos),
88   fV_v0(v0->V_ca),
89   fV0_birth(v0->V0_birth),
90   fBeta_neg(tNeg->beta),
91   fBeta_pos(tPos->beta),
92   fLabel_neg(v0->d_label[0]),
93   fLabel_pos(v0->d_label[1]),
94   fPathMarksNeg(),
95   fPathMarksPos(),
96   fRnrStyle(rs),
97   fPolyLineNeg(),
98   fPolyLinePos(),
99   fPolyLineV0(),
100   fESDIndex(-1),
101   fDaughterDCA(999),
102   fCosPointingAngle(999),
103   fDecayLength(0)
104 {
105   fPolyLineV0.SetLineColor(fMarkerColor);
106   fPolyLinePos.SetLineColor(2);  // red
107   fPolyLineNeg.SetLineColor(7);  // light blue
108
109   fMainColorPtr = &fMarkerColor;
110   fMarkerStyle = 20;
111   fMarkerColor = 5;
112   fMarkerSize = 0.3;
113 }
114
115
116 V0::~V0()
117 {
118   for (vpPathMark_i i=fPathMarksNeg.begin(); i!=fPathMarksNeg.end(); ++i)
119     delete *i;
120   for (vpPathMark_i i=fPathMarksPos.begin(); i!=fPathMarksPos.end(); ++i)
121     delete *i;
122 }
123
124
125 void V0::Reset(TPolyLine3D* polyLine) {
126   //polyLine->SetPolyLine(n_points);
127   polyLine->SetPolyLine(0);
128 }
129
130 //______________________________________________________________________
131 void V0::SetDecayLength(Float_t primx, Float_t primy, Float_t primz) {
132
133   Float_t dx = fV_v0.x-primx;
134   Float_t dy = fV_v0.y-primy;
135   Float_t dz = fV_v0.z-primz;
136
137   fDecayLength = sqrt(dx*dx+dy*dy+dz*dz);
138
139   // This is probably wrong but I can only do this for now
140   Float_t distNorm = fDecayLength/GetMomentum();
141   fV0_birth.x = fV_v0.x - distNorm*GetPx();
142   fV0_birth.y = fV_v0.y - distNorm*GetPy();
143   fV0_birth.z = fV_v0.z - distNorm*GetPz();
144 }
145
146
147
148 //______________________________________________________________________
149 void V0::MakeTrack(vpPathMark_t& pathMark, Reve::Vector& vtx,  Reve::Vector& p,
150                    Int_t charge, Float_t beta, TPolyLine3D& polyLine) {
151
152   TrackRnrStyle& RS((fRnrStyle != 0) ? *fRnrStyle : TrackRnrStyle::fgDefStyle);
153
154   Float_t px = p.x, py = p.y, pz = p.z;  
155
156   MCVertex  mc_v0;
157   mc_v0.x = vtx.x;
158   mc_v0.y = vtx.y; 
159   mc_v0.z = vtx.z; 
160   mc_v0.t = 0;
161
162   std::vector<MCVertex> track_points;
163   Bool_t decay = kFALSE;
164
165   if ((TMath::Abs(vtx.z) > RS.fMaxZ) || (vtx.x*vtx.x + vtx.y*vtx.y > RS.fMaxR*RS.fMaxR)) 
166     goto make_polyline;
167   
168   if (TMath::Abs(RS.fMagField) > 1e-5) {
169
170     // Charged particle in magnetic field
171
172     Float_t a = RS.fgkB2C * RS.fMagField * charge;
173    
174     MCHelix helix(fRnrStyle, &mc_v0, TMath::C()*beta, &track_points, a); //m->cm
175     helix.Init(TMath::Sqrt(px*px+py*py), pz);
176    
177     if(!pathMark.empty()){
178       for(std::vector<Reve::PathMark*>::iterator i=pathMark.begin();
179           i!=pathMark.end(); ++i) {
180
181         Reve::PathMark* pm = *i;
182         
183         if(RS.fFitDaughters &&  pm->type == Reve::PathMark::Daughter){
184           if(TMath::Abs(pm->V.z) > RS.fMaxZ 
185              || TMath::Sqrt(pm->V.x*pm->V.x + pm->V.y*pm->V.y) > RS.fMaxR )
186             goto helix_bounds;
187
188           //printf("%s fit daughter  \n", fName.Data()); 
189           helix.LoopToVertex(p.x, p.y, p.z, pm->V.x, pm->V.y, pm->V.z);
190           p.x -=  pm->P.x;
191           p.y -=  pm->P.y;
192           p.z -=  pm->P.z;
193         }
194         if(RS.fFitDecay &&  pm->type == Reve::PathMark::Decay){
195           
196           if(TMath::Abs(pm->V.z) > RS.fMaxZ 
197              || TMath::Sqrt(pm->V.x*pm->V.x + pm->V.y*pm->V.y) > RS.fMaxR )
198             goto helix_bounds;
199           helix.LoopToVertex(p.x, p.y, p.z, pm->V.x, pm->V.y, pm->V.z);
200           decay = true;
201           break;
202         }
203       }
204     }
205   helix_bounds:
206     //go to bounds
207     if(!decay || RS.fFitDecay == kFALSE){
208       helix.LoopToBounds(px,py,pz);
209       // printf("%s loop to bounds  \n",fName.Data() );
210     }
211
212   } else {
213
214     // Neutral particle or no field
215
216     MCLine line(fRnrStyle, &mc_v0, TMath::C()*beta, &track_points);
217    
218     if(!pathMark.empty()) {
219       for(std::vector<Reve::PathMark*>::iterator i=pathMark.begin();
220           i!=pathMark.end(); ++i) {
221         Reve::PathMark* pm = *i;
222
223         if(RS.fFitDaughters &&  pm->type == Reve::PathMark::Daughter){
224           if(TMath::Abs(pm->V.z) > RS.fMaxZ 
225              || TMath::Sqrt(pm->V.x*pm->V.x + pm->V.y*pm->V.y) > RS.fMaxR )
226             goto line_bounds;
227           line.GotoVertex(pm->V.x, pm->V.y, pm->V.z);
228           p.x -=  pm->P.x;
229           p.y -=  pm->P.y;
230           p.z -=  pm->P.z;
231         }
232
233         if(RS.fFitDecay &&  pm->type == Reve::PathMark::Decay){
234           if(TMath::Abs(pm->V.z) > RS.fMaxZ 
235              || TMath::Sqrt(pm->V.x*pm->V.x + pm->V.y*pm->V.y) > RS.fMaxR )
236             goto line_bounds;
237           line.GotoVertex(pm->V.x, pm->V.y, pm->V.z);
238           decay = true;
239           break;
240         }
241       }
242     }
243
244   line_bounds:
245     if(!decay || RS.fFitDecay == kFALSE)
246       line.GotoBounds(px,py,pz);
247
248   }
249 make_polyline:
250   Reset(&polyLine);
251   for(std::vector<MCVertex>::iterator i=track_points.begin();
252       i!=track_points.end(); ++i) {
253     polyLine.SetNextPoint(i->x,i->y, i->z);
254   }
255
256 }
257
258
259 //______________________________________________________________________
260 void V0::MakeV0path() {
261   
262   MCVertex  mc_v0;
263   mc_v0.x = fV_v0.x;
264   mc_v0.y = fV_v0.y; 
265   mc_v0.z = fV_v0.z; 
266   mc_v0.t = 0;
267
268   std::vector<MCVertex> track_points;
269   MCLine line(fRnrStyle, &mc_v0, TMath::C()*0.99, &track_points);
270
271  line.GotoVertex(fV0_birth.x,fV0_birth.y,fV0_birth.z);
272
273   Reset(&fPolyLineV0);
274   for(std::vector<MCVertex>::iterator i=track_points.begin();
275       i!=track_points.end(); ++i) {
276     fPolyLineV0.SetNextPoint(i->x,i->y, i->z);
277   }
278 }
279
280
281 //______________________________________________________________________
282 void V0::MakeV0()
283 {
284   SetNextPoint(fV_v0.x, fV_v0.y, fV_v0.z);
285   MakeTrack(fPathMarksNeg, fV_neg, fP_neg, -1, fBeta_neg, fPolyLineNeg);
286   MakeTrack(fPathMarksPos, fV_pos, fP_pos,  1, fBeta_pos, fPolyLinePos);
287   MakeV0path();
288   //fN = fPolyLineNeg.GetN();
289 }
290
291
292 //______________________________________________________________________
293 Float_t V0::GetAlphaArmenteros() const
294 {
295   Float_t  posXv0 = fP_pos.x*GetPx() + fP_pos.y*GetPy() + fP_pos.z*GetPz();
296   Float_t  negXv0 = fP_neg.x*GetPx() + fP_neg.y*GetPy() + fP_neg.z*GetPz();
297
298   if (posXv0+negXv0 > 1.e-39)
299     return (posXv0-negXv0)/(posXv0+negXv0);
300   else return -999;
301 }
302
303 //______________________________________________________________________
304 Float_t V0::GetPtArmenteros() const
305 {
306   Float_t  posXv0 = fP_pos.x*GetPx() + fP_pos.y*GetPy() + fP_pos.z*GetPz();
307   Float_t  v0mom2  = GetP2();
308
309   if (v0mom2 > 1.e-39)
310     return  TMath::Sqrt( GetPosP2() - posXv0*posXv0/v0mom2 ) ;
311   else return -999;
312 }
313
314
315
316 /***********************************************************************
317 *
318 *  V0List class
319 *
320 ************************************************************************/
321
322 ClassImp(Alieve::V0List)
323
324 //______________________________________________________________________
325 V0List::V0List() :
326   RenderElementList(),
327   fTitle(),
328   fRnrStyle(0),
329   fRnrDaughters(kTRUE),
330   fRnrV0vtx(kTRUE),
331   fRnrV0path(kTRUE),
332   fNegColor(0),
333   fPosColor(0)
334 {
335   fChildClass = V0::Class(); // override member from base RenderElementList
336
337   for (Int_t i=0; i<fgkNcutVar; i++)
338     fHist[i] = 0;
339   for (Int_t i=0; i<fgkNcutVar2D; i++)
340     fHist2D[i] = 0;
341 }
342
343 //______________________________________________________________________
344 V0List::V0List(TrackRnrStyle* rs) :
345   RenderElementList(),
346   fTitle(),
347   fRnrStyle(rs),
348   fRnrDaughters(kTRUE),
349   fRnrV0vtx(kTRUE),
350   fRnrV0path(kTRUE),
351   fNegColor(0),
352   fPosColor(0)
353 {
354   fChildClass = V0::Class(); // override member from base RenderElementList
355
356   Init();
357 }
358
359 //______________________________________________________________________
360 V0List::V0List(const Text_t* name, TrackRnrStyle* rs) :
361   RenderElementList(),
362   fTitle(),
363   fRnrStyle(rs),
364   fRnrDaughters(kTRUE),
365   fRnrV0vtx(kTRUE),
366   fRnrV0path(kTRUE),
367   fNegColor(0),
368   fPosColor(0)
369 {
370   fChildClass = V0::Class(); // override member from base RenderElementList
371
372   Init();
373   SetName(name);
374 }
375
376 //______________________________________________________________________
377 void V0List::Init()
378 {
379   if (fRnrStyle== 0) fRnrStyle = new TrackRnrStyle;
380
381   fMin[0]  = 0;   fMax[0] = 10; // pt
382   fMin[1]  = 0;   fMax[1] = 10; // K0s mass
383   fMin[2]  = 0;   fMax[2] = 10; // lambda mass
384   fMin[3]  = 0;   fMax[3] = 10; // anti-lambda mass
385   fMin[4]  = 0;   fMax[4] = 10; // daughter DCA
386   fMin[5]  = 0.8; fMax[5] = 1;  // cos of pointing angle
387
388   fMin[6]  =  0;   fMax[6] = 200;  // radius
389   fMin[7]  = -2;   fMax[7] =   2;  // PseudoRapidity
390   fMin[8]  =  0;   fMax[8] =  10;  // NegPt
391   fMin[9]  = -2;   fMax[9] =   2;  // NegPseudoRapidity
392   fMin[10] =  0;   fMax[10] = 10;  // PosPt
393   fMin[11] = -2;   fMax[11] =  2;  // PosPseudoRapidity
394   fMin[12] =  0;   fMax[12] =  1e5;  // index, would be good to be able to adjust it
395
396   char *ch = "ptV0";
397   fHist[0] = new TH1F(ch,ch, 100, fMin[0], fMax[0]);
398   ch = "K0sMass";
399   fHist[1] = new TH1F(ch,ch, 100, fMin[1], fMax[1]);
400   ch = "LambdaMass";
401   fHist[2] = new TH1F(ch,ch, 100, fMin[2], fMax[2]);
402   ch = "AntiLambdaMass";
403   fHist[3] = new TH1F(ch,ch, 100, fMin[3], fMax[3]);
404   ch = "daughterDCA";
405   fHist[4] = new TH1F(ch,ch, 100, fMin[4], fMax[4]);
406   ch = "cosPointingAngle";
407   fHist[5] = new TH1F(ch,ch, 100, fMin[5], fMax[5]);
408   ch = "radius";
409   fHist[6] = new TH1F(ch,ch, 100, fMin[6], fMax[6]);
410   ch = "PseudoRapidity";
411   fHist[7] = new TH1F(ch,ch, 100, fMin[7], fMax[7]);
412
413   ch = "NegPt";
414   fHist[8] = new TH1F(ch,ch, 100, fMin[8], fMax[8]);
415   ch = "NegPseudoRapidity";
416   fHist[9] = new TH1F(ch,ch, 100, fMin[9], fMax[9]);
417   ch = "PosPt";
418   fHist[10] = new TH1F(ch,ch, 100, fMin[10], fMax[10]);
419   ch = "PosPseudoRapidity";
420   fHist[11] = new TH1F(ch,ch, 100, fMin[11], fMax[11]);
421   ch = "v0Index";
422   fHist[12] = new TH1F(ch,ch, 100, fMin[12], fMax[12]);
423
424
425   fMinX[0] = -1.2;
426   fMaxX[0] = 1.2;
427   fMinY[0] = 0;
428   fMaxY[0] = 0.4;
429   ch = "ArmenterosPodolansky";
430   fHist2D[0] = new TH2F(ch,ch, 70, fMinX[0], fMaxX[0], 70,
431                         fMinY[0], fMaxY[0]);
432
433   for (Int_t i=0; i<fgkNcutVar; i++) {
434     fHist[i]->GetXaxis()->SetLabelSize(0.07);
435     fHist[i]->GetYaxis()->SetLabelSize(0.07);
436     fHist[i]->SetStats(0);
437   }
438   for (Int_t i=0; i<fgkNcutVar2D; i++) {
439     fHist2D[i]->GetXaxis()->SetLabelSize(0.07);
440     fHist2D[i]->GetYaxis()->SetLabelSize(0.07);
441     fHist2D[i]->SetStats(0);
442   }
443 }
444
445 //______________________________________________________________________
446 V0List::~V0List() {
447
448   for (Int_t i=0; i<fgkNcutVar; i++)
449     if (fHist[i]) delete fHist[i];
450   for (Int_t i=0; i<fgkNcutVar2D; i++)
451     if (fHist2D[i]) delete fHist2D[i];
452
453 }
454
455 //______________________________________________________________________
456 void V0List::Paint(Option_t* option) {
457   if(fRnrSelf) {
458
459     if(fRnrV0vtx) {
460       for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
461         if((*i)->GetRnrSelf()) {
462           ((V0*)(*i))->Paint(option);
463         }
464       }
465     }
466
467     if(fRnrDaughters) {
468       for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
469         if((*i)->GetRnrSelf()) {
470           ((V0*)(*i))->PaintDaughters(option);
471         }
472       }
473     }
474
475     if(fRnrV0path) {
476       for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
477         if((*i)->GetRnrSelf()) {
478           ((V0*)(*i))->PaintPath(option);
479         }
480       }
481     }
482   }
483 }
484
485
486 //______________________________________________________________________
487
488 void V0List::SetRnrV0vtx(Bool_t rnr) {
489   fRnrV0vtx = rnr;
490   gReve->Redraw3D();
491 }
492
493 void V0List::SetRnrV0path(Bool_t rnr) {
494   fRnrV0path = rnr;
495   gReve->Redraw3D();
496 }
497
498 void V0List::SetRnrDaughters(Bool_t rnr) {
499   fRnrDaughters = rnr;
500   gReve->Redraw3D();
501 }
502
503
504 //______________________________________________________________________
505
506 void V0List::MakeV0s() {
507   for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
508     ((V0*)(*i))->MakeV0();
509   }
510   gReve->Redraw3D();
511 }
512
513
514 void V0List::MakeMarkers() {
515   gReve->Redraw3D();
516 }
517
518
519 //_________________________________________________________________________
520 void V0List::AdjustHist(Int_t iHist) {
521
522   if ((iHist<0)||(iHist>=fgkNcutVar)) return;
523   if (! fHist[iHist]) return;
524   
525   TString name = fHist[iHist]->GetName();
526   Int_t nBin = fHist[iHist]->GetXaxis()->GetNbins();
527   delete fHist[iHist];
528   fHist[iHist] = new TH1F(name.Data(), name.Data(), nBin, GetMin(iHist),
529                           GetMax(iHist));
530   fHist[iHist]->GetXaxis()->SetLabelSize(0.07);
531   fHist[iHist]->GetYaxis()->SetLabelSize(0.07);
532   fHist[iHist]->SetStats(0);
533
534 }
535
536
537 //______________________________________________________________________
538 void V0List::UnFill(V0* v0) {
539
540     Int_t bin = fHist[0]->GetXaxis()->FindBin(v0->GetPt());
541     fHist[0]->SetBinContent( bin, fHist[0]->GetBinContent(bin)-1 );
542
543     bin = fHist[1]->GetXaxis()->FindBin( v0->GetK0mass() );
544     fHist[1]->SetBinContent( bin, fHist[1]->GetBinContent(bin)-1 );
545
546     bin = fHist[2]->GetXaxis()->FindBin( v0->GetLamMass() );
547     fHist[2]->SetBinContent( bin, fHist[2]->GetBinContent(bin)-1 );
548
549     bin = fHist[3]->GetXaxis()->FindBin( v0->GetAntiLamMass() );
550     fHist[3]->SetBinContent( bin, fHist[3]->GetBinContent(bin)-1 );
551
552     bin = fHist[4]->GetXaxis()->FindBin( v0->GetDaughterDCA() );
553     fHist[4]->SetBinContent( bin, fHist[4]->GetBinContent(bin)-1 );
554
555     bin = fHist[5]->GetXaxis()->FindBin( v0->GetCosPointingAngle() );
556     fHist[5]->SetBinContent( bin, fHist[5]->GetBinContent(bin)-1 );
557
558
559     bin = fHist[6]->GetXaxis()->FindBin( v0->GetRadius() );// radius
560     fHist[6]->SetBinContent( bin, fHist[6]->GetBinContent(bin)-1 );
561
562     bin = fHist[7]->GetXaxis()->FindBin( v0->GetPseudoRapidity() ); // PseudoRapidity
563     fHist[7]->SetBinContent( bin, fHist[7]->GetBinContent(bin)-1 );
564
565     bin = fHist[8]->GetXaxis()->FindBin( v0->GetNegPt() );// NegPt
566     fHist[8]->SetBinContent( bin, fHist[8]->GetBinContent(bin)-1 );
567
568     bin = fHist[9]->GetXaxis()->FindBin( v0->GetNegPseudoRapidity() );// NegPseudoRapidity
569     fHist[9]->SetBinContent( bin, fHist[9]->GetBinContent(bin)-1 );
570
571     bin = fHist[10]->GetXaxis()->FindBin( v0->GetPosPt() ); // PosPt
572     fHist[10]->SetBinContent( bin, fHist[10]->GetBinContent(bin)-1 );
573
574     bin = fHist[11]->GetXaxis()->FindBin( v0->GetPosPseudoRapidity() ); // PosPseudoRapidity
575     fHist[11]->SetBinContent( bin, fHist[11]->GetBinContent(bin)-1 );
576
577     bin = fHist[12]->GetXaxis()->FindBin( v0->GetESDIndex() ); // ESD index
578     fHist[12]->SetBinContent( bin, fHist[12]->GetBinContent(bin)-1 );
579
580     //---
581           bin  = fHist2D[0]->GetXaxis()->FindBin( v0->GetAlphaArmenteros() );
582     Int_t binY = fHist2D[0]->GetYaxis()->FindBin( v0->GetPtArmenteros() );
583     fHist2D[0]->SetBinContent( bin, binY, fHist2D[0]->GetBinContent(bin,binY)-1 );
584 }
585
586
587 //______________________________________________________________________
588 void V0List::Filter(V0* v0) {
589
590   Float_t pt = v0->GetPt();
591   if ((pt<fMin[0])||(pt>fMax[0])) return;
592
593   Float_t k0sMass = v0->GetK0mass();
594   if ( (k0sMass<fMin[1])||(k0sMass>fMax[1]) ) return;
595
596   Float_t lamMass = v0->GetLamMass();
597   if ( (lamMass<fMin[2])||(lamMass>fMax[2]) ) return;
598
599   Float_t alamMass = v0->GetAntiLamMass();
600   if ( (alamMass<fMin[3])||(alamMass>fMax[3]) ) return;
601
602   Float_t daughtDCA = v0->GetDaughterDCA();
603   if ( (daughtDCA<fMin[4])||(daughtDCA>fMax[4]) ) return;
604
605   Float_t cosPointing = v0->GetCosPointingAngle();
606   if ( (cosPointing<fMin[5])||(cosPointing>fMax[5]) ) return;
607
608
609   Float_t radius = v0->GetRadius();
610   if ( (radius<fMin[6])||(radius>fMax[6]) ) return;
611
612   Float_t pseudoRapidity = v0->GetPseudoRapidity();
613   if ( (pseudoRapidity<fMin[7])||(pseudoRapidity>fMax[7]) ) return;
614
615   Float_t negPt = v0->GetNegPt();
616   if ( (negPt<fMin[8])||(negPt>fMax[8]) ) return;
617
618   Float_t negPseudoRapidity = v0->GetNegPseudoRapidity();
619   if ( (negPseudoRapidity<fMin[9])||(negPseudoRapidity>fMax[9]) ) return;
620
621   Float_t posPt = v0->GetPosPt();
622   if ( (posPt<fMin[10])||(posPt>fMax[10]) ) return;
623
624   Float_t posPseudoRapidity = v0->GetPosPseudoRapidity();
625   if ( (posPseudoRapidity<fMin[11])||(posPseudoRapidity>fMax[11]) ) return;
626
627    Float_t esdIndex = v0->GetESDIndex();
628    if ( (esdIndex<fMin[12])||(esdIndex>fMax[12]) ) return;
629
630    Float_t alphaArm = v0->GetAlphaArmenteros();
631 //   if ( (alphaArm<fMinX[0])||(alphaArm>fMaxX[0]) ) return;
632
633    Float_t ptArm = v0->GetPtArmenteros();
634 //   if ( (ptArm<fMinY[0])||(ptArm>fMaxY[0]) ) return;
635
636   v0->SetRnrSelf(kTRUE);
637   fHist[0]->Fill(pt);
638   fHist[1]->Fill(k0sMass);
639   fHist[2]->Fill(lamMass);
640   fHist[3]->Fill(alamMass);
641   fHist[4]->Fill(daughtDCA);
642   fHist[5]->Fill(cosPointing);
643   fHist[6]->Fill(radius);
644   fHist[7]->Fill(pseudoRapidity);
645   fHist[8]->Fill(negPt);
646   fHist[9]->Fill(negPseudoRapidity);
647   fHist[10]->Fill(posPt);
648   fHist[11]->Fill(posPseudoRapidity);
649   fHist[12]->Fill(esdIndex);
650   fHist2D[0]->Fill(alphaArm, ptArm);
651 }
652
653 //______________________________________________________________________
654 void V0List::FilterAll() {
655
656   for (Int_t i=0; i<fgkNcutVar; i++)
657     fHist[i]->Reset();
658
659   for (Int_t i=0; i<fgkNcutVar2D; i++)
660     fHist2D[i]->Reset();
661   
662   V0* myV0;
663   for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
664
665     myV0 = (V0*)(*i);
666     Filter(myV0);
667   }
668 }
669
670
671 //______________________________________________________________________
672 void V0List::GetV0IndexRange(Int_t &imin, Int_t &imax) {
673
674   Int_t index;
675   V0* myV0;
676   List_i i=fChildren.begin();
677   myV0 = (V0*)(*i);
678   index = myV0->GetESDIndex();
679   imin = index;
680   imax = index;
681
682   for(; i!=fChildren.end(); ++i) {
683
684     myV0 = (V0*)(*i);
685     index = myV0->GetESDIndex();
686     if (index<imin) imin = index;
687     if (index>imax) imax = index;
688   }
689 }
690
691
692 //______________________________________________________________________
693 void V0List::PtFilter(Float_t min, Float_t max) {
694
695   fMin[0] = min;
696   fMax[0] = max;
697
698   Float_t val;
699   Bool_t wasSelected;
700   Bool_t isSelected;
701   V0* myV0;
702
703   for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
704
705     myV0 = (V0*)(*i);
706     val = myV0->GetPt();
707     wasSelected = myV0->GetRnrSelf();
708     isSelected = ( (val>=min) && (val<=max) );
709
710     if (wasSelected) {
711       if (! isSelected) {
712         UnFill(myV0);
713         myV0->SetRnrSelf(isSelected);
714       }
715     } else {
716       if (isSelected) Filter(myV0);
717     }
718   }
719 }
720
721
722 //______________________________________________________________________
723 void V0List::K0sMFilter(Float_t min, Float_t max) {
724
725   fMin[1] = min;
726   fMax[1] = max;
727
728   Float_t val;
729   Bool_t wasSelected;
730   Bool_t isSelected;
731   V0* myV0;
732
733   for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
734
735     myV0 = (V0*)(*i);
736     val = myV0->GetK0mass();
737     wasSelected = myV0->GetRnrSelf();
738     isSelected = ( (val>=min) && (val<=max) );
739
740     if (wasSelected) {
741       if (! isSelected) {
742         UnFill(myV0);
743         myV0->SetRnrSelf(isSelected);
744       }
745     } else {
746       if (isSelected) Filter(myV0);
747     }
748   }
749 }
750
751 //______________________________________________________________________
752 void V0List::LamMFilter(Float_t min, Float_t max) {
753
754   fMin[2] = min;
755   fMax[2] = max;
756
757   Float_t val;
758   Bool_t wasSelected;
759   Bool_t isSelected;
760   V0* myV0;
761
762   for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
763
764     myV0 = (V0*)(*i);
765     val = myV0->GetLamMass();
766     wasSelected = myV0->GetRnrSelf();
767     isSelected = ( (val>=min) && (val<=max) );
768
769     if (wasSelected) {
770       if (! isSelected) {
771         UnFill(myV0);
772         myV0->SetRnrSelf(isSelected);
773       }
774     } else {
775       if (isSelected) Filter(myV0);
776     }
777   }
778 }
779
780
781
782 //______________________________________________________________________
783 void V0List::ALamMFilter(Float_t min, Float_t max) {
784
785   fMin[3] = min;
786   fMax[3] = max;
787
788   Float_t val;
789   Bool_t wasSelected;
790   Bool_t isSelected;
791   V0* myV0;
792
793   for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
794
795     myV0 = (V0*)(*i);
796     val = myV0->GetAntiLamMass();
797     wasSelected = myV0->GetRnrSelf();
798     isSelected = ( (val>=min) && (val<=max) );
799
800     if (wasSelected) {
801       if (! isSelected) {
802         UnFill(myV0);
803         myV0->SetRnrSelf(isSelected);
804       }
805     } else {
806       if (isSelected) Filter(myV0);
807     }
808   }
809 }
810
811 //______________________________________________________________________
812 void V0List::CosPointingFilter(Float_t min, Float_t max) {
813
814   fMin[5] = min;
815   fMax[5] = max;
816
817   Float_t val;
818   Bool_t wasSelected;
819   Bool_t isSelected;
820   V0* myV0;
821
822   for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
823
824     myV0 = (V0*)(*i);
825     val = myV0->GetCosPointingAngle();
826     wasSelected = myV0->GetRnrSelf();
827     isSelected = ( (val>=min) && (val<=max) );
828
829     if (wasSelected) {
830       if (! isSelected) {
831         UnFill(myV0);
832         myV0->SetRnrSelf(isSelected);
833       }
834     } else {
835       if (isSelected) Filter(myV0);
836     }
837   }
838 }
839
840 //______________________________________________________________________
841 void V0List::DaughterDCAFilter(Float_t min, Float_t max) {
842
843   fMin[4] = min;
844   fMax[4] = max;
845
846   Float_t val;
847   Bool_t wasSelected;
848   Bool_t isSelected;
849   V0* myV0;
850
851   for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
852
853     myV0 = (V0*)(*i);
854     val = myV0->GetDaughterDCA();
855     wasSelected = myV0->GetRnrSelf();
856     isSelected = ( (val>=min) && (val<=max) );
857
858     if (wasSelected) {
859       if (! isSelected) {
860         UnFill(myV0);
861         myV0->SetRnrSelf(isSelected);
862       }
863     } else {
864       if (isSelected) Filter(myV0);
865     }
866   }
867 }
868
869 //______________________________________________________________________
870 void V0List::RadiusFilter(Float_t min, Float_t max) {
871
872   fMin[6] = min;
873   fMax[6] = max;
874
875   Float_t val;
876   Bool_t wasSelected;
877   Bool_t isSelected;
878   V0* myV0;
879
880   for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
881
882     myV0 = (V0*)(*i);
883     val = myV0->GetRadius();
884     wasSelected = myV0->GetRnrSelf();
885     isSelected = ( (val>=min) && (val<=max) );
886
887     if (wasSelected) {
888       if (! isSelected) {
889         UnFill(myV0);
890         myV0->SetRnrSelf(isSelected);
891       }
892     } else {
893       if (isSelected) Filter(myV0);
894     }
895   }
896 }
897
898 //______________________________________________________________________
899 void V0List::EtaFilter(Float_t min, Float_t max) {
900
901   fMin[7] = min;
902   fMax[7] = max;
903
904   Float_t val;
905   Bool_t wasSelected;
906   Bool_t isSelected;
907   V0* myV0;
908
909   for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
910
911     myV0 = (V0*)(*i);
912     val = myV0->GetPseudoRapidity();
913     wasSelected = myV0->GetRnrSelf();
914     isSelected = ( (val>=min) && (val<=max) );
915
916     if (wasSelected) {
917       if (! isSelected) {
918         UnFill(myV0);
919         myV0->SetRnrSelf(isSelected);
920       }
921     } else {
922       if (isSelected) Filter(myV0);
923     }
924   }
925 }
926
927 //______________________________________________________________________
928 void V0List::NegPtFilter(Float_t min, Float_t max) {
929
930   fMin[8] = min;
931   fMax[8] = max;
932
933   Float_t val;
934   Bool_t wasSelected;
935   Bool_t isSelected;
936   V0* myV0;
937
938   for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
939
940     myV0 = (V0*)(*i);
941     val = myV0->GetNegPt();
942     wasSelected = myV0->GetRnrSelf();
943     isSelected = ( (val>=min) && (val<=max) );
944
945     if (wasSelected) {
946       if (! isSelected) {
947         UnFill(myV0);
948         myV0->SetRnrSelf(isSelected);
949       }
950     } else {
951       if (isSelected) Filter(myV0);
952     }
953   }
954 }
955
956 //______________________________________________________________________
957 void V0List::NegEtaFilter(Float_t min, Float_t max) {
958
959   fMin[9] = min;
960   fMax[9] = max;
961
962   Float_t val;
963   Bool_t wasSelected;
964   Bool_t isSelected;
965   V0* myV0;
966
967   for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
968
969     myV0 = (V0*)(*i);
970     val = myV0->GetNegPseudoRapidity();
971     wasSelected = myV0->GetRnrSelf();
972     isSelected = ( (val>=min) && (val<=max) );
973
974     if (wasSelected) {
975       if (! isSelected) {
976         UnFill(myV0);
977         myV0->SetRnrSelf(isSelected);
978       }
979     } else {
980       if (isSelected) Filter(myV0);
981     }
982   }
983 }
984
985 //______________________________________________________________________
986 void V0List::PosPtFilter(Float_t min, Float_t max) {
987
988   fMin[10] = min;
989   fMax[10] = max;
990
991   Float_t val;
992   Bool_t wasSelected;
993   Bool_t isSelected;
994   V0* myV0;
995
996   for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
997
998     myV0 = (V0*)(*i);
999     val = myV0->GetPosPt();
1000     wasSelected = myV0->GetRnrSelf();
1001     isSelected = ( (val>=min) && (val<=max) );
1002
1003     if (wasSelected) {
1004       if (! isSelected) {
1005         UnFill(myV0);
1006         myV0->SetRnrSelf(isSelected);
1007       }
1008     } else {
1009       if (isSelected) Filter(myV0);
1010     }
1011   }
1012 }
1013
1014 //______________________________________________________________________
1015 void V0List::PosEtaFilter(Float_t min, Float_t max) {
1016
1017   fMin[11] = min;
1018   fMax[11] = max;
1019
1020   Float_t val;
1021   Bool_t wasSelected;
1022   Bool_t isSelected;
1023   V0* myV0;
1024
1025   for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
1026
1027     myV0 = (V0*)(*i);
1028     val = myV0->GetPosPseudoRapidity();
1029     wasSelected = myV0->GetRnrSelf();
1030     isSelected = ( (val>=min) && (val<=max) );
1031
1032     if (wasSelected) {
1033       if (! isSelected) {
1034         UnFill(myV0);
1035         myV0->SetRnrSelf(isSelected);
1036       }
1037     } else {
1038       if (isSelected) Filter(myV0);
1039     }
1040   }
1041 }
1042
1043 //______________________________________________________________________
1044 void V0List::IndexFilter(Float_t min, Float_t max) {
1045
1046   fMin[12] = min;
1047   fMax[12] = max;
1048
1049   Float_t val;
1050   Bool_t wasSelected;
1051   Bool_t isSelected;
1052   V0* myV0;
1053
1054   for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
1055
1056     myV0 = (V0*)(*i);
1057     val = myV0->GetESDIndex();
1058     wasSelected = myV0->GetRnrSelf();
1059     isSelected = ( (val>=min) && (val<=max) );
1060
1061     if (wasSelected) {
1062       if (! isSelected) {
1063         UnFill(myV0);
1064         myV0->SetRnrSelf(isSelected);
1065       }
1066     } else {
1067       if (isSelected) Filter(myV0);
1068     }
1069   }
1070 }
1071