]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVE/Reve/CascadeEditors.cxx
Two new items for context menu: PrintPathMarks and ImportDaughters.
[u/mrichter/AliRoot.git] / EVE / Reve / CascadeEditors.cxx
1
2 /***********************************************************************
3   This editor appears in the Reve window when v0 are visualize.
4 It allows to select the v0 as a function of some useful parameters.
5
6 Ludovic Gaudichet (gaudichet@to.infn.it)
7 ************************************************************************/
8
9
10
11
12 #include "CascadeEditors.h"
13 #include <Reve/Cascade.h>
14
15 #include <Reve/RGValuators.h>
16
17 #include <TVirtualPad.h>
18 #include <TColor.h>
19
20 #include <TGLabel.h>
21 #include <TGButton.h>
22 #include <TGNumberEntry.h>
23 #include <TGColorSelect.h>
24 #include <TGDoubleSlider.h>
25
26 #include <TGTab.h>
27 #include <TRootEmbeddedCanvas.h>
28 #include <TCanvas.h>
29 #include <TH1.h>
30 #include <TH1F.h>
31 #include <TH2F.h>
32
33
34 using namespace Reve;
35
36 //______________________________________________________________________
37 // CascadeListEditor
38 //
39
40 ClassImp(CascadeListEditor)
41
42 CascadeListEditor::CascadeListEditor(const TGWindow *p,
43                                  Int_t width, Int_t height,
44                                  UInt_t options, Pixel_t back) :
45   TGedFrame(p, width, height, options | kVerticalFrame, back),
46   fMList(0),
47   fRnrV0Daughters(0),
48   fRnrV0path(0),
49   fRnrVtx(0),
50   fRnrBach(0),
51   fRnrCasPath(0),
52   fMainTabA(0),
53   fMainTabB(0)
54 {
55   MakeTitle("CascadeList");
56  
57   //TGHorizontalFrame* frame = new TGHorizontalFrame(this);
58  
59   // --- Rendering control
60
61   fRnrVtx = new TGCheckButton(this, "Render v0 and cascade vertices");
62   AddFrame(fRnrVtx, new TGLayoutHints(kLHintsTop, 3, 1, 1, 0));
63   fRnrVtx->Connect("Toggled(Bool_t)",
64                      "Reve::CascadeListEditor", this, "DoRnrVtx()");
65
66   fRnrV0path = new TGCheckButton(this, "Render v0 path");
67   AddFrame(fRnrV0path, new TGLayoutHints(kLHintsTop, 3, 1, 1, 0));
68   fRnrV0path->Connect("Toggled(Bool_t)",
69                      "Reve::CascadeListEditor", this, "DoRnrV0path()");  
70
71   fRnrCasPath = new TGCheckButton(this, "Render cascade path");
72   AddFrame(fRnrCasPath, new TGLayoutHints(kLHintsTop, 3, 1, 1, 0));
73   fRnrCasPath->Connect("Toggled(Bool_t)",
74                        "Reve::CascadeListEditor", this, "DoRnrCasPath()");  
75
76   fRnrV0Daughters = new TGCheckButton(this, "Render v0 daughter tracks");
77   AddFrame(fRnrV0Daughters, new TGLayoutHints(kLHintsTop, 3, 1, 1, 0));
78   fRnrV0Daughters->Connect("Toggled(Bool_t)",
79                         "Reve::CascadeListEditor", this, "DoRnrV0Daughters()");
80
81   fRnrBach = new TGCheckButton(this, "Render bachelor track");
82   AddFrame(fRnrBach, new TGLayoutHints(kLHintsTop, 3, 1, 1, 0));
83   fRnrBach->Connect("Toggled(Bool_t)",
84                         "Reve::CascadeListEditor", this, "DoRnrBach()");
85     
86   for (Int_t i=0; i<fgkNRange; i++) fRange[i]=0;
87   for (Int_t i=0; i<fgkNCanvas; i++) fCanvasA[i]=0;
88   for (Int_t i=0; i<fgkNCanvas; i++) fCanvasB[i]=0;
89
90   AddSelectTab();
91   AddSeeTab();
92
93   TGTextButton* resetCutsButton = new TGTextButton(this, "Reset all cuts", 40);
94   resetCutsButton->Connect("Clicked()", "Reve::CascadeListEditor", this, "ResetCuts()");
95   AddFrame(resetCutsButton, new TGLayoutHints(kLHintsTop, 3, 1, 1, 0));
96 }
97
98 CascadeListEditor::~CascadeListEditor()
99 {}
100
101
102 //_________________________________________________________________________
103 void CascadeListEditor::SetModel(TObject* obj)
104 {
105   fMList = dynamic_cast<CascadeList*>(obj);
106
107   for (Int_t i=0; i<fgkNRange; i++)
108     if (fRange[i])
109       fRange[i]->SetValues( fMList->GetMin(i), fMList->GetMax(i) );
110
111   fRnrV0Daughters->SetState(fMList->GetRnrV0Daughters() ? kButtonDown : kButtonUp);
112   fRnrV0path->SetState(fMList->GetRnrV0path() ? kButtonDown : kButtonUp);
113   fRnrVtx->SetState(fMList->GetRnrCasVtx() ? kButtonDown : kButtonUp);
114   fRnrBach->SetState(fMList->GetRnrBachelor() ? kButtonDown : kButtonUp);
115   fRnrCasPath->SetState(fMList->GetRnrCasPath() ? kButtonDown : kButtonUp);
116
117   FillCanvas();
118 }
119
120
121 //_________________________________________________________________________
122 TGCompositeFrame* CascadeListEditor::AddTab(TGTab* tab, Int_t i, Int_t can,
123                                        char *name) {
124
125   TGCompositeFrame* frameTab = tab->AddTab(name);
126   frameTab->SetLayoutManager(new TGVerticalLayout( frameTab ));
127
128   TRootEmbeddedCanvas** embeddedCanvas = 0;
129
130   if (can==1) embeddedCanvas = &fCanvasA[i];
131   else if (can==2) embeddedCanvas = &fCanvasB[i];
132
133   *embeddedCanvas = new TRootEmbeddedCanvas("EmbeddedCanvas", frameTab, 200, 200);
134   TCanvas *c1 = (*embeddedCanvas)->GetCanvas();
135   c1->SetBorderMode(0);
136   c1->SetTicks(1,0);
137   c1->SetGrid();
138   c1->SetBorderMode(0);
139  
140   frameTab->AddFrame(*embeddedCanvas, new TGLayoutHints(kLHintsTop|kLHintsExpandX,
141                                                        0, 0, 0, 0));
142   return frameTab;
143 }
144
145
146 //_________________________________________________________________________
147 TGCompositeFrame** CascadeListEditor::CreateTab(TGTab **pMainTab, TGTab **ptab, Int_t can) {
148
149   //------
150   // tab widget
151   pMainTab[0] = new TGTab(this,0,0);
152   pMainTab[0]->Connect("Selected(Int_t)", "Reve::CascadeListEditor", this, "UpdateSelectedTab()");
153   this->AddFrame(pMainTab[0], new TGLayoutHints( kLHintsTop | kLHintsExpandX,2,2,2,2));
154
155   //------
156   // container of "Tab1"
157   TGCompositeFrame *frameTab1 = pMainTab[0]->AddTab("ident.");
158   frameTab1->SetLayoutManager(new TGVerticalLayout(frameTab1));
159   
160   // tab widget
161   ptab[0] = new TGTab(frameTab1,2,2);
162   ptab[0]->Resize(ptab[0]->GetDefaultSize());
163   // The following is for updating the canvas of a tab if this one is selected
164   // (it updates every canvas)
165   ptab[0]->Connect("Selected(Int_t)", "Reve::CascadeListEditor", this, "UpdateSelectedTab()");
166   frameTab1->AddFrame(ptab[0], new TGLayoutHints(kLHintsLeft| kLHintsExpandX,0,0,0,0));
167   
168   //------
169   // container of "Tab2"
170   TGCompositeFrame *frameTab2 = pMainTab[0]->AddTab("cascade");
171   frameTab2->SetLayoutManager(new TGVerticalLayout(frameTab2));
172   
173   // tab widget
174   ptab[1] = new TGTab(frameTab2,440,299);
175   ptab[1]->Resize(ptab[1]->GetDefaultSize());
176   ptab[1]->Connect("Selected(Int_t)", "Reve::CascadeListEditor", this, "UpdateSelectedTab()");
177   frameTab2->AddFrame(ptab[1], new TGLayoutHints(kLHintsLeft| kLHintsExpandX ,0,0,0,0));
178   
179   //------
180   // container of "Tab3"
181   TGCompositeFrame *frameTab3 = pMainTab[0]->AddTab("V0daugh/bach.");
182   frameTab3->SetLayoutManager(new TGVerticalLayout(frameTab3));
183   
184   // tab widget
185   ptab[2] = new TGTab(frameTab3,440,299);
186   ptab[2]->Resize(ptab[2]->GetDefaultSize());
187   ptab[2]->Connect("Selected(Int_t)", "Reve::CascadeListEditor", this, "UpdateSelectedTab()");
188   frameTab3->AddFrame(ptab[2], new TGLayoutHints(kLHintsLeft| kLHintsExpandX ,0,0,0,0));
189
190   //------
191   TGCompositeFrame **frameTab = new TGCompositeFrame*[fgkNCanvas];
192   
193   frameTab[0] = AddTab(ptab[0], 0, can, "Xi");
194   frameTab[1] = AddTab(ptab[0], 1, can, "Omega");
195   frameTab[2] = AddTab(ptab[0], 2, can, "Arm.Podo.");
196   frameTab[3] = AddTab(ptab[0], 3, can, "Index");
197
198   frameTab[4] = AddTab(ptab[1], 4, can, "cosPointing");
199   frameTab[5] = AddTab(ptab[1], 5, can, "bachV0DCA");
200   frameTab[6] = AddTab(ptab[1], 6, can, "radius");
201   frameTab[7] = AddTab(ptab[1], 7, can, "Pt");
202   frameTab[8] = AddTab(ptab[1], 8, can, "eta");
203
204   frameTab[9] = AddTab(ptab[2], 9, can, "neg Pt");
205   frameTab[10] = AddTab(ptab[2], 10, can, "neg eta");
206   frameTab[11] = AddTab(ptab[2], 11, can, "pos Pt");
207   frameTab[12] = AddTab(ptab[2], 12, can, "pos eta");
208   frameTab[13] = AddTab(ptab[2], 13, can, "bach Pt");
209   frameTab[14] = AddTab(ptab[2], 14, can, "bach eta");
210
211   pMainTab[0]->SetTab(0);
212   ptab[0]->SetTab(0);
213   ptab[1]->SetTab(0);
214   ptab[2]->SetTab(0);
215
216   Pixel_t darkgrey;
217   gClient->GetColorByName("grey50", darkgrey);
218   ptab[0]->SetBackgroundColor(darkgrey);
219   ptab[1]->SetBackgroundColor(darkgrey);
220   ptab[2]->SetBackgroundColor(darkgrey);
221
222   return frameTab;
223 }
224
225
226 //_________________________________________________________________________
227 void CascadeListEditor::AddValuator(TGCompositeFrame* frame, char *name,
228                                Float_t min, Float_t max, Int_t pres,
229                                char *func, Int_t iHist) {
230
231   TGCompositeFrame* downFrame = new TGCompositeFrame(frame,
232                                     60, 60, kHorizontalFrame);
233
234    // --- Selectors
235   fRange[iHist] = new RGDoubleValuator(downFrame, name, 200, 0);
236   fRange[iHist]->SetNELength(6);
237   fRange[iHist]->Build();
238   fRange[iHist]->GetSlider()->SetWidth(200);
239
240   fRange[iHist]->SetLimits(min, max, TGNumberFormat::kNESRealTwo);
241   if (pres==0)
242     fRange[iHist]->SetLimits(min, max, TGNumberFormat::kNESInteger);
243   else if (pres==3)
244     fRange[iHist]->SetLimits(min, max, TGNumberFormat::kNESRealThree);
245   else if (pres==4)
246     fRange[iHist]->SetLimits(min, max, TGNumberFormat::kNESRealFour);
247   else if (pres==5)
248     fRange[iHist]->SetLimits(min, max, TGNumberFormat::kNESReal);
249
250   fRange[iHist]->Connect("ValueSet()",
251                          "Reve::CascadeListEditor", this, func);
252   downFrame->AddFrame(fRange[iHist], new TGLayoutHints(kLHintsLeft,
253                                                       0, 0, 0, 0));
254
255   TGTextButton* adjustButton = new TGTextButton(downFrame, "Adjust Hist", 40);
256
257   char ch[40];
258   sprintf(ch,"AdjustHist(=%i)",iHist);
259   adjustButton->Connect("Clicked()", "Reve::CascadeListEditor", this, ch);
260   downFrame->AddFrame(adjustButton, new TGLayoutHints(kLHintsTop, 0, 0, 0, 0));
261
262   frame->AddFrame(downFrame, new TGLayoutHints(kLHintsTop| kLHintsExpandY,
263                                                 0, 0, 0, 0));
264 }
265
266
267 //_________________________________________________________________________
268 void CascadeListEditor::AddSelectTab() {
269  
270   TGCompositeFrame** tab = CreateTab(&fMainTabA, fTabA, 1);
271
272   AddValuator(tab[0],  "mass Xi",       0,  5, 3, "MassXiRange()",     0);
273   AddValuator(tab[1],  "mass Omega",    0,  5, 3, "MassOmegaRange()",  1);
274
275   AddValuator(tab[3],  "Index",                  0,   1e5, 0, "IndexRange()",       2);
276   AddValuator(tab[4],  "cos pointing angle",     0.8,   1, 5, "CosPointingRange()", 3);
277   AddValuator(tab[5],  "bach-V0 DCA",            0,     5, 3, "BachV0DCARange()",   4);
278   AddValuator(tab[6],  "radius",                 0,   100, 2, "RadiusRange()",      5);
279   AddValuator(tab[7],  "Pt",                     0,    10, 2, "PtRange()",          6);
280   AddValuator(tab[8],  "Pseudo-rapidity",       -2,     2, 2, "PseudoRapRange()",   7);
281
282   AddValuator(tab[9],  "neg Pt",                 0,    10, 2, "NegPtRange()",       8);
283   AddValuator(tab[10], "neg pseudo-rapidity",   -2,     2, 2, "NegEtaRange()",      9);
284   AddValuator(tab[11], "pos Pt",                 0,    10, 2, "PosPtRange()",      10);
285   AddValuator(tab[12], "pos pseudo-rapidity",   -2,     2, 2, "PosEtaRange()",     11);
286   AddValuator(tab[13], "bach. Pt",               0,    10, 2, "BachPtRange()",     12);
287   AddValuator(tab[14], "bach. pseudo-rapidity", -2,     2, 2, "BachEtaRange()",    13);
288
289   delete[] tab;
290 }
291
292
293 //_________________________________________________________________________
294 void CascadeListEditor::AddSeeTab() {
295
296   TGCompositeFrame** tab = CreateTab(&fMainTabB, fTabB, 2);
297   delete[] tab;
298 }
299
300
301 //_________________________________________________________________________
302 void CascadeListEditor::AdjustHist(Int_t iHist) {
303
304   if (! fMList) return;
305   fMList->AdjustHist(iHist);
306
307   FillCanvas();
308 }
309
310 //_________________________________________________________________________
311 void CascadeListEditor::ResetCuts() {
312
313   if (! fMList) return;
314
315   Float_t min,max;
316
317   for (Int_t i=0; i<fgkNRange;i++) {
318     
319     if (i==2) continue;
320     min = fRange[i]->GetLimitMin();
321     max = fRange[i]->GetLimitMax();
322     fMList->SetMin(i, min);
323     fMList->SetMax(i, max);
324     fRange[i]->SetValues(min, max);
325     fMList->AdjustHist(i);
326   }
327
328   // for the Index we scan its actual range
329   Int_t imin, imax;
330   fMList->GetCasIndexRange(imin, imax);
331   if (imin<imax) {
332     Int_t minH = imin-(imax-imin)/20;
333     Int_t maxH = imax+(imax-imin)/20;
334     fMList->SetMin(2, minH);
335     fMList->SetMax(2, maxH);
336     fRange[2]->SetLimits(minH, maxH, TGNumberFormat::kNESInteger);
337     fRange[2]->SetValues(minH, maxH);
338     fMList->AdjustHist(2);
339
340   }
341   FillCanvas();
342   Update();
343 }
344
345
346 //_________________________________________________________________________
347 void CascadeListEditor::FillCanvas() {
348
349   fMList->FilterAll();
350
351   TCanvas *c1, *c1b;
352   TH1F *hist = 0;
353   TH2F *hist2D = 0;
354   Bool_t is2D;
355
356   Int_t canvasMap[fgkNCanvas]={0,1,1000,2,3,4,5,6,7,8,9,10,11,12,13};
357
358   for (Int_t i=0; i<fgkNCanvas; i++)
359     if (fCanvasA[i]) {
360
361       is2D = canvasMap[i]>999;
362
363       if (is2D) hist2D = fMList->GetHist2D(canvasMap[i]-1000);
364       else hist = fMList->GetHist(canvasMap[i]);
365
366       c1 = fCanvasA[i]->GetCanvas();
367       c1->cd();
368       if (is2D) hist2D->Draw("colz"); else hist->Draw();
369       c1->Modified();
370       c1->Update();
371       
372       c1b = fCanvasB[i]->GetCanvas();
373       c1b->cd();
374       if (is2D) hist2D->Draw("colz"); else hist->Draw();
375       c1b->Modified();
376       c1b->Update();
377   }
378   UpdateSelectedTab();
379 }
380 //_________________________________________________________________________
381 void CascadeListEditor::UpdateSelectedTab() {
382
383   Pixel_t yellow;
384   Pixel_t grey;
385   gClient->GetColorByName("yellow", yellow);
386   gClient->GetColorByName("grey", grey);
387
388   TGTabElement *tabElem; 
389   for (Int_t i=0; i<fMainTabA->GetNumberOfTabs(); i++) {
390
391     tabElem = fMainTabA->GetTabTab(i);
392     tabElem->ChangeBackground(grey);
393
394     for (Int_t j=0; j<fTabA[i]->GetNumberOfTabs();j++) {
395       tabElem = fTabA[i]->GetTabTab(j);
396       tabElem->ChangeBackground(grey);
397     }
398   }
399
400   Int_t currentTab = fMainTabA->GetCurrent();
401   Int_t currentSubTab = fTabA[currentTab]->GetCurrent();
402   tabElem = fMainTabA->GetTabTab(currentTab);
403   tabElem->ChangeBackground(yellow);
404   tabElem = fTabA[currentTab]->GetTabTab(currentSubTab);
405   tabElem->ChangeBackground(yellow);
406
407   TCanvas *c1;
408   Int_t iCan = 0;
409   Int_t i=0;
410
411   while (currentTab>i) {
412     iCan += fTabA[i]->GetNumberOfTabs();
413     i++;
414   }
415   iCan += currentSubTab;
416   c1 = fCanvasA[iCan]->GetCanvas();
417   c1->GetCanvas()->Modified();
418   c1->GetCanvas()->Update();
419
420   //---
421
422   for (Int_t i=0; i<fMainTabB->GetNumberOfTabs(); i++) {
423
424     tabElem = fMainTabB->GetTabTab(i);
425     tabElem->ChangeBackground(grey);
426
427     for (Int_t j=0; j<fTabB[i]->GetNumberOfTabs();j++) {
428       tabElem = fTabB[i]->GetTabTab(j);
429       tabElem->ChangeBackground(grey);
430     }
431   }
432
433   currentTab = fMainTabB->GetCurrent();
434   currentSubTab = fTabB[currentTab]->GetCurrent();
435   tabElem = fMainTabB->GetTabTab(currentTab);
436   tabElem->ChangeBackground(yellow);
437   tabElem = fTabB[currentTab]->GetTabTab(currentSubTab);
438   tabElem->ChangeBackground(yellow);
439
440   iCan = 0;
441   i=0;
442
443   while (currentTab>i) {
444     iCan += fTabB[i]->GetNumberOfTabs();
445     i++;
446   }
447   iCan += currentSubTab;
448   c1 = fCanvasB[iCan]->GetCanvas();
449   c1->GetCanvas()->Modified();
450   c1->GetCanvas()->Update();
451 }
452
453
454 //_________________________________________________________________________
455 void CascadeListEditor::UpdateAll(Int_t iCanA) {
456
457   TCanvas *c1 = fCanvasA[iCanA]->GetCanvas();
458   c1->Modified();
459   c1->Update();
460
461   static Int_t iCan, i;
462   iCan=0;
463   i=0;
464   while (fMainTabB->GetCurrent()>i) {
465     iCan += fTabB[i]->GetNumberOfTabs();
466     i++;
467   }
468   iCan += fTabB[i]->GetCurrent();
469   c1 = fCanvasB[iCan]->GetCanvas();
470   c1->GetCanvas()->Modified();
471   c1->GetCanvas()->Update();
472
473   Update();
474 }
475
476
477 //_________________________________________________________________________
478
479 void CascadeListEditor::DoRnrVtx()
480 {
481   fMList->SetRnrCasVtx(fRnrVtx->IsOn());
482   Update();
483 }
484
485 void CascadeListEditor::DoRnrV0path()
486 {
487   fMList->SetRnrV0path(fRnrV0path->IsOn());
488   Update();
489 }
490
491 void CascadeListEditor::DoRnrV0Daughters()
492 {
493   fMList->SetRnrV0Daughters(fRnrV0Daughters->IsOn());
494   Update();
495 }
496
497 void CascadeListEditor::DoRnrBach()
498 {
499   fMList->SetRnrBachelor(fRnrBach->IsOn());
500   Update();
501 }
502
503 void CascadeListEditor::DoRnrCasPath()
504 {
505   fMList->SetRnrCasPath(fRnrCasPath->IsOn());
506   Update();
507 }
508
509
510 //_________________________________________________________________________
511 void CascadeListEditor::MassXiRange() {
512
513   fMList->XiMassFilter(fRange[0]->GetMin(), fRange[0]->GetMax());
514   UpdateAll(0);
515 }
516
517 //_________________________________________________________________________
518   void CascadeListEditor::MassOmegaRange() {
519   fMList->OmegaMassFilter(fRange[1]->GetMin(), fRange[1]->GetMax());
520   UpdateAll(1);
521 }
522
523 //_________________________________________________________________________
524   void CascadeListEditor::IndexRange() {
525   fMList->IndexFilter(fRange[2]->GetMin(), fRange[2]->GetMax());
526   UpdateAll(3);
527 }
528 //_________________________________________________________________________
529   void CascadeListEditor::CosPointingRange() {
530   fMList->CosPointingFilter(fRange[3]->GetMin(), fRange[3]->GetMax());
531   UpdateAll(4);
532 }
533
534 //_________________________________________________________________________
535   void CascadeListEditor::BachV0DCARange() {
536   fMList->BachV0DCAFilter(fRange[4]->GetMin(), fRange[4]->GetMax());
537   UpdateAll(5);
538 }
539
540 //_________________________________________________________________________
541   void CascadeListEditor::RadiusRange() {
542   fMList->RadiusFilter(fRange[5]->GetMin(), fRange[5]->GetMax());
543   UpdateAll(6);
544 }
545
546 //_________________________________________________________________________
547   void CascadeListEditor::PtRange() {
548   fMList->PtFilter(fRange[6]->GetMin(), fRange[6]->GetMax());
549   UpdateAll(7);
550 }
551
552 //_________________________________________________________________________
553   void CascadeListEditor::PseudoRapRange() {
554   fMList->PseudoRapFilter(fRange[7]->GetMin(), fRange[7]->GetMax());
555   UpdateAll(8);
556 }
557
558 //_________________________________________________________________________
559   void CascadeListEditor::NegPtRange() {
560   fMList->NegPtFilter(fRange[8]->GetMin(), fRange[8]->GetMax());
561   UpdateAll(9);
562 }
563
564 //_________________________________________________________________________
565   void CascadeListEditor::NegEtaRange() {
566   fMList->NegEtaFilter(fRange[9]->GetMin(), fRange[9]->GetMax());
567   UpdateAll(10);
568 }
569
570 //_________________________________________________________________________
571   void CascadeListEditor::PosPtRange() {
572   fMList->PosPtFilter(fRange[10]->GetMin(), fRange[10]->GetMax());
573   UpdateAll(11);
574 }
575
576 //_________________________________________________________________________
577   void CascadeListEditor::PosEtaRange() {
578   fMList->PosEtaFilter(fRange[11]->GetMin(), fRange[11]->GetMax());
579   UpdateAll(12);
580 }
581
582 //_________________________________________________________________________
583   void CascadeListEditor::BachPtRange() {
584   fMList->BachPtFilter(fRange[12]->GetMin(), fRange[12]->GetMax());
585   UpdateAll(13);
586 }
587
588 //_________________________________________________________________________
589   void CascadeListEditor::BachEtaRange() {
590   fMList->BachEtaFilter(fRange[13]->GetMin(), fRange[13]->GetMax());
591   UpdateAll(14);
592 }
593
594