]> git.uio.no Git - u/mrichter/AliRoot.git/blob - test/gun/test.C
added plots for vertex position and extended axis ranges for 7TeV data
[u/mrichter/AliRoot.git] / test / gun / test.C
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 "AliESDEvent.h"
19 #include "AliRunLoader.h"
20 #include "AliRun.h"
21 #include "AliStack.h"
22
23 #endif
24
25 void test(const char * sdir =".") {
26
27   TStopwatch timer;
28   timer.Start();
29
30   TString name;
31
32   // Signal file and tree
33   name = sdir;
34   name += "/AliESDs.root";
35   TFile * fSig = TFile::Open(name.Data());
36   TTree * tSig = (TTree*)fSig->Get("esdTree");
37
38   AliESDEvent * esdSig = new AliESDEvent();// The signal ESD object is put here
39   esdSig->ReadFromTree(tSig);
40
41   // Run loader (signal events)
42   name = sdir;
43   name += "/galice.root";
44   AliRunLoader* rlSig = AliRunLoader::Open(name.Data());
45
46   // gAlice
47   rlSig->LoadgAlice();
48   gAlice = rlSig->GetAliRun();
49
50   // Now load kinematics and event header
51   rlSig->LoadKinematics();
52   rlSig->LoadHeader();
53
54   // Loop on events: check that MC and data contain the same number of events
55   Long64_t nevSig = rlSig->GetNumberOfEvents();
56
57   cout << nevSig << " signal events" << endl;
58
59   Int_t lab[3]; // Labels from TOF
60   Double_t mom[3]; // Track momentum
61
62   for (Int_t iev=0; iev<nevSig; iev++) {
63     cout << "---------- Signal event ----------" << iev << endl;
64
65     // Get signal ESD
66     tSig->GetEntry(iev);
67
68     // Particle stack
69     rlSig->GetEvent(iev);
70     AliStack * stackSig = rlSig->Stack();
71     stackSig->DumpPStack();
72     Int_t nPartSig = stackSig->GetNtrack();
73
74     Int_t nrec = esdSig->GetNumberOfTracks();
75     cout << nrec << " reconstructed tracks" << endl;
76     for(Int_t irec=0; irec<nrec; irec++) {
77       AliESDtrack * track = esdSig->GetTrack(irec);
78       cout << "Labels:" << endl;
79       cout << "Global: "<< track->GetLabel() << endl;
80       cout << "ITS: "<< track->GetITSLabel() << endl;
81       cout << "TPC: "<< track->GetTPCLabel() << endl;
82       cout << "TRD: "<< track->GetTRDLabel() << endl;
83       track->GetTOFLabel(lab);
84       cout << "TOF: "<< lab[0] <<" "<< lab[1] <<" "<< lab[2] << endl;
85       UInt_t label = TMath::Abs(track->GetLabel());
86       if (label>=10000000) {
87         // Underlying event. 10000000 is the
88         // value of fkMASKSTEP in AliRunDigitizer
89         cout <<"Strange, there should be no underlying event"<<endl;
90       }
91       else {
92         if (label>=nPartSig) {
93           cout <<"Strange, label outside the range"<< endl;
94           continue;
95         }
96         TParticle * part = stackSig->Particle(label);
97         if(part) part->Print();
98         track->GetPxPyPz(mom);
99         cout <<"Momentum: "<< mom[0] <<" "<< mom[1] <<" "<< mom[2] <<endl; 
100       }
101
102     }
103
104   }
105
106   fSig->Close();
107
108   timer.Stop();
109   timer.Print();
110 }