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 **************************************************************************/
19 /***********************************************************************
20 * This code defines the reconstructed cascades visualized with EVE
22 * Ludovic Gaudichet (gaudichet@to.infn.it)
23 ************************************************************************/
27 #include <Reve/Track.h>
28 #include <Reve/MCHelixLine.hi>
30 #include <TPolyLine3D.h>
31 #include <TPolyMarker3D.h>
35 #include <Reve/ReveManager.h>
44 using namespace Alieve;
47 /***********************************************************************
51 ************************************************************************/
53 const Float_t Cascade::fgkMassPion2 = 0.13956995*0.13956995;
54 const Float_t Cascade::fgkMassKaon2 = 0.493677*0.493677;
55 const Float_t Cascade::fgkMassProton2 = 0.93827231*0.93827231;
56 const Float_t Cascade::fgkMassLambda2 = 1.115683*1.115683;
58 ClassImp(Alieve::Cascade)
86 fCasCosPointingAngle(999),
89 fPolyLinePos.SetLineColor(2); // red
90 fPolyLineNeg.SetLineColor(7); // light blue
91 fPolyLineBach.SetLineColor(4); // blue
100 Cascade::Cascade(TrackRnrStyle* rs) :
125 fCasCosPointingAngle(999),
128 fMarkerColor = fRnrStyle->fFVAtt.GetMarkerColor();
129 fPolyLineV0.SetLineColor(fMarkerColor);
130 fPolyLinePos.SetLineColor(2); // red
131 fPolyLineNeg.SetLineColor(7); // light blue
132 fPolyLineBach.SetLineColor(4); // blue
134 fMainColorPtr = &fMarkerColor;
143 for (vpPathMark_i i=fPathMarksNeg.begin(); i!=fPathMarksNeg.end(); ++i)
145 for (vpPathMark_i i=fPathMarksPos.begin(); i!=fPathMarksPos.end(); ++i)
147 for (vpPathMark_i i=fPathMarksBach.begin(); i!=fPathMarksBach.end(); ++i)
150 void Cascade::Reset(TPolyLine3D* polyLine) {
151 //polyLine->SetPolyLine(n_points);
152 polyLine->SetPolyLine(0);
156 //______________________________________________________________________
157 void Cascade::SetDecayLength(Float_t primx, Float_t primy, Float_t primz) {
160 Float_t dx = fV_decay.x-primx;
161 Float_t dy = fV_decay.y-primy;
162 Float_t dz = fV_decay.z-primz;
164 fCasDecayLength = sqrt(dx*dx+dy*dy+dz*dz);
165 // This is probably wrong but I can only do this for now
166 Float_t distNorm = fCasDecayLength/GetMomentum();
167 fV_birth.x = fV_decay.x - distNorm*GetPx();
168 fV_birth.y = fV_decay.y - distNorm*GetPy();
169 fV_birth.z = fV_decay.z - distNorm*GetPz();
171 fV_bach.x = fV_decay.x;
172 fV_bach.y = fV_decay.y;
173 fV_bach.z = fV_decay.z;
177 //______________________________________________________________________
178 void Cascade::MakeTrack(vpPathMark_t& pathMark, Reve::Vector& vtx, Reve::Vector& p,
179 Int_t charge, Float_t beta, TPolyLine3D& polyLine) {
181 TrackRnrStyle& RS((fRnrStyle != 0) ? *fRnrStyle : TrackRnrStyle::fgDefStyle);
183 Float_t px = p.x, py = p.y, pz = p.z;
191 std::vector<MCVertex> track_points;
192 Bool_t decay = kFALSE;
194 if ((TMath::Abs(vtx.z) > RS.fMaxZ) || (vtx.x*vtx.x + vtx.y*vtx.y > RS.fMaxR*RS.fMaxR))
197 if (TMath::Abs(RS.fMagField) > 1e-5) {
199 // Charged particle in magnetic field
201 Float_t a = RS.fgkB2C * RS.fMagField * charge;
203 MCHelix helix(fRnrStyle, &mc_v0, TMath::C()*beta, &track_points, a); //m->cm
204 helix.Init(TMath::Sqrt(px*px+py*py), pz);
206 if(!pathMark.empty()){
207 for(std::vector<Reve::PathMark*>::iterator i=pathMark.begin();
208 i!=pathMark.end(); ++i) {
210 Reve::PathMark* pm = *i;
212 if(RS.fFitDaughters && pm->type == Reve::PathMark::Daughter){
213 if(TMath::Abs(pm->V.z) > RS.fMaxZ
214 || TMath::Sqrt(pm->V.x*pm->V.x + pm->V.y*pm->V.y) > RS.fMaxR )
217 //printf("%s fit daughter \n", fName.Data());
218 helix.LoopToVertex(p.x, p.y, p.z, pm->V.x, pm->V.y, pm->V.z);
223 if(RS.fFitDecay && pm->type == Reve::PathMark::Decay){
225 if(TMath::Abs(pm->V.z) > RS.fMaxZ
226 || TMath::Sqrt(pm->V.x*pm->V.x + pm->V.y*pm->V.y) > RS.fMaxR )
228 helix.LoopToVertex(p.x, p.y, p.z, pm->V.x, pm->V.y, pm->V.z);
236 if(!decay || RS.fFitDecay == kFALSE){
237 helix.LoopToBounds(px,py,pz);
238 // printf("%s loop to bounds \n",fName.Data() );
243 // Neutral particle or no field
245 MCLine line(fRnrStyle, &mc_v0, TMath::C()*beta, &track_points);
247 if(!pathMark.empty()) {
248 for(std::vector<Reve::PathMark*>::iterator i=pathMark.begin();
249 i!=pathMark.end(); ++i) {
250 Reve::PathMark* pm = *i;
252 if(RS.fFitDaughters && pm->type == Reve::PathMark::Daughter){
253 if(TMath::Abs(pm->V.z) > RS.fMaxZ
254 || TMath::Sqrt(pm->V.x*pm->V.x + pm->V.y*pm->V.y) > RS.fMaxR )
256 line.GotoVertex(pm->V.x, pm->V.y, pm->V.z);
262 if(RS.fFitDecay && pm->type == Reve::PathMark::Decay){
263 if(TMath::Abs(pm->V.z) > RS.fMaxZ
264 || TMath::Sqrt(pm->V.x*pm->V.x + pm->V.y*pm->V.y) > RS.fMaxR )
266 line.GotoVertex(pm->V.x, pm->V.y, pm->V.z);
274 if(!decay || RS.fFitDecay == kFALSE)
275 line.GotoBounds(px,py,pz);
280 for(std::vector<MCVertex>::iterator i=track_points.begin();
281 i!=track_points.end(); ++i) {
282 polyLine.SetNextPoint(i->x,i->y, i->z);
287 //______________________________________________________________________
288 void Cascade::MakeV0path() {
291 mc_v0.x = (fV_neg.x+fV_pos.x)/2;
292 mc_v0.y = (fV_neg.y+fV_pos.y)/2;
293 mc_v0.z = (fV_neg.z+fV_pos.z)/2;
296 std::vector<MCVertex> track_points;
297 MCLine line(fRnrStyle, &mc_v0, TMath::C()*0.99, &track_points);
299 line.GotoVertex(fV_decay.x,fV_decay.y,fV_decay.z);
302 for(std::vector<MCVertex>::iterator i=track_points.begin();
303 i!=track_points.end(); ++i) {
304 fPolyLineV0.SetNextPoint(i->x,i->y, i->z);
309 //______________________________________________________________________
310 void Cascade::MakeCasPath() {
313 mc_v0.x = fV_birth.x;
314 mc_v0.y = fV_birth.y;
315 mc_v0.z = fV_birth.z;
318 std::vector<MCVertex> track_points;
319 MCLine line(fRnrStyle, &mc_v0, TMath::C()*0.99, &track_points);
321 line.GotoVertex(fV_decay.x,fV_decay.y,fV_decay.z);
323 Reset(&fPolyLineCas);
324 for(std::vector<MCVertex>::iterator i=track_points.begin();
325 i!=track_points.end(); ++i) {
326 fPolyLineCas.SetNextPoint(i->x,i->y, i->z);
331 //______________________________________________________________________
332 void Cascade::MakeCascade()
334 SetNextPoint(fV_neg.x, fV_neg.y, fV_neg.z);
335 SetNextPoint(fV_decay.x, fV_decay.y, fV_decay.z);
337 MakeTrack(fPathMarksNeg, fV_neg, fP_neg, -1, fBeta_neg, fPolyLineNeg);
338 MakeTrack(fPathMarksPos, fV_pos, fP_pos, 1, fBeta_pos, fPolyLinePos);
340 MakeTrack(fPathMarksBach, fV_bach, fP_bach, 1, fBeta_bach, fPolyLineBach);
342 MakeTrack(fPathMarksBach, fV_bach, fP_bach, -1, -fBeta_bach, fPolyLineBach);
350 //______________________________________________________________________
351 Float_t Cascade::GetCasAlphaArmenteros() const
353 Float_t px = GetPx(), py = GetPy(), pz = GetPz();
354 Float_t posXcas, negXcas;
357 posXcas = fP_bach.x*px + fP_bach.y*py + fP_bach.z*pz;
358 negXcas = (fP_neg.x+fP_pos.x)*px + (fP_neg.y+fP_pos.y)*py + (fP_neg.z+fP_pos.z)*pz;
360 posXcas = (fP_neg.x+fP_pos.x)*px + (fP_neg.y+fP_pos.y)*py + (fP_neg.z+fP_pos.z)*pz;
361 negXcas = fP_bach.x*px + fP_bach.y*py + fP_bach.z*pz;
364 if (posXcas + negXcas > 1.e-39)
365 return (posXcas - negXcas)/(posXcas + negXcas);
370 //______________________________________________________________________
371 Float_t Cascade::GetCasPtArmenteros() const
373 Float_t px = GetPx(), py = GetPy(), pz = GetPz();
374 Float_t p2 = px*px + py*py + pz*pz;
375 if (p2 < 1.e-39) return -999;
377 Float_t posXcas, posP2;
380 posXcas = fP_bach.x*px + fP_bach.y*py + fP_bach.z*pz;
383 posXcas = (fP_neg.x+fP_pos.x)*px + (fP_neg.y+fP_pos.y)*py + (fP_neg.z+fP_pos.z)*pz;
386 return sqrt( posP2 - posXcas*posXcas/p2 );
391 /***********************************************************************
395 ************************************************************************/
397 ClassImp(Alieve::CascadeList)
400 //______________________________________________________________________
401 CascadeList::CascadeList(TrackRnrStyle* rs) :
406 fRnrV0Daughters(kTRUE),
415 fChildClass = Cascade::Class(); // override member from base RenderElementList
421 //______________________________________________________________________
422 CascadeList::CascadeList(const Text_t* name, TrackRnrStyle* rs) :
427 fRnrV0Daughters(kTRUE),
436 fChildClass = Cascade::Class(); // override member from base RenderElementList
443 //______________________________________________________________________
444 void CascadeList::Init()
446 if (fRnrStyle== 0) fRnrStyle = new TrackRnrStyle;
448 fMin[0] = 0; fMax[0] = 5; // Xi mass
449 fMin[1] = 0; fMax[1] = 5; // Omega mass
450 fMin[2] = 0; fMax[2] = 1e5; // Index
451 fMin[3] = 0.8; fMax[3] = 1; // cosPointingAngle
452 fMin[4] = 0; fMax[4] = 5; // bachV0DCA
453 fMin[5] = 0; fMax[5] = 100; // radius
454 fMin[6] = 0; fMax[6] = 10; // Pt
455 fMin[7] = -2; fMax[7] = 2; // PseudoRapidity
456 fMin[8] = 0; fMax[8] = 10; // negPt
457 fMin[9] = -2; fMax[9] = 2; // negEta
458 fMin[10] = 0; fMax[10] = 10; // posPt
459 fMin[11] = -2; fMax[11] = 2; // posEta
460 fMin[12] = 0; fMax[12] = 10; // bachPt
461 fMin[13] = -2; fMax[13] = 2; // backEta
464 fHist[0] = new TH1F(ch,ch, 100, fMin[0], fMax[0]);
466 fHist[1] = new TH1F(ch,ch, 100, fMin[1], fMax[1]);
468 fHist[2] = new TH1F(ch,ch, 100, fMin[2], fMax[2]);
470 ch = "cosPointingAngle";
471 fHist[3] = new TH1F(ch,ch, 100, fMin[3], fMax[3]);
473 fHist[4] = new TH1F(ch,ch, 100, fMin[4], fMax[4]);
475 fHist[5] = new TH1F(ch,ch, 100, fMin[5], fMax[5]);
477 fHist[6] = new TH1F(ch,ch, 100, fMin[6], fMax[6]);
478 ch = "PseudoRapidity";
479 fHist[7] = new TH1F(ch,ch, 100, fMin[7], fMax[7]);
482 fHist[8] = new TH1F(ch,ch, 100, fMin[8], fMax[8]);
484 fHist[9] = new TH1F(ch,ch, 100, fMin[9], fMax[9]);
486 fHist[10] = new TH1F(ch,ch, 100, fMin[10], fMax[10]);
488 fHist[11] = new TH1F(ch,ch, 100, fMin[11], fMax[11]);
490 fHist[12] = new TH1F(ch,ch, 100, fMin[12], fMax[12]);
492 fHist[13] = new TH1F(ch,ch, 100, fMin[13], fMax[13]);
498 ch = "ArmenterosPodolansky";
499 fHist2D[0] = new TH2F(ch,ch, 70, fMinX[0], fMaxX[0], 70,
502 for (Int_t i=0; i<fgkNcutVar; i++) {
503 fHist[i]->GetXaxis()->SetLabelSize(0.07);
504 fHist[i]->GetYaxis()->SetLabelSize(0.07);
505 fHist[i]->SetStats(0);
507 for (Int_t i=0; i<fgkNcutVar2D; i++) {
508 fHist2D[i]->GetXaxis()->SetLabelSize(0.07);
509 fHist2D[i]->GetYaxis()->SetLabelSize(0.07);
510 fHist2D[i]->SetStats(0);
515 //______________________________________________________________________
516 void CascadeList::Paint(Option_t* option) {
520 for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
521 if((*i)->GetRnrSelf()) {
522 ((Cascade*)(*i))->PaintBachelor(option);
527 if(fRnrV0Daughters) {
528 for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
529 if((*i)->GetRnrSelf()) {
530 ((Cascade*)(*i))->PaintV0Daughters(option);
536 for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
537 if((*i)->GetRnrSelf()) {
538 ((Cascade*)(*i))->PaintV0Path(option);
544 for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
545 if((*i)->GetRnrSelf()) {
546 ((Cascade*)(*i))->Paint(option);
552 for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
553 if((*i)->GetRnrSelf()) {
554 ((Cascade*)(*i))->PaintCasPath(option);
559 } // end if(fRnrSelf)
563 //______________________________________________________________________
564 void CascadeList::SetRnrV0vtx(Bool_t rnr)
570 void CascadeList::SetRnrV0path(Bool_t rnr)
576 void CascadeList::SetRnrV0Daughters(Bool_t rnr)
578 fRnrV0Daughters = rnr;
583 void CascadeList::SetRnrCasPath(Bool_t rnr)
589 void CascadeList::SetRnrCasVtx(Bool_t rnr)
595 void CascadeList::SetRnrBachelor(Bool_t rnr)
602 //______________________________________________________________________
604 void CascadeList::MakeCascades()
606 for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
607 ((Cascade*)(*i))->MakeCascade();
612 //_________________________________________________________________________
613 void CascadeList::AdjustHist(Int_t iHist) {
615 if ((iHist<0)||(iHist>=fgkNcutVar)) return;
616 if (! fHist[iHist]) return;
618 TString name = fHist[iHist]->GetName();
619 Int_t nBin = fHist[iHist]->GetXaxis()->GetNbins();
621 fHist[iHist] = new TH1F(name.Data(), name.Data(), nBin, GetMin(iHist),
623 fHist[iHist]->GetXaxis()->SetLabelSize(0.07);
624 fHist[iHist]->GetYaxis()->SetLabelSize(0.07);
625 fHist[iHist]->SetStats(0);
629 //______________________________________________________________________
630 void CascadeList::UnFill(Cascade* cas) {
633 Int_t bin = fHist[0]->GetXaxis()->FindBin(cas->GetXiMass());
634 fHist[0]->SetBinContent( bin, fHist[0]->GetBinContent(bin)-1 );
636 bin = fHist[1]->GetXaxis()->FindBin( cas->GetOmegaMass() );
637 fHist[1]->SetBinContent( bin, fHist[1]->GetBinContent(bin)-1 );
640 bin = fHist[2]->GetXaxis()->FindBin( cas->GetESDIndex() );
641 fHist[2]->SetBinContent( bin, fHist[2]->GetBinContent(bin)-1 );
643 bin = fHist[3]->GetXaxis()->FindBin( cas->GetCasCosPointingAngle() );
644 fHist[3]->SetBinContent( bin, fHist[3]->GetBinContent(bin)-1 );
646 bin = fHist[4]->GetXaxis()->FindBin( cas->GetDCA_v0_Bach() );
647 fHist[4]->SetBinContent( bin, fHist[4]->GetBinContent(bin)-1 );
649 bin = fHist[5]->GetXaxis()->FindBin( cas->GetRadius() );
650 fHist[5]->SetBinContent( bin, fHist[5]->GetBinContent(bin)-1 );
652 bin = fHist[6]->GetXaxis()->FindBin( cas->GetPt() );
653 fHist[6]->SetBinContent( bin, fHist[6]->GetBinContent(bin)-1 );
655 bin = fHist[7]->GetXaxis()->FindBin( cas->GetPseudoRapidity() );
656 fHist[7]->SetBinContent( bin, fHist[7]->GetBinContent(bin)-1 );
659 bin = fHist[8]->GetXaxis()->FindBin( cas->GetNegPt() );
660 fHist[8]->SetBinContent( bin, fHist[8]->GetBinContent(bin)-1 );
662 bin = fHist[9]->GetXaxis()->FindBin( cas->GetNegPseudoRapidity() );
663 fHist[9]->SetBinContent( bin, fHist[9]->GetBinContent(bin)-1 );
665 bin = fHist[10]->GetXaxis()->FindBin( cas->GetPosPt() );
666 fHist[10]->SetBinContent( bin, fHist[10]->GetBinContent(bin)-1 );
668 bin = fHist[11]->GetXaxis()->FindBin( cas->GetPosPseudoRapidity() );
669 fHist[11]->SetBinContent( bin, fHist[11]->GetBinContent(bin)-1 );
671 bin = fHist[12]->GetXaxis()->FindBin( cas->GetBachPt() );
672 fHist[12]->SetBinContent( bin, fHist[12]->GetBinContent(bin)-1 );
674 bin = fHist[13]->GetXaxis()->FindBin( cas->GetBachPseudoRapidity() );
675 fHist[13]->SetBinContent( bin, fHist[13]->GetBinContent(bin)-1 );
678 bin = fHist2D[0]->GetXaxis()->FindBin( cas->GetCasAlphaArmenteros() );
679 Int_t binY = fHist2D[0]->GetYaxis()->FindBin( cas->GetCasPtArmenteros() );
680 fHist2D[0]->SetBinContent( bin, binY, fHist2D[0]->GetBinContent(bin,binY)-1 );
684 //______________________________________________________________________
685 void CascadeList::Filter(Cascade* cas) {
687 Float_t xiMass = cas->GetXiMass();
688 if ((xiMass<fMin[0])||(xiMass>fMax[0])) return;
690 Float_t omegaMass = cas->GetOmegaMass();
691 if ( (omegaMass<fMin[1])||(omegaMass>fMax[1]) ) return;
694 Float_t index = cas->GetESDIndex();
695 if ( (index<fMin[2])||(index>fMax[2]) ) return;
697 Float_t cosPointingAngle = cas->GetCasCosPointingAngle();
698 if ( (cosPointingAngle<fMin[3])||(cosPointingAngle>fMax[3]) ) return;
700 Float_t bachV0DCA = cas->GetDCA_v0_Bach();
701 if ( (bachV0DCA<fMin[4])||(bachV0DCA>fMax[4]) ) return;
703 Float_t radius = cas->GetRadius();
704 if ( (radius<fMin[5])||(radius>fMax[5]) ) return;
706 Float_t pt = cas->GetPt();
707 if ( (pt<fMin[6])||(pt>fMax[6]) ) return;
709 Float_t pseudoRapidity = cas->GetPseudoRapidity();
710 if ( (pseudoRapidity<fMin[7])||(pseudoRapidity>fMax[7]) ) return;
712 Float_t negPt = cas->GetNegPt();
713 if ( (negPt<fMin[8])||(negPt>fMax[8]) ) return;
715 Float_t negEta = cas->GetNegPseudoRapidity();
716 if ( (negEta<fMin[9])||(negEta>fMax[9]) ) return;
718 Float_t posPt = cas->GetPosPt();
719 if ( (posPt<fMin[10])||(posPt>fMax[10]) ) return;
721 Float_t posEta = cas->GetPosPseudoRapidity();
722 if ( (posEta<fMin[11])||(posEta>fMax[11]) ) return;
724 Float_t bachPt = cas->GetBachPt();
725 if ( (bachPt<fMin[12])||(bachPt>fMax[12]) ) return;
727 Float_t bachEta = cas->GetBachPseudoRapidity();
728 if ( (bachEta<fMin[13])||(bachEta>fMax[13]) ) return;
730 cas->SetRnrSelf(kTRUE);
731 fHist[0]->Fill(xiMass);
732 fHist[1]->Fill(omegaMass);
733 fHist[2]->Fill(index);
734 fHist[3]->Fill(cosPointingAngle);
735 fHist[4]->Fill(bachV0DCA);
736 fHist[5]->Fill(radius);
738 fHist[7]->Fill(pseudoRapidity);
739 fHist[8]->Fill(negPt);
740 fHist[9]->Fill(negEta);
741 fHist[10]->Fill(posPt);
742 fHist[11]->Fill(posEta);
743 fHist[12]->Fill(bachPt);
744 fHist[13]->Fill(bachEta);
746 fHist2D[0]->Fill(cas->GetCasAlphaArmenteros(), cas->GetCasPtArmenteros() );
749 //______________________________________________________________________
750 void CascadeList::FilterAll() {
752 for (Int_t i=0; i<fgkNcutVar; i++)
755 for (Int_t i=0; i<fgkNcutVar2D; i++)
759 for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
761 myCas = (Cascade*)(*i);
767 //______________________________________________________________________
768 void CascadeList::GetCasIndexRange(Int_t &imin, Int_t &imax) {
772 List_i i = fChildren.begin();
773 myCas = (Cascade*)(*i);
774 index = myCas->GetESDIndex();
778 for(; i!=fChildren.end(); ++i) {
780 myCas = (Cascade*)(*i);
781 index = myCas->GetESDIndex();
782 if (index<imin) imin = index;
783 if (index>imax) imax = index;
787 //______________________________________________________________________
788 void CascadeList::XiMassFilter(Float_t min, Float_t max) {
798 for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
800 myCas = (Cascade*)(*i);
801 val = myCas->GetXiMass();
802 wasSelected = myCas->GetRnrSelf();
803 isSelected = ( (val>=min) && (val<=max) );
808 myCas->SetRnrSelf(isSelected);
811 if (isSelected) Filter(myCas);
816 //______________________________________________________________________
817 void CascadeList::OmegaMassFilter(Float_t min, Float_t max) {
827 for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
829 myCas = (Cascade*)(*i);
830 val = myCas->GetOmegaMass();
831 wasSelected = myCas->GetRnrSelf();
832 isSelected = ( (val>=min) && (val<=max) );
837 myCas->SetRnrSelf(isSelected);
840 if (isSelected) Filter(myCas);
845 //______________________________________________________________________
846 void CascadeList::IndexFilter(Float_t min, Float_t max) {
856 for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
858 myCas = (Cascade*)(*i);
859 val = myCas->GetESDIndex();
860 wasSelected = myCas->GetRnrSelf();
861 isSelected = ( (val>=min) && (val<=max) );
866 myCas->SetRnrSelf(isSelected);
869 if (isSelected) Filter(myCas);
874 //______________________________________________________________________
875 void CascadeList::CosPointingFilter(Float_t min, Float_t max) {
885 for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
887 myCas = (Cascade*)(*i);
888 val = myCas->GetCasCosPointingAngle();
889 wasSelected = myCas->GetRnrSelf();
890 isSelected = ( (val>=min) && (val<=max) );
895 myCas->SetRnrSelf(isSelected);
898 if (isSelected) Filter(myCas);
904 //______________________________________________________________________
905 void CascadeList::BachV0DCAFilter(Float_t min, Float_t max) {
915 for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
917 myCas = (Cascade*)(*i);
918 val = myCas->GetDCA_v0_Bach();
919 wasSelected = myCas->GetRnrSelf();
920 isSelected = ( (val>=min) && (val<=max) );
925 myCas->SetRnrSelf(isSelected);
928 if (isSelected) Filter(myCas);
934 //______________________________________________________________________
935 void CascadeList::RadiusFilter(Float_t min, Float_t max) {
945 for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
947 myCas = (Cascade*)(*i);
948 val = myCas->GetRadius();
949 wasSelected = myCas->GetRnrSelf();
950 isSelected = ( (val>=min) && (val<=max) );
955 myCas->SetRnrSelf(isSelected);
958 if (isSelected) Filter(myCas);
964 //______________________________________________________________________
965 void CascadeList::PtFilter(Float_t min, Float_t max) {
975 for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
977 myCas = (Cascade*)(*i);
978 val = myCas->GetPt();
979 wasSelected = myCas->GetRnrSelf();
980 isSelected = ( (val>=min) && (val<=max) );
985 myCas->SetRnrSelf(isSelected);
988 if (isSelected) Filter(myCas);
994 //______________________________________________________________________
995 void CascadeList::PseudoRapFilter(Float_t min, Float_t max) {
1005 for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
1007 myCas = (Cascade*)(*i);
1008 val = myCas->GetPseudoRapidity();
1009 wasSelected = myCas->GetRnrSelf();
1010 isSelected = ( (val>=min) && (val<=max) );
1015 myCas->SetRnrSelf(isSelected);
1018 if (isSelected) Filter(myCas);
1024 //______________________________________________________________________
1025 void CascadeList::NegPtFilter(Float_t min, Float_t max) {
1035 for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
1037 myCas = (Cascade*)(*i);
1038 val = myCas->GetNegPt();
1039 wasSelected = myCas->GetRnrSelf();
1040 isSelected = ( (val>=min) && (val<=max) );
1045 myCas->SetRnrSelf(isSelected);
1048 if (isSelected) Filter(myCas);
1054 //______________________________________________________________________
1055 void CascadeList::NegEtaFilter(Float_t min, Float_t max) {
1065 for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
1067 myCas = (Cascade*)(*i);
1068 val = myCas->GetNegPseudoRapidity();
1069 wasSelected = myCas->GetRnrSelf();
1070 isSelected = ( (val>=min) && (val<=max) );
1075 myCas->SetRnrSelf(isSelected);
1078 if (isSelected) Filter(myCas);
1084 //______________________________________________________________________
1085 void CascadeList::PosPtFilter(Float_t min, Float_t max) {
1095 for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
1097 myCas = (Cascade*)(*i);
1098 val = myCas->GetPosPt();
1099 wasSelected = myCas->GetRnrSelf();
1100 isSelected = ( (val>=min) && (val<=max) );
1105 myCas->SetRnrSelf(isSelected);
1108 if (isSelected) Filter(myCas);
1113 //______________________________________________________________________
1114 void CascadeList::PosEtaFilter(Float_t min, Float_t max) {
1124 for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
1126 myCas = (Cascade*)(*i);
1127 val = myCas->GetPosPseudoRapidity();
1128 wasSelected = myCas->GetRnrSelf();
1129 isSelected = ( (val>=min) && (val<=max) );
1134 myCas->SetRnrSelf(isSelected);
1137 if (isSelected) Filter(myCas);
1143 //______________________________________________________________________
1144 void CascadeList::BachPtFilter(Float_t min, Float_t max) {
1154 for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
1156 myCas = (Cascade*)(*i);
1157 val = myCas->GetBachPt();
1158 wasSelected = myCas->GetRnrSelf();
1159 isSelected = ( (val>=min) && (val<=max) );
1164 myCas->SetRnrSelf(isSelected);
1167 if (isSelected) Filter(myCas);
1172 //______________________________________________________________________
1173 void CascadeList::BachEtaFilter(Float_t min, Float_t max) {
1183 for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
1185 myCas = (Cascade*)(*i);
1186 val = myCas->GetBachPseudoRapidity();
1187 wasSelected = myCas->GetRnrSelf();
1188 isSelected = ( (val>=min) && (val<=max) );
1193 myCas->SetRnrSelf(isSelected);
1196 if (isSelected) Filter(myCas);