]> git.uio.no Git - u/mrichter/AliRoot.git/blame - test/merge/test.C
Using AliTOFv6T0 instead of AliTOFv5T0
[u/mrichter/AliRoot.git] / test / merge / test.C
CommitLineData
3928b038 1// Usage in compiled mode
2// gSystem->SetIncludePath("-I$ROOTSYS/include -I$ALICE_ROOT/include");
3// gROOT->LoadMacro("test.C+");
4// test()
5
6#if !defined(__CINT__) || defined(__MAKECINT__)
7
8// Root include files
9#include <Riostream.h>
10#include <TFile.h>
11#include <TTree.h>
12#include <TBranch.h>
13#include <TStopwatch.h>
14#include <TObject.h>
15#include <TParticle.h>
16
17// AliRoot include files
18#include "AliESD.h"
19#include "AliRunLoader.h"
20#include "AliRun.h"
21#include "AliStack.h"
22
23#endif
24
25void test(const char * sdir ="signal",
26 const char * bdir ="backgr") {
27
28 TStopwatch timer;
29 timer.Start();
30
31 TString name;
32
33 // Signal file, tree, and branch
34 name = sdir;
35 name += "/AliESDs.root";
36 TFile * fSig = TFile::Open(name.Data());
37 TTree * tSig = (TTree*)fSig->Get("esdTree");
38 TBranch * bSig = tSig->GetBranch("ESD");
39
40 AliESD * esdSig = 0; // The signal ESD object is put here
41 bSig->SetAddress(&esdSig);
42
43 // Run loader (signal events)
44 name = sdir;
45 name += "/galice.root";
46 AliRunLoader* rlSig = AliRunLoader::Open(name.Data());
47
48 // Run loader (underlying events)
49 name = bdir;
50 name += "/galice.root";
51 AliRunLoader* rlUnd = AliRunLoader::Open(name.Data(),"Underlying");
52
53 // gAlice
54 rlSig->LoadgAlice();
55 rlUnd->LoadgAlice();
56 gAlice = rlSig->GetAliRun();
57
58 // Now load kinematics and event header
59 rlSig->LoadKinematics();
60 rlSig->LoadHeader();
61 rlUnd->LoadKinematics();
62 rlUnd->LoadHeader();
63
64 // Loop on events: check that MC and data contain the same number of events
65 Long64_t nevSig = rlSig->GetNumberOfEvents();
66 Long64_t nevUnd = rlUnd->GetNumberOfEvents();
67 Long64_t nSigPerUnd = nevSig/nevUnd;
68
69 cout << nevSig << " signal events" << endl;
70 cout << nevUnd << " underlying events" << endl;
71 cout << nSigPerUnd << " signal events per one underlying" << endl;
72
73 for (Int_t iev=0; iev<nevSig; iev++) {
74 cout << "Signal event " << iev << endl;
75 Int_t ievUnd = iev/nSigPerUnd;
76 cout << "Underlying event " << ievUnd << endl;
77
78 // Get signal ESD
79 bSig->GetEntry(iev);
80 // Get signal kinematics
81 rlSig->GetEvent(iev);
82 // Get underlying kinematics
83 rlUnd->GetEvent(ievUnd);
84
85 // Particle stack
86 AliStack * stackSig = rlSig->Stack();
87 Int_t nPartSig = stackSig->GetNtrack();
88 AliStack * stackUnd = rlUnd->Stack();
89 Int_t nPartUnd = stackUnd->GetNtrack();
90
91 Int_t nrec = esdSig->GetNumberOfTracks();
92 cout << nrec << " reconstructed tracks" << endl;
93 for(Int_t irec=0; irec<nrec; irec++) {
94 AliESDtrack * track = esdSig->GetTrack(irec);
95 UInt_t label = TMath::Abs(track->GetTPCLabel());
96 if (label>=10000000) {
97 // Underlying event. 10000000 is the
98 // value of fkMASKSTEP in AliRunDigitizer
99// cout << " Track from the underlying event" << endl;
100 label %=10000000;
101 if (label>=nPartUnd) continue;
102 TParticle * part = stackUnd->Particle(label);
103 if(part) part->Print();
104 }
105 else {
106 cout << " Track " << label << " from the signal event" << endl;
107 if (label>=nPartSig) {
108 cout <<"Strange, label outside the range "<< endl;
109 continue;
110 }
111 TParticle * part = stackSig->Particle(label);
112 if(part) part->Print();
113 }
114
115 }
116
117 }
118
119 fSig->Close();
120
121 timer.Stop();
122 timer.Print();
123}