]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TRD/AliTRDanalyzeHits.C
Update of the tracking by Sergei
[u/mrichter/AliRoot.git] / TRD / AliTRDanalyzeHits.C
CommitLineData
793ff80c 1Int_t AliTRDanalyzeHits()
2{
3 //
4 // Analyzes the hits and fills QA-histograms
5 //
6
7 Int_t rc = 0;
8
9 if (!gAlice) {
10 cout << "<AliTRDanalyzeHits> No AliRun object found" << endl;
11 rc = 1;
12 return rc;
13 }
14 gAlice->GetEvent(0);
15
16 // Get the pointer to the TRD detector
17 AliDetector *TRD = gAlice->GetDetector("TRD");
18 if (!TRD) {
19 cout << "<AliTRDanalyzeHits> No TRD detector found" << endl;
20 rc = 2;
21 return rc;
22 }
23
24 // Define the histograms
25 TH1F *hQdedx = new TH1F("hQdedx","Charge dedx-hits",100,0.0,1000.0);
26 TH1F *hQtr = new TH1F("hQtr" ,"Charge TR-hits" ,100,0.0,1000.0);
27 TH2F *hZY = new TH2F("hZY" ,"Y vs Z",50,-100.0,100.0,40,290.0,370.0);
28 TH2F *hXZ = new TH2F("hXZ" ,"Z vs X",50,-100.0,100.0,50,-100.0,100.0);
29
30 // Get the pointer hit tree
31 TTree *HitTree = gAlice->TreeH();
32 if (!HitTree) {
33 cout << "<AliTRDanalyzeHits> No hit tree found" << endl;
34 rc = 3;
35 return rc;
36 }
37
38 Int_t countHits = 0;
39 Int_t nBytes = 0;
40
41 // Get the number of entries in the hit tree
42 // (Number of primary particles creating a hit somewhere)
43 Int_t nTrack = (Int_t) HitTree->GetEntries();
44 cout << "<AliTRDanalyzeHits> Found " << nTrack
45 << " primary particles with hits" << endl;
46
47 // Loop through all entries in the tree
48 for (Int_t iTrack = 0; iTrack < nTrack; iTrack++) {
49
50 gAlice->ResetHits();
51 nBytes += HitTree->GetEvent(iTrack);
52
53 // Get the number of hits in the TRD created by this particle
54 Int_t nHit = TRD->Hits()->GetEntriesFast();
55 //cout << "<AliTRDanalyzeHits> Found " << nHit
56 // << " hits for primary particle " << iTrack << endl;
57
58 // Loop through the TRD hits
59 for (Int_t iHit = 0; iHit < nHit; iHit++) {
60
61 countHits++;
62
63 AliTRDhit *hit = (AliTRDhit *) TRD->Hits()->UncheckedAt(iHit);
64
65 Float_t x = hit->X();
66 Float_t y = hit->Y();
67 Float_t z = hit->Z();
68 Float_t q = hit->GetCharge();
69
70 if (q > 0)
71 hQdedx->Fill(q);
72 else
73 hQtr->Fill(TMath::Abs(q));
74
75 hZY->Fill(z,y);
76 hXZ->Fill(x,z);
77
78 }
79
80 }
81
82 cout << "<AliTRDanalyzeHits> Found " << countHits << " hits in total" << endl;
83
84 TCanvas *cHits = new TCanvas("cHits","AliTRDanalyzeHits",50,50,600,600);
85 cHits->Divide(2,2);
86 cHits->cd(1);
87 hXZ->Draw("COL");
88 cHits->cd(2);
89 hZY->Draw("COL");
90 cHits->cd(3);
91 gPad->SetLogy();
92 hQdedx->Draw();
93 cHits->cd(4);
94 gPad->SetLogy();
95 hQtr->Draw();
96
97 return rc;
98
99}