]>
Commit | Line | Data |
---|---|---|
23abf4db | 1 | #include <TH1F.h> |
2 | #include <TMath.h> | |
3 | #include <TObjArray.h> | |
4 | #include <TProfile.h> | |
5 | ||
6 | #include "AliTRDcluster.h" | |
7 | #include "AliTRDseedV1.h" | |
8 | #include "AliTRDtrackV1.h" | |
9 | #include "TTreeStream.h" | |
10 | ||
11 | #include "AliTRDtrackInfo/AliTRDtrackInfo.h" | |
12 | #include "AliTRDcheckDetector.h" | |
13 | ||
14 | //////////////////////////////////////////////////////////////////////////// | |
15 | // // | |
16 | // Reconstruction QA // | |
17 | // // | |
18 | // Task doing basic checks for tracking and detector performance // | |
19 | // // | |
20 | // Authors: // | |
21 | // Anton Andronic <A.Andronic@gsi.de> // | |
22 | // Markus Fasel <M.Fasel@gsi.de> // | |
23 | // // | |
24 | //////////////////////////////////////////////////////////////////////////// | |
25 | ||
26 | //_______________________________________________________ | |
27 | AliTRDcheckDetector::AliTRDcheckDetector(): | |
2b468513 | 28 | AliTRDrecoTask("detChecker", "Utility to check TRD detector signal") |
29 | ,fPHSdetector(0x0) | |
30 | ,fPHSsector(0x0) | |
31 | ,fQCLdetector(0x0) | |
32 | ,fQCLsector(0x0) | |
33 | ,fQTdetector(0x0) | |
34 | ,fQTsector(0x0) | |
23abf4db | 35 | { |
2b468513 | 36 | // |
37 | // Default constructor | |
38 | // | |
23abf4db | 39 | } |
40 | ||
41 | //_______________________________________________________ | |
42 | AliTRDcheckDetector::~AliTRDcheckDetector(){ | |
2b468513 | 43 | // |
44 | // Destructor | |
45 | // | |
46 | if(fPHSdetector) delete fPHSdetector; | |
47 | if(fPHSsector) delete fPHSsector; | |
48 | if(fQCLdetector) delete fQCLdetector; | |
49 | if(fQCLsector) delete fQCLsector; | |
50 | if(fQTdetector) delete fQTdetector; | |
51 | if(fQTsector) delete fQTsector; | |
23abf4db | 52 | } |
53 | ||
54 | //_______________________________________________________ | |
55 | void AliTRDcheckDetector::CreateOutputObjects(){ | |
2b468513 | 56 | // |
57 | // Create Output Objects | |
58 | // | |
59 | fContainer = new TObjArray(); | |
60 | TH1F *hNtrks = new TH1F("hNtrks", "Number of Tracks per event", 100, 0, 100); | |
61 | TH1F *hNcls = new TH1F("hNcls", "Nr. of clusters per tracklet", 30, 0, 30); | |
23abf4db | 62 | TH1F *hNtls = new TH1F("hNtls", "Nr. tracklets per track", 10, 0, 10); |
63 | TH1F *hNclTls = new TH1F("hNclTls","Mean Number of clusters per tracklet", 30, 0, 30); | |
64 | TH1F *hChi2 = new TH1F("hChi2", "Chi2", 200, 0, 20); | |
65 | TH1F *hChi2n = new TH1F("hChi2n", "Norm. Chi2 (tracklets)", 50, 0, 5); | |
66 | TH1F *hSM = new TH1F("hSM", "Track Counts in Supermodule", 18, 0., 18.); | |
67 | TH1F *hQcl = new TH1F("hQcl1", "Cluster charge", 200, -1200, 1200); | |
68 | TProfile *hPH = new TProfile("hPH", "Average PH", 31, -0.5, 30.5); | |
69 | TH1F * hQT = new TH1F("hQT", "Total Charge Deposit", 6000, 0, 6000); | |
70 | // Register Histograms | |
71 | fContainer->Add(hNtrks); | |
2b468513 | 72 | fContainer->Add(hNcls); |
73 | fContainer->Add(hNtls); | |
74 | fContainer->Add(hNclTls); | |
75 | fContainer->Add(hChi2); | |
76 | fContainer->Add(hChi2n); | |
77 | fContainer->Add(hSM); | |
78 | fContainer->Add(hPH); | |
79 | fContainer->Add(hQcl); | |
80 | fContainer->Add(hQT); | |
23abf4db | 81 | |
2b468513 | 82 | // Initialise the PHS, cluster charge and total charge deposit histograms for each detector |
83 | fPHSdetector = new TObjArray(); | |
84 | fQCLdetector = new TObjArray(); | |
85 | fQTdetector = new TObjArray(); | |
86 | for(Int_t idet = 0; idet < kNDetectors; idet++){ | |
87 | fPHSdetector->Add(new TProfile(Form("hPHSdet%d", idet), Form("Average Pulse Height in Detector %d", idet), 31, -0.5, 30.5)); | |
88 | fQCLdetector->Add(new TH1F(Form("hQd%d",idet),Form("Qd%d",idet), 100,0,5000)); | |
89 | fQTdetector->Add(new TH1F(Form("hQTd%d", idet), Form("Qtotal%d", idet), 6000, 0, 6000)); | |
90 | } | |
91 | fContainer->Add(fPHSdetector); | |
92 | fContainer->Add(fQCLdetector); | |
93 | fContainer->Add(fQTdetector); | |
94 | ||
95 | // Initialise the PHS, cluster charge and total charge deposit histograms for each sector | |
96 | fPHSsector = new TObjArray(); | |
97 | fQCLsector = new TObjArray(); | |
98 | fQTsector = new TObjArray(); | |
99 | for(Int_t isec = 0; isec < kNSectors; isec++){ | |
100 | fPHSsector->Add(new TProfile(Form("hPHSsec%d", isec), Form("Average Pulse Height in Sector %d", isec), 31, -0.5, 30.5)); | |
101 | fQCLsector->Add(new TH1F(Form("hQd%d",isec),Form("Qd%d",isec), 100,0,5000)); | |
102 | fQTsector->Add(new TH1F(Form("hQTd%d", isec), Form("Qtotal%d", isec), 6000, 0, 6000)); | |
103 | } | |
104 | fContainer->Add(fPHSsector); | |
105 | fContainer->Add(fQCLsector); | |
106 | fContainer->Add(fQTsector); | |
23abf4db | 107 | } |
108 | ||
109 | //_______________________________________________________ | |
110 | void AliTRDcheckDetector::Exec(Option_t *){ | |
2b468513 | 111 | // |
112 | // Execution function | |
113 | // Filling TRD quality histos | |
114 | // | |
115 | Int_t nTracks = 0; // Count the number of tracks per event | |
116 | AliTRDtrackInfo *fTrackInfo = 0x0; | |
117 | AliTRDtrackV1 *fTRDtrack = 0x0; | |
118 | AliTRDseedV1 *fTracklet = 0x0; | |
119 | AliTRDcluster *fTRDcluster = 0x0; | |
120 | for(Int_t iti = 0; iti < fTracks->GetEntriesFast(); iti++){ | |
121 | fTrackInfo = dynamic_cast<AliTRDtrackInfo *>(fTracks->UncheckedAt(iti)); | |
122 | if(!fTrackInfo || !(fTRDtrack = fTrackInfo->GetTRDtrack())) continue; | |
123 | Int_t nclusters = fTRDtrack->GetNumberOfClusters(); | |
124 | Int_t ntracklets = fTRDtrack->GetNumberOfTracklets(); | |
125 | Float_t chi2 = fTRDtrack->GetChi2(); | |
126 | // Fill Histograms | |
127 | dynamic_cast<TH1F *>(fContainer->UncheckedAt(kNclustersHist))->Fill(nclusters); | |
128 | dynamic_cast<TH1F *>(fContainer->UncheckedAt(kNtrackletsHist))->Fill(ntracklets); | |
129 | dynamic_cast<TH1F *>(fContainer->UncheckedAt(kChi2))->Fill(chi2); | |
130 | dynamic_cast<TH1F *>(fContainer->UncheckedAt(kChi2Normalized))->Fill(chi2/static_cast<Float_t>(ntracklets)); | |
131 | // now loop over single tracklets | |
132 | if(ntracklets > 2){ | |
133 | Int_t sector = -1; | |
134 | for(Int_t ilayer = 0; ilayer < kNLayers; ilayer++){ | |
135 | if(!(fTracklet = fTRDtrack->GetTracklet(ilayer)) || !fTracklet->IsOK()) continue; | |
136 | Float_t nClustersTracklet = fTracklet->GetN(); | |
137 | dynamic_cast<TH1F *>(fContainer->UncheckedAt(kNclusterTrackletHist))->Fill(nClustersTracklet); | |
138 | if(nClustersTracklet){ | |
139 | Float_t Qtot = 0; | |
140 | Int_t detector = -1; | |
141 | for(Int_t itb = 0; itb < kNTimebins; itb++){ | |
142 | if(!(fTRDcluster = fTracklet->GetClusters(itb))) continue; | |
143 | Int_t localtime = fTRDcluster->GetLocalTimeBin(); | |
144 | Double_t clusterCharge = fTRDcluster->GetQ(); | |
145 | detector = fTRDcluster->GetDetector(); | |
146 | sector = static_cast<Int_t>(detector/kNDetectorsSector); | |
147 | Double_t absolute_charge = TMath::Abs(clusterCharge); | |
148 | Qtot += absolute_charge; | |
149 | // Fill Histograms | |
150 | // 1st. PHS | |
151 | dynamic_cast<TProfile *>(fContainer->UncheckedAt(kPulseHeight))->Fill(localtime, absolute_charge); | |
152 | dynamic_cast<TProfile *>(fPHSdetector->UncheckedAt(detector))->Fill(localtime, absolute_charge); | |
153 | dynamic_cast<TProfile *>(fPHSsector->UncheckedAt(sector))->Fill(localtime, absolute_charge); | |
154 | // 2nd. cluster charge | |
155 | dynamic_cast<TH1F *>(fContainer->UncheckedAt(kClusterCharge))->Fill(absolute_charge); | |
156 | dynamic_cast<TH1F *>(fQCLdetector->UncheckedAt(detector))->Fill(absolute_charge); | |
157 | dynamic_cast<TH1F *>(fQCLsector->UncheckedAt(sector))->Fill(absolute_charge); | |
158 | } | |
159 | // Fill total cluster charge histograms | |
160 | dynamic_cast<TH1F *>(fContainer->UncheckedAt(kChargeDeposit))->Fill(Qtot); | |
161 | dynamic_cast<TH1F *>(fQTdetector->UncheckedAt(detector))->Fill(Qtot); | |
162 | dynamic_cast<TH1F *>(fQTsector->UncheckedAt(sector))->Fill(Qtot); | |
163 | } | |
164 | } | |
165 | dynamic_cast<TH1F *>(fContainer->UncheckedAt(kNTracksSectorHist))->Fill(sector); | |
166 | } | |
167 | nTracks++; | |
168 | } | |
169 | dynamic_cast<TH1F *>(fContainer->UncheckedAt(kNTracksEventHist))->Fill(nTracks); | |
170 | PostData(0, fContainer); | |
23abf4db | 171 | } |
172 | ||
173 | //_______________________________________________________ | |
174 | void AliTRDcheckDetector::Terminate(Option_t *){ | |
2b468513 | 175 | // |
176 | // Terminate function | |
177 | // | |
23abf4db | 178 | } |
179 |