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