AliPhysics  b6a3523 (b6a3523)
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
TupleSelector.C
Go to the documentation of this file.
1 
38 #ifndef SELECTOR_C
39 #define SELECTOR_C
40 
41 #include <TSelector.h>
42 #include <TQObject.h>
43 #ifndef __CINT__
44 # include <TH1.h>
45 # include <TString.h>
46 # include <TCanvas.h>
47 # include <TLegend.h>
48 # include <TLegendEntry.h>
49 # include <TError.h>
50 # include <TMath.h>
51 # include <THStack.h>
52 # include <TTree.h>
53 # include <TClonesArray.h>
54 # include <TChain.h>
55 # include <TSystem.h>
56 # include <TOutputListSelectorDataMap.h>
57 # include <TLatex.h>
58 # include <TROOT.h>
59 # include <TFile.h>
60 # include <TDirectory.h>
61 # include <TSystemDirectory.h>
62 # include <TRegexp.h>
63 # include <TKey.h>
64 # include <TFileCollection.h>
65 # include <THashList.h>
66 # include <TDSet.h>
67 # include <TQConnection.h>
68 # include <iostream>
69 # include <TProof.h>
70 # include <TStopwatch.h>
71 # include "AliFMDMCTrackELoss.h"
72 #else
73 // Forward declarations of types in the interface, and to trigger
74 // autoloading of some libraries
75 class TTree;
76 class TChain;
77 class TCanvas;
78 class TH1;
79 class THStack;
80 class TLegend;
81 class TCanvas;
82 class TString;
83 class TDirectory;
84 class TSystemDirectory;
85 class TRegexp;
86 class TDSet;
87 class TProof;
88 class AliFMDMCTrackELoss;
90 #endif
91 
92 //====================================================================
96 struct Spectra : public TObject
97 {
108  : TObject(), fName(""), fPrimary(0), fSecondary(0)
109  {}
118  Spectra(const char* name,
119  const char* title,
120  Color_t color,
121  const TArrayD& bins)
122  : TObject(),
123  fName(name),
124  fPrimary(0),
125  fSecondary(0)
126  {
127  fPrimary = new TH1D(Form("primary%s", name),
128  title, bins.GetSize()-1,
129  bins.GetArray());
130  fPrimary->SetXTitle("#beta#gamma");
131  fPrimary->SetMarkerColor(color);
132  fPrimary->SetLineColor(color);
133  fPrimary->SetMarkerStyle(20);
134  fPrimary->SetFillStyle(0);
135  fPrimary->SetDirectory(0);
136  fPrimary->Sumw2();
137 
138  fSecondary = static_cast<TH1*>(fPrimary->Clone(Form("secondary%s",name)));
139  fSecondary->SetMarkerStyle(24);
140  fSecondary->SetMarkerSize(1.2);
141  fSecondary->SetDirectory(0);
142  }
148  Spectra(const Spectra& o)
149  : TObject(o),
150  fName(o.fName),
151  fPrimary(0),
152  fSecondary(0)
153  {
154  if (o.fPrimary)
155  fPrimary = static_cast<TH1*>(o.fPrimary->Clone());
156  if (o.fSecondary)
157  fSecondary = static_cast<TH1*>(o.fSecondary->Clone());
158  }
162  virtual ~Spectra()
163  {
164  if (fPrimary) delete fPrimary;
165  if (fSecondary) delete fSecondary;
166  }
175  {
176  if (&o == this) return *this;
177 
178  TObject::operator=(o);
179  fName = o.fName;
180  if (fPrimary) { delete fPrimary; fPrimary = 0; }
181  if (fSecondary) { delete fSecondary; fSecondary = 0; }
182  if (o.fPrimary)
183  fPrimary = static_cast<TH1*>(o.fPrimary->Clone());
184  if (o.fSecondary)
185  fSecondary = static_cast<TH1*>(o.fSecondary->Clone());
186 
187  return *this;
188  }
194  const char* GetName() const { return fName.Data(); }
201  void Fill(Bool_t primary, Double_t x)
202  {
203  if (primary) fPrimary->Fill(x);
204  else fSecondary->Fill(x);
205  }
214  TH1* Hist(Bool_t primary) const
215  {
216  return (primary ? fPrimary : fSecondary);
217  }
224  virtual void Stack(THStack* stack, TLegend* leg)
225  {
226  TH1* hists[] = { fPrimary, fSecondary, 0 };
227  for (Int_t i = 0; i < 2; i++) {
228  if (!hists[i]) continue;
229 
230  Int_t n = hists[i]->GetEntries();
231  if (n <= 0) continue;
232 
233  if (stack) {
234  hists[i]->Scale(1. / n, "width");
235  stack->Add(hists[i]);
236 
237  if (i != 0) continue;
238  if (!leg) continue;
239 
240  leg->AddEntry(hists[i], hists[i]->GetTitle(), "p");
241  }
242  else if (leg) {
243  TLegendEntry* e = leg->AddEntry("",
244  (i == 0 ? "Primary" : "Secondary"),
245  "p");
246  e->SetMarkerStyle(hists[i]->GetMarkerStyle());
247  e->SetFillStyle(hists[i]->GetFillStyle());
248  e->SetFillColor(kBlack);
249  }
250  }
251  }
260  {
261  TIter nxt(list);
262  TObject* obj = 0;
263  Long64_t cnt = 0;
264  while ((obj = nxt())) {
265  if (!obj->IsA()->InheritsFrom(this->IsA())) {
266  Warning("Merge", "Will not merge a Spectra with a %s",
267  obj->IsA()->GetName());
268  continue;
269  }
270  Spectra* oth = static_cast<Spectra*>(obj);
271  if (Add(oth)) cnt++;
272  }
273  // Info("Merge", "%s merged %lld entries", fName.Data(), cnt);
274  return cnt;
275  }
276  Bool_t Add(const Spectra* oth)
277  {
278  if (!oth->fName.EqualTo(fName)) {
279  Warning("Merge", "Will not merge %s with %s",
280  oth->fName.Data(), fName.Data());
281  return false;
282  }
283  // Int_t nPrmOld = fPrimary->GetEntries();
284  // Int_t nSecOld = fPrimary->GetEntries();
285 
286  fPrimary ->Add(oth->fPrimary);
287  fSecondary->Add(oth->fSecondary);
288 
289  // Int_t nPrmNow = fPrimary->GetEntries();
290  // Int_t nSecNow = fPrimary->GetEntries();
291 
292  // Info("Add", "%s: %d/%d merged to %d/%d",
293  // fName.Data(), nPrmOld, nSecOld, nPrmNow, nSecNow);
294  return true;
295  }
301  void ls(Option_t* option="") const
302  {
303  TObject::ls(option);
304  gROOT->IncreaseDirLevel();
305  TH1* hists[] = { fPrimary,
306  fSecondary,
307  0 };
308  for (Int_t i = 0; i < 2; i++) if (hists[i]) hists[i]->ls();
309  gROOT->DecreaseDirLevel();
310  }
316  void Store(TDirectory* dir)
317  {
318  TDirectory* out = dir->mkdir(fName);
319  out->cd();
320  fPrimary->Clone("primary")->Write();
321  fSecondary->Clone("secondary")->Write();
322  dir->cd();
323  }
324  ClassDef(Spectra,1);
325 };
326 //====================================================================
330 struct Type : public Spectra
331 {
335  Type() : Spectra() {}
343  Type(const char* name,
344  const char* title,
345  Color_t color)
346  : Spectra() // Do not use the base class ctor
347  {
348  fName = name;
349  fPrimary = new TH1D(Form("primary%s", name),
350  title, 6, .5, 6.5);
351  fPrimary->SetXTitle("particle type");
352  fPrimary->SetMarkerColor(color);
353  fPrimary->SetMarkerStyle(20);
354  fPrimary->SetFillColor(color);
355  fPrimary->SetFillStyle(3001);
356  fPrimary->SetLineColor(color);
357  fPrimary->SetDirectory(0);
358  fPrimary->Sumw2();
359  fPrimary->GetXaxis()->SetBinLabel(1, "e^{#pm}");
360  fPrimary->GetXaxis()->SetBinLabel(2, "#pi^{#pm}");
361  fPrimary->GetXaxis()->SetBinLabel(3, "K^{#pm}");
362  fPrimary->GetXaxis()->SetBinLabel(4, "p/#bar{p}");
363  fPrimary->GetXaxis()->SetBinLabel(5, "strange");
364  fPrimary->GetXaxis()->SetBinLabel(6, "other");
365 
366  fSecondary =
367  static_cast<TH1*>(fPrimary->Clone(Form("secondary%s",name)));
368  fSecondary->SetMarkerStyle(24);
369  fSecondary->SetMarkerSize(1.2);
370  fSecondary->SetFillStyle(3002);
371  fSecondary->SetDirectory(0);
372  }
373  ClassDef(Type,1);
374 };
375 
376 //====================================================================
377 struct Ring : public TObject
378 {
379 
388 
390  enum {
396  };
397 
401  Ring()
402  : TObject(),
403  fD(0),
404  fR('\0'),
405  fName(),
406  fSpectra(0)
407  {}
415  : fD(d),
416  fR(r),
417  fName(Form("FMD%d%c", d, r)),
418  fSpectra(0)
419  {
420  Color_t color = AliForwardUtil::RingColor(fD, fR);
421  TArrayD bgArray(601);
422  AliForwardUtil::MakeLogScale(300, -2, 5, bgArray);
423 
425  fSpectra->SetOwner(true);
426  fSpectra->AddAt(new Spectra("BetaGamma", fName, color, bgArray),
427  kBetaGamma);
428  fSpectra->AddAt(new Spectra("BetaGammaPi", fName, color, bgArray),
429  kBetaGammaPi);
430  fSpectra->AddAt(new Spectra("BetaGammaElectron", fName, color, bgArray),
432  fSpectra->AddAt(new Type("Type", fName, color), kTypes);
433  // fSpectra->ls();
434  }
440  Ring(const Ring& r)
441  : TObject(r),
442  fD(r.fD),
443  fR(r.fR),
444  fName(r.fName),
445  fSpectra(0)
446  {
447  if (r.fSpectra) fSpectra = new TObjArray(*r.fSpectra);
448  }
456  Ring& operator=(const Ring& r)
457  {
458  if (&r == this) return *this;
459 
460  fD = r.fD;
461  fR = r.fR;
462  fName = r.fName;
463  if (fSpectra) { delete fSpectra; fSpectra = 0; }
464  if (r.fSpectra) fSpectra = new TObjArray(*r.fSpectra);
465 
466  return *this;
467  }
468  const char* GetName() const { return fName.Data(); }
476  Spectra* Get(UInt_t which)
477  {
478  if (!fSpectra) return 0;
479  Int_t i = which;
480  if (i > fSpectra->GetEntriesFast()) return 0;
481  return static_cast<Spectra*>(fSpectra->At(i));
482  }
490  const Spectra* Get(UInt_t which) const
491  {
492  if (!fSpectra) return 0;
493  Int_t i = which;
494  if (i > fSpectra->GetEntriesFast()) return 0;
495  return static_cast<Spectra*>(fSpectra->At(i));
496  }
503  {
504  Bool_t prim = hit->IsPrimary();
505  Int_t apdg = hit->AbsPdg();
506  Int_t mfl = Int_t(apdg/TMath::Power(10,Int_t(TMath::Log10(apdg))));
507  Int_t type = (hit->IsElectron() ? 1 :
508  hit->IsPion() ? 2 :
509  hit->IsKaon() ? 3 :
510  hit->IsProton() ? 4 :
511  mfl == 3 ? 5 : 6);
512 
513  Get(kBetaGamma) ->Fill(prim, hit->BetaGamma());
514  Get(kTypes) ->Fill(prim, type);
515  if (hit->IsPion() || hit->IsProton() || hit->IsKaon() || apdg > 100)
516  Get(kBetaGammaPi) ->Fill(prim, hit->BetaGamma());
517  if (hit->IsElectron())
518  Get(kBetaGammaElectron)->Fill(prim, hit->BetaGamma());
519 
520  }
529  TH1* Hist(Bool_t primary, UShort_t which) const
530  {
531  const Spectra* spe = Get(which);
532  if (!spe) return 0;
533  return spe->Hist(primary);
534  }
543  {
544  TIter nxt(list);
545  TObject* obj = 0;
546  Long64_t cnt = 0;
547  while ((obj = nxt())) {
548  if (!obj->IsA()->InheritsFrom(this->IsA())) {
549  Warning("Merge", "Will not merge a Ring with a %s",
550  obj->IsA()->GetName());
551  continue;
552  }
553  Ring* oth = static_cast<Ring*>(obj);
554  if (oth->fD != fD || oth->fR != fR) {
555  Warning("Merge", "Will not merge FMD%d%c with FMD%d%c",
556  fD, fR, oth->fD, oth->fR);
557  continue;
558  }
559 
560  if (fSpectra) {
561  // fSpectra->Merge(oth->fSpectra);
562  TIter thsNext(fSpectra);
563  TIter othNext(oth->fSpectra);
564  Spectra* thsSpec = 0;
565  Spectra* othSpec = 0;
566 
567  while ((thsSpec = static_cast<Spectra*>(thsNext())) &&
568  (othSpec = static_cast<Spectra*>(othNext())))
569  thsSpec->Add(othSpec);
570 
571  }
572  cnt++;
573  }
574  // Info("Merge", "FMD%d%c merged %lld entries", fD, fR, cnt);
575  return cnt;
576  }
577  void Draw(Option_t* opt="")
578  {
579  Info("Draw", "Should draw rings here");
580  ls(opt);
581  }
587  void ls(Option_t* option="") const
588  {
589  TObject::ls(option);
590  gROOT->IncreaseDirLevel();
591  if (fSpectra) fSpectra->ls(option);
592  gROOT->DecreaseDirLevel();
593  }
599  void Store(TDirectory* dir)
600  {
601  TDirectory* out = dir->mkdir(fName);
602  out->cd();
603  TIter next(fSpectra);
604  Spectra* spec = 0;
605  while ((spec = static_cast<Spectra*>(next())))
606  spec->Store(out);
607  dir->cd();
608  }
609 
610  ClassDef(Ring,2);
611 };
612 
613 //====================================================================
614 struct Monitor : public TObject, public TQObject
615 {
617  {
618  if (!gProof) return;
619 
620  fName = gProof->GetSessionTag();
621  gDirectory->Add(this);
622  // _must_ specify signal _exactly_ like below, or we won't be called
623  Bool_t ret = gProof->Connect("Feedback(TList *objs)", "Monitor", this,
624  "Feedback(TList *objs)");
625  if (!ret)
626  Warning("Monitor", "Failed to connect to Proof");
627  }
628  virtual ~Monitor()
629  {
630  if (!gProof) return;
631  gProof->Disconnect("Feedback(TList *objs)",this,
632  "Feedback(TList* objs)");
633  }
634  void SetName(const char* name) { fName = name; }
635  const char* GetName() const { return fName.Data(); }
636  void Feedback(TList* objs)
637  {
638  Info("Feedback", "Got a list of objects (%p)", objs);
639  if (!objs) return;
640  objs->ls();
641  }
643  ClassDef(Monitor,1);
644 };
645 
646 
647 //====================================================================
648 struct TupleSelector : public TSelector
649 {
652  TClonesArray* fHits;
655 
661  TupleSelector(const char* name="") :
662  fTitle(name), fTree(0), fHits(0), fI(0), fRings(0)
663  {
664  if (fTitle.IsNull())
665  fTitle = gSystem->BaseName(gSystem->WorkingDirectory());
666  // SetOption("fb=rings");
667 
668  }
669  const char* GetName() const { return fTitle.Data(); }
683  {
684  Bool_t inner = (r == 'i' || r == 'I');
685  switch (d) {
686  case 1: return 0;
687  case 2: return (inner ? 1 : 2);
688  case 3: return (inner ? 4 : 3);
689  }
690  ::Warning("", "Unknown ring FMD%d%c", d, r);
691  return 0xFFFF;
692  }
693  void Init(TTree* tree)
694  {
695  Info("Init", "Got a tree: %p", tree);
696 
697  if (!tree) {
698  Warning("Init", "No tree passed");
699  return;
700  }
701 
702  fTree = tree;
703  fHits = new TClonesArray("AliFMDMCTrackELoss::Hit");
704  fTree->SetBranchAddress("hits", &fHits);
705  }
706  void Begin(TTree* tree)
707  {
708  Info("Begin", "Called w/tree @ %p", tree);
709  // if (tree) SlaveBegin(tree);
710  new Monitor;
711  }
717  void SlaveBegin(TTree* tree)
718  {
719  Info("SlaveBegin", "Got a tree: %p", tree);
720  fRings = new TObjArray(5);
721  fRings->SetName("rings");
722  fRings->SetOwner(true);
723  fRings->AddAt(new Ring(1,'I'), Index(1,'I'));
724  fRings->AddAt(new Ring(2,'I'), Index(2,'I'));
725  fRings->AddAt(new Ring(2,'O'), Index(2,'O'));
726  fRings->AddAt(new Ring(3,'O'), Index(3,'O'));
727  fRings->AddAt(new Ring(3,'I'), Index(3,'I'));
728 
729  fOutput->Add(fRings);
730  if (tree) {
731  // In case of local mode, we get the tree here,
732  // and init isn't called
733  Init(tree);
734  }
735  }
736  /* @} */
737 
743  {
744  TFile* file = (fTree ? fTree->GetCurrentFile() : 0);
745  Info("Notify","processing file: %p (%s)",
746  file, (file ? file->GetName() : "nil"));
747  return true;
748  }
757  {
758  if (!fTree) {
759  Warning("Process", "No tree");
760  return false;
761  }
762  if (!fTree->GetTree()) {
763  Warning("Process", "No real tree");
764  return false;
765  }
766  fHits->Clear();
767  fTree->GetTree()->GetEntry(entry);
768  fI++;
769  // if (fI % 100 == 0) {
770  // printf("Event # %6d (%6d hits)\r", fI, fHits->GetEntries());
771  // fflush(stdout);
772  // }
773  // Printf("Event # %7d, %6d hits", fI, fHits->GetEntries());
774 
775  TIter next(fHits);
776  AliFMDMCTrackELoss::Hit* hit = 0;
777  while ((hit = static_cast<AliFMDMCTrackELoss::Hit*>(next())))
778  Fill(hit);
779  return true;
780  }
790  {
791  return static_cast<Ring*>(fRings->At(Index(d,r)));
792  }
801  const Ring* Get(UShort_t d, Char_t r) const
802  {
803  return static_cast<Ring*>(fRings->At(Index(d,r)));
804  }
811  {
812  Ring* r = Get(hit->D(), hit->R());
813  r->Fill(hit);
814  }
815  /* @} */
816 
822  {
823  }
828  void Terminate()
829  {
830  Printf("\nDone (rings=%p)", fRings);
831 
832  if (!fRings)
833  fRings = static_cast<TObjArray*>(fOutput->FindObject("rings"));
834  if (!fRings) {
835  Error("Terminate", "No rings in output array");
836  return;
837  }
838  // fRings->ls();
839 
840  TFile* out = TFile::Open("tuple_summary.root", "RECREATE");
841  Store(out);
842 
843  PlotOne(Ring::kBetaGamma, "#beta#gamma - all", out);
844  PlotOne(Ring::kBetaGammaPi, "#beta#gamma - Hadrons", out);
845  PlotOne(Ring::kBetaGammaElectron,"#beta#gamma - e^{#pm}", out);
846  PlotOne(Ring::kTypes, "Relative abundance", out);
847  }
856  void PlotOne(UShort_t which, const char* title,
857  TDirectory* dir=0, Option_t* opt="")
858  {
859  static Int_t cnt = 1;
860  const char* name = Form("c%02d", cnt++);
861 
862  TCanvas* c = new TCanvas(name, title, 1000, 1000);
863  c->SetFillColor(0);
864  c->SetFillStyle(0);
865  c->SetTopMargin(0.01);
866  c->SetRightMargin(0.03);
867 
868  TLegend* leg = new TLegend(0.65, 0.65, 0.975, 0.975);
869  leg->SetFillStyle(0);
870  leg->SetBorderSize(0);
871 
872  TString xTitle = "";
873  THStack* stack = new THStack(name, title);
874  for (Int_t i = 0; i < 5; i++) {
875  Ring* ring = static_cast<Ring*>(fRings->At(i));
876  Spectra* spec = ring->Get(which);
877  if (!spec) {
878  Warning("PlotOne", "No spectra for %d", which);
879  continue;
880  }
881 
882  // Add to our stack
883  spec->Stack(stack, leg);
884 
885  // Get X title
886  if (xTitle.IsNull() && spec->Hist(true))
887  xTitle = spec->Hist(true)->GetXaxis()->GetTitle();
888 
889  // If we're not at the last one continue
890  if (i < 4) continue;
891 
892  // Add primary/secondary entry to legend
893  spec->Stack(0, leg);
894 
895  }
896  c->cd();
897  if (!stack->GetHists() ||
898  stack->GetHists()->GetEntries() < 0) {
899  Warning("PlotOne", "No histograms added");
900  return;
901  }
902  stack->Draw(Form("nostack %s", opt));
903  stack->GetHistogram()->SetXTitle(xTitle);
904  stack->GetHistogram()->GetListOfFunctions()->Add(leg);
905 
906  leg->Draw();
907 
908  if (which != Ring::kTypes) {
909  c->SetLogx();
910  c->SetLogy();
911  }
912 
913  TLatex* ltx = new TLatex(0.02, 0.02, fTitle);
914  ltx->SetNDC();
915  ltx->SetTextColor(kRed+2);
916  ltx->SetTextAlign(11);
917  ltx->SetTextSize(0.02);
918  ltx->Draw();
919  stack->GetHistogram()->GetListOfFunctions()->Add(ltx);
920 
921  if (dir) {
922  dir->cd();
923  stack->Clone()->Write();
924  }
925 
926 
927  c->Modified();
928  c->Update();
929  c->cd();
930 
931  c->Print(Form("%s.pdf", name));
932  }
938  void Store(TDirectory* dir)
939  {
940  TDirectory* out = dir->mkdir(fTitle);
941  out->cd();
942  TIter next(fRings);
943  Ring* ring = 0;
944  while ((ring = static_cast<Ring*>(next())))
945  ring->Store(out);
946  dir->cd();
947  }
948  /* @} */
949 
959  Int_t Version() const { return 1; }
965  Long64_t GetStatus() const { return fI; }
966 
967  /* @} */
968 
969  //------------------------------------------------------------------
974  //------------------------------------------------------------------
982  static TObject* GetChainOrDataSet(const char* file)
983  {
984  if (gSystem->AccessPathName(file)) return 0;
985 
986  TFile* in = TFile::Open(file, "READ");
987  if (!in) return 0;
988 
989  TObject* ret = in->Get("tree");
990  if (!ret) {
991  in->Close();
992  return 0;
993  }
994 
995  if (ret->IsA()->InheritsFrom(TChain::Class()))
996  static_cast<TChain*>(ret)->SetDirectory(0);
997  else if (ret->IsA()->InheritsFrom(TDSet::Class()))
998  static_cast<TDSet*>(ret)->SetDirectory(0);
999  else {
1000  ::Warning("GetChainOrDataSet", "Found object is a %s",
1001  ret->IsA()->GetName());
1002  ret = 0;
1003  }
1004  in->Close();
1005  return ret;
1006  }
1007  //------------------------------------------------------------------
1017  static TDSet* MakeDataSet(const TString& src=".",
1018  Bool_t recursive=false,
1019  Bool_t verbose=false)
1020  {
1021  TString dsFile(Form("%s/dataset.root", src.Data()));
1022  TDSet* dataset = static_cast<TDSet*>(GetChainOrDataSet(dsFile));
1023  if (dataset) {
1025  return dataset;
1026  }
1027 
1028  TChain* c = DoMakeChain(src, recursive, verbose);
1029  if (!c) return 0;
1030 
1031  dataset = new TDSet(*c, false);
1032  dataset->SetName("tree");
1033  dataset->SetLookedUp();
1034  dataset->Validate();
1035 
1036  delete c;
1037  if (dataset) {
1038  TFile* out = TFile::Open(dsFile, "RECREATE");
1039  dataset->Write();
1040  dataset->SetDirectory(0);
1041  out->Close();
1042  }
1043 
1044  return dataset;
1045 
1046  }
1047  //------------------------------------------------------------------
1057  static TChain* MakeChain(const TString& src=".",
1058  Bool_t recursive=false,
1059  Bool_t verbose=false)
1060  {
1061  TString chFile(Form("%s/chain.root", src.Data()));
1062  if (!gSystem->AccessPathName(chFile)) {
1063  TFile* in = TFile::Open(chFile, "READ");
1064  if (in) {
1065  TChain* ret = static_cast<TChain*>(in->Get("tree"));
1066  if (ret) {
1067  ret->SetDirectory(0);
1068  // ret->GetListOfFiles()->ls();
1069  }
1070  in->Close();
1071  if (ret) return ret;
1072  }
1073  }
1074 
1075  TChain* chain = DoMakeChain(src, recursive, verbose);
1076  if (chain) {
1077  TFile* out = TFile::Open(chFile, "RECREATE");
1078  chain->Write();
1079  chain->SetDirectory(0);
1080  out->Close();
1081  }
1082 
1083  return chain;
1084  }
1085  //------------------------------------------------------------------
1095  static TChain* DoMakeChain(const TString& src=".",
1096  Bool_t recursive=false,
1097  Bool_t verbose=false)
1098  {
1099  TChain* chain = new TChain("tree");
1100  TString savdir(gSystem->WorkingDirectory());
1101  TSystemDirectory d(gSystem->BaseName(src.Data()), src.Data());
1102  if (!ScanDirectory(&d, chain, "forward_tree_*", recursive, verbose)) {
1103  delete chain;
1104  chain = 0;
1105  }
1106  else if (!chain->GetListOfFiles() ||
1107  chain->GetListOfFiles()->GetEntries() <= 0) {
1108  delete chain;
1109  chain = 0;
1110  }
1111  return chain;
1112  }
1113 
1114  //------------------------------------------------------------------
1127  static Bool_t ScanDirectory(TSystemDirectory* dir,
1128  TChain* chain,
1129  const TString& pattern,
1130  Bool_t recursive,
1131  Bool_t verbose)
1132  {
1133  if (!dir) {
1134  ::Error("ScanDirectory", "No diretory passed");
1135  return false;
1136  }
1137  if (verbose) ::Info("ScanDirectory", "Scanning %s", dir->GetName());
1138 
1139  Bool_t ret = false;
1140  TRegexp wild(pattern, true);
1141  TString oldDir(gSystem->WorkingDirectory());
1142  TList* files = dir->GetListOfFiles();
1143 
1144  if (!gSystem->ChangeDirectory(oldDir)) {
1145  ::Error("ScanDirectory", "Failed to go back to %s",
1146  oldDir.Data());
1147  return false;
1148  }
1149  if (!files) {
1150  ::Warning("ScanDirectory", "No files");
1151  return false;
1152  }
1153 
1154  TList toAdd;
1155  toAdd.SetOwner();
1156 
1157  // Sort list of files and check if we should add it
1158  files->Sort();
1159  TIter next(files);
1160  TSystemFile* file = 0;
1161  while ((file = static_cast<TSystemFile*>(next()))) {
1162  TString name(file->GetName());
1163  TString title(file->GetTitle());
1164  TString full(gSystem->ConcatFileName(file->GetTitle(), name.Data()));
1165  if (file->IsA()->InheritsFrom(TSystemDirectory::Class())) full = title;
1166 
1167 
1168  if (name.EqualTo(".")||name.EqualTo("..")) {
1169  // Ignore special links
1170  if (verbose) ::Info("ScanDirectory", "Ignoring special %s",
1171  name.Data());
1172  continue;
1173  }
1174 
1175  if (verbose) ::Info("ScanDirectory", "Looking at %s", full.Data());
1176 
1177  FileStat_t fs;
1178  if (gSystem->GetPathInfo(full.Data(), fs)) {
1179  ::Warning("ScanDirectory", "Cannot stat %s (%s)",
1180  full.Data(), gSystem->WorkingDirectory());
1181  continue;
1182  }
1183  // Check if this is a directory
1184  if (file->IsDirectory(full)) {
1185  if (verbose) ::Info("ScanDirectory", "Got a directory: %s",
1186  full.Data());
1187  if (recursive) {
1188  // if (title[0] == '/')
1189  TSystemDirectory* d = new TSystemDirectory(file->GetName(),
1190  full.Data());
1191  if (ScanDirectory(d, chain, pattern, recursive, verbose))
1192  ret = true;
1193  delete d;
1194  }
1195  continue;
1196  }
1197 
1198  // If this is not a root file, ignore
1199  if (!name.EndsWith(".root") &&
1200  !name.EndsWith(".zip",TString::kIgnoreCase)) {
1201  if (verbose) ::Info("ScanDirectory", "Ignoring non-ROOT or ZIP %s",
1202  name.Data());
1203  continue;
1204  }
1205 
1206  // If this file does not contain AliESDs, ignore
1207  if (!name.Contains(wild)) {
1208  if (verbose) ::Info("ScanDirectory", "%s does not match %s",
1209  name.Data(), pattern.Data());
1210  continue;
1211  }
1212  // Add
1213  if (verbose)
1214  ::Info("::ScanDirectory", "Candidate %s", full.Data());
1215  toAdd.Add(new TObjString(full));
1216  }
1217  TIter nextAdd(&toAdd);
1218  TObjString* s = 0;
1219  Int_t added = 0;
1220  while ((s = static_cast<TObjString*>(nextAdd()))) {
1221  // Info("ChainBuilder::ScanDirectory",
1222  // "Adding %s", s->GetString().Data());
1223  TString fn = s->GetString();
1224  if (!CheckFile(fn, chain, true/*verbose*/)) continue;
1225 
1226  added++;
1227  }
1228  if (added > 0) ret = true;
1229  if (verbose) ::Info("ScanDirectory", "Added %d files", added);
1230 
1231  gSystem->ChangeDirectory(oldDir);
1232  return ret;
1233  }
1234  //------------------------------------------------------------------
1244  static Bool_t CheckFile(const TString& path,
1245  TChain* chain,
1246  Bool_t verbose)
1247  {
1248  // Info("", "Checking %s", path.Data());
1249  TString fn = path;
1250 
1251  gSystem->RedirectOutput("/dev/null", "w");
1252  TFile* test = TFile::Open(fn, "READ");
1253  gSystem->RedirectOutput(0);
1254  if (!test) {
1255  ::Warning("CheckFile", "Failed to open %s", fn.Data());
1256  return false;
1257  }
1258 
1259  Bool_t ok = false;
1260  TObject* o = test->Get(chain->GetName());
1261  TTree* t = dynamic_cast<TTree*>(o);
1262  TFileCollection* c = dynamic_cast<TFileCollection*>(o);
1263  if (t) {
1264  Int_t n = t->GetEntries();
1265  test->Close();
1266  chain->Add(fn, n);
1267  ok = true;
1268  if (verbose) ::Info("CheckFile", "Added file %s (%d)", fn.Data(), n);
1269  } else if (c) {
1270  chain->AddFileInfoList(c->GetList());
1271  ok = true;
1272  if (verbose) ::Info("CheckFile", "Added file collection %s", fn.Data());
1273  } else {
1274  // Let's try to find a TFileCollection
1275  TList* l = test->GetListOfKeys();
1276  TIter next(l);
1277  TKey* k = 0;
1278  while ((k = static_cast<TKey*>(next()))) {
1279  TString cl(k->GetClassName());
1280  if (!cl.EqualTo("TFileCollection")) continue;
1281  c = dynamic_cast<TFileCollection*>(k->ReadObj());
1282  if (!c) {
1283  ::Warning("CheckFile", "Returned collection invalid");
1284  continue;
1285  }
1286  // Info("", "Adding file collection");
1287  chain->AddFileInfoList(c->GetList());
1288  ok = true;
1289  if (verbose) ::Info("CheckFile", "Added file collection %s", fn.Data());
1290  }
1291  test->Close();
1292  }
1293 
1294  if (!ok)
1295  ::Warning("CheckFile",
1296  "The file %s does not contain the tree %s "
1297  "or a file collection",
1298  path.Data(), chain->GetName());
1299  return ok;
1300  }
1301 
1310  static Bool_t Run(Long64_t maxEvents=-1,
1311  const char* title="")
1312  {
1313  TStopwatch timer;
1314  timer.Start();
1315  TChain* chain = MakeChain("tuple");
1316  if (!chain) {
1317  ::Error("Run", "No chain!");
1318  return false;
1319  }
1320 
1321  TupleSelector* s = new TupleSelector(title);
1322  Int_t status= chain->Process(s, "", maxEvents);
1323  timer.Print();
1324  return status >= 0;
1325  }
1335  static Bool_t Proof(Long64_t maxEvents,
1336  const char* opt="",
1337  const char* title="")
1338  {
1339  TStopwatch timer;
1340  timer.Start();
1341  TProof::Reset("lite:///?workers=8");
1342  TProof::Open("lite:///?workers=8");
1343  gProof->ClearCache();
1344  TString ali = gSystem->ExpandPathName("$(ALICE_PHYSICS)");
1345  TString fwd = ali + "/PWGLF/FORWARD/analysis2";
1346  gProof->AddIncludePath(Form("%s/include", ali.Data()));
1347  gProof->Load(Form("%s/scripts/LoadLibs.C",fwd.Data()), true);
1348  gProof->Exec("LoadLibs()");
1349  // gROOT->ProcessLine("gProof->SetLogLevel(5);");
1350  gProof->Load(Form("%s/scripts/TupleSelector.C+%s", fwd.Data(), opt),true);
1351 
1352  TDSet* dataset = MakeDataSet("tuple");
1353  if (!dataset) {
1354  ::Error("Proof", "No dataset");
1355  return false;
1356  }
1357 
1358  TupleSelector* s = new TupleSelector(title);
1359  gProof->AddFeedback("rings");
1360  gProof->Process(dataset, s, "", maxEvents);
1361  timer.Print();
1362  return true; // status >= 0;
1363  }
1364 
1365  /* @} */
1367 
1368 private:
1372  TupleSelector(const TupleSelector&); // {}
1378  TupleSelector& operator=(const TupleSelector&); // { return *this; }
1379 
1380 };
1381 
1382 
1383 
1384 #endif
1385 // Local Variables:
1386 // mode: C++
1387 // End:
1388 //
1389 // EOF
1390 //
Int_t color[]
Spectra(const char *name, const char *title, Color_t color, const TArrayD &bins)
void Init(TTree *tree)
virtual ~Spectra()
TString fName
void Store(TDirectory *dir)
double Double_t
Definition: External.C:58
Double_t BetaGamma() const
const char * title
Definition: MakeQAPdf.C:26
UShort_t fD
const char * GetName() const
TH1 * Hist(Bool_t primary) const
Spectra(const Spectra &o)
TClonesArray * fHits
Must not be persistent.
long long Long64_t
Definition: External.C:43
Long64_t Merge(TCollection *list)
ClassDef(Monitor, 1)
UShort_t Index(UShort_t d, Char_t r) const
void Fill(AliFMDMCTrackELoss::Hit *hit)
TObjArray * fSpectra
TSystem * gSystem
Ring * Get(UShort_t d, Char_t r)
TTree * fTree
Must not be persistent.
char Char_t
Definition: External.C:18
void ls(Option_t *option="") const
TList * list
void Store(TDirectory *dir)
TCanvas * c
Definition: TestFitELoss.C:172
static TChain * MakeChain(const TString &src=".", Bool_t recursive=false, Bool_t verbose=false)
void SlaveBegin(TTree *tree)
static TDSet * MakeDataSet(const TString &src=".", Bool_t recursive=false, Bool_t verbose=false)
const char * GetName() const
ClassDef(Ring, 2)
const Spectra * Get(UInt_t which) const
static Bool_t ScanDirectory(TSystemDirectory *dir, TChain *chain, const TString &pattern, Bool_t recursive, Bool_t verbose)
int Int_t
Definition: External.C:63
static Color_t RingColor(UShort_t d, Char_t r)
Long64_t Merge(TCollection *list)
unsigned int UInt_t
Definition: External.C:33
TH1 * Hist(Bool_t primary, UShort_t which) const
void SlaveTerminate()
const char * GetName() const
Definition: External.C:212
Bool_t Add(const Spectra *oth)
virtual void Stack(THStack *stack, TLegend *leg)
TString fName
Definition: TupleSelector.C:99
static TChain * DoMakeChain(const TString &src=".", Bool_t recursive=false, Bool_t verbose=false)
Spectra & operator=(const Spectra &o)
ClassDef(TupleSelector, 2)
static TObject * GetChainOrDataSet(const char *file)
ClassDef(Spectra, 1)
static Bool_t Run(Long64_t maxEvents=-1, const char *title="")
static Bool_t Proof(Long64_t maxEvents, const char *opt="", const char *title="")
void Begin(TTree *tree)
const char * fwd
void Fill(Bool_t primary, Double_t x)
virtual ~Monitor()
ClassDef(Type, 1)
TupleSelector(const char *name="")
Not persistent.
void Fill(AliFMDMCTrackELoss::Hit *hit)
TString fName
Ring(const Ring &r)
Type(const char *name, const char *title, Color_t color)
const char * GetName() const
TupleSelector & operator=(const TupleSelector &)
TFile * file
void Draw(Option_t *opt="")
TH1 * fSecondary
static void MakeLogScale(Int_t nBins, Int_t minOrder, Int_t maxOrder, TArrayD &bins)
TObjArray * fRings
Not persistent.
unsigned short UShort_t
Definition: External.C:28
Int_t Version() const
const Ring * Get(UShort_t d, Char_t r) const
const char Option_t
Definition: External.C:48
Bool_t Process(Long64_t entry)
Ring(UShort_t d, Char_t r)
Int_t fI
Must not be persistent.
void PlotOne(UShort_t which, const char *title, TDirectory *dir=0, Option_t *opt="")
TH1 * fPrimary
void Feedback(TList *objs)
bool Bool_t
Definition: External.C:53
Ring & operator=(const Ring &r)
Bool_t Notify()
void Store(TDirectory *dir)
static Bool_t CheckFile(const TString &path, TChain *chain, Bool_t verbose)
Long64_t GetStatus() const
void ls(Option_t *option="") const
void SetName(const char *name)
Spectra * Get(UInt_t which)
Definition: External.C:196
Char_t fR