]>
Commit | Line | Data |
---|---|---|
1 | /************************************************************************** | |
2 | * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. * | |
3 | * * | |
4 | * Author: The ALICE Off-line Project. * | |
5 | * Contributors are mentioned in the code where appropriate. * | |
6 | * * | |
7 | * Permission to use, copy, modify and distribute this software and its * | |
8 | * documentation strictly for non-commercial purposes is hereby granted * | |
9 | * without fee, provided that the above copyright notice appears in all * | |
10 | * copies and that both the copyright notice and this permission notice * | |
11 | * appear in the supporting documentation. The authors make no claims * | |
12 | * about the suitability of this software for any purpose. It is * | |
13 | * provided "as is" without express or implied warranty. * | |
14 | **************************************************************************/ | |
15 | ||
16 | /* $Id$ */ | |
17 | ||
18 | /////////////////////////////////////////////////////////////////////////////// | |
19 | // // | |
20 | // This class creates and fills the monitor histograms for the ITS // | |
21 | // // | |
22 | /////////////////////////////////////////////////////////////////////////////// | |
23 | ||
24 | ||
25 | #include "AliMonitorITS.h" | |
26 | #include "AliMonitorHisto.h" | |
27 | #include "AliITSgeom.h" | |
28 | #include "AliITSRawStreamSDD.h" | |
29 | #include "AliITSRawStreamSSD.h" | |
30 | #include "AliITSclusterV2.h" | |
31 | #include "AliITStrackV2.h" | |
32 | #include "AliRunLoader.h" | |
33 | #include <TFolder.h> | |
34 | #include <TTree.h> | |
35 | ||
36 | ||
37 | ClassImp(AliMonitorITS) | |
38 | ||
39 | ||
40 | //_____________________________________________________________________________ | |
41 | AliMonitorITS::AliMonitorITS(AliITSgeom* geom) | |
42 | { | |
43 | // create a ITS monitor object with the given geometry | |
44 | ||
45 | fGeom = geom; | |
46 | } | |
47 | ||
48 | ||
49 | //_____________________________________________________________________________ | |
50 | AliMonitorITS::AliMonitorITS(const AliMonitorITS& monitor) : | |
51 | AliMonitor(monitor) | |
52 | { | |
53 | Fatal("AliMonitorITS", "copy constructor not implemented"); | |
54 | } | |
55 | ||
56 | //_____________________________________________________________________________ | |
57 | AliMonitorITS& AliMonitorITS::operator = (const AliMonitorITS& /*monitor*/) | |
58 | { | |
59 | Fatal("operator =", "assignment operator not implemented"); | |
60 | return *this; | |
61 | } | |
62 | ||
63 | //_____________________________________________________________________________ | |
64 | void AliMonitorITS::CreateHistos(TFolder* folder) | |
65 | { | |
66 | // create the ITS monitor histograms | |
67 | ||
68 | fFolder = folder->AddFolder("ITS", "ITS"); | |
69 | ||
70 | fSDDDigitsCharge = CreateHisto1("SDDDigitsCharge", | |
71 | "charge distribution of SDD digits", | |
72 | 128, 0, 256, "charge", "#Delta N/N", | |
73 | AliMonitorHisto::kNormEvents); | |
74 | ||
75 | fSSDDigitsCharge = CreateHisto1("SSDDigitsCharge", | |
76 | "charge distribution of SSD digits", | |
77 | 100, 0, 200, "charge", "#Delta N/N", | |
78 | AliMonitorHisto::kNormEvents); | |
79 | ||
80 | fSDDClustersCharge = CreateHisto1("SDDClustersCharge", | |
81 | "charge distribution of SDD clusters", | |
82 | 100, 0, 200, "charge", "#Delta N/N", | |
83 | AliMonitorHisto::kNormEvents); | |
84 | ||
85 | fSSDClustersCharge = CreateHisto1("SSDClustersCharge", | |
86 | "charge distribution of SSD clusters", | |
87 | 100, 0, 400, "charge", "#Delta N/N", | |
88 | AliMonitorHisto::kNormEvents); | |
89 | ||
90 | Int_t nModulesSPD = fGeom->GetLastSPD() - fGeom->GetStartSPD() + 1; | |
91 | fSPDNClustersVsModule = CreateHisto1("SPDNClustersVsModule", | |
92 | "mean number of clusters per SPD module", | |
93 | nModulesSPD, -0.5, nModulesSPD-0.5, | |
94 | "module", "<N_{clusters}>", | |
95 | AliMonitorHisto::kNormEvents); | |
96 | ||
97 | Int_t nModulesSDD = fGeom->GetLastSDD() - fGeom->GetStartSDD() + 1; | |
98 | fSDDNClustersVsModule = CreateHisto1("SDDNClustersVsModule", | |
99 | "mean number of clusters per SDD module", | |
100 | nModulesSDD, -0.5, nModulesSDD-0.5, | |
101 | "module", "<N_{clusters}>", | |
102 | AliMonitorHisto::kNormEvents); | |
103 | ||
104 | Int_t nModulesSSD = fGeom->GetLastSSD() - fGeom->GetStartSSD() + 1; | |
105 | fSSDNClustersVsModule = CreateHisto1("SSDNClustersVsModule", | |
106 | "mean number of clusters per SSD module", | |
107 | nModulesSSD, -0.5, nModulesSSD-0.5, | |
108 | "module", "<N_{clusters}>", | |
109 | AliMonitorHisto::kNormEvents); | |
110 | ||
111 | Int_t nLayer = fGeom->GetNlayers(); | |
112 | fNClustersVsLayer = CreateHisto1("NClustersVsLayer", | |
113 | "mean number of clusters per layer", | |
114 | nLayer, 0.5, nLayer+0.5, | |
115 | "layer", "<N_{clusters}>", | |
116 | AliMonitorHisto::kNormEvents); | |
117 | ||
118 | fNTracks = CreateHisto1("NTracks", "number of tracks per event", | |
119 | 300, 0, 30000, "N_{tracks}", "#Delta N/N", | |
120 | AliMonitorHisto::kNormEvents); | |
121 | ||
122 | fNTracksITSTPC = CreateHisto2("NTracksTPCITS", | |
123 | "correlation of number of ITS and TPC tracks per event", | |
124 | 30, 0, 30000, 30, 0, 30000, | |
125 | "N_{tracks}(ITS)", "N_{tracks}(TPC)", | |
126 | "#Delta N/N", AliMonitorHisto::kNormEvents); | |
127 | ||
128 | fTrackPt = CreateHisto1("TrackPt", "pt distribution of tracks", | |
129 | 90, 0, 3, "p_{t} [GeV/c]", "#Delta N/N", | |
130 | AliMonitorHisto::kNormEntries); | |
131 | ||
132 | fTrackEta = CreateHisto1("TrackEta", "eta distribution of tracks", | |
133 | 100, -2, 2, "#eta", "#Delta N/N", | |
134 | AliMonitorHisto::kNormEntries); | |
135 | ||
136 | fTrackPhi = CreateHisto1("TrackPhi", "phi distribution of tracks", | |
137 | 120, 0, 360, "#phi [#circ]", "#Delta N/N", | |
138 | AliMonitorHisto::kNormEntries); | |
139 | ||
140 | fTrackDEdxVsP = CreateHisto2("TrackDEdxVsP", "dE/dx of tracks", | |
141 | 100, 0, 3, 100, 0, 150, | |
142 | "p [GeV/c]", "dE/dx", "#Delta N/N", | |
143 | AliMonitorHisto::kNormEntries); | |
144 | } | |
145 | ||
146 | ||
147 | //_____________________________________________________________________________ | |
148 | void AliMonitorITS::FillHistos(AliRunLoader* runLoader, | |
149 | AliRawReader* rawReader) | |
150 | { | |
151 | // fill the ITS monitor histogrms | |
152 | ||
153 | rawReader->Reset(); | |
154 | AliITSRawStreamSDD inputSDD(rawReader); | |
155 | while (inputSDD.Next()) { | |
156 | fSDDDigitsCharge->Fill(inputSDD.GetSignal()); | |
157 | } | |
158 | ||
159 | rawReader->Reset(); | |
160 | AliITSRawStreamSSD inputSSD(rawReader); | |
161 | while (inputSSD.Next()) { | |
162 | fSSDDigitsCharge->Fill(inputSSD.GetSignal()); | |
163 | } | |
164 | ||
165 | AliLoader* itsLoader = runLoader->GetLoader("ITSLoader"); | |
166 | if (!itsLoader) return; | |
167 | ||
168 | ||
169 | itsLoader->LoadRecPoints(); | |
170 | TTree* clustersTree = itsLoader->TreeR(); | |
171 | if (!clustersTree) return; | |
172 | TClonesArray* clusters = new TClonesArray("AliITSclusterV2"); | |
173 | clustersTree->SetBranchAddress("Clusters", &clusters); | |
174 | ||
175 | for (Int_t iModule = 0; iModule < clustersTree->GetEntries(); iModule++) { | |
176 | clustersTree->GetEntry(iModule); | |
177 | Int_t iLayer, iLadder, iDetector; | |
178 | fGeom->GetModuleId(iModule, iLayer, iLadder, iDetector); | |
179 | ||
180 | for (Int_t j = 0; j < clusters->GetEntriesFast(); j++) { | |
181 | AliITSclusterV2* cluster = (AliITSclusterV2*) clusters->At(j); | |
182 | switch (iLayer) { | |
183 | case 1: case 2: { | |
184 | fSPDNClustersVsModule->Fill(iModule - fGeom->GetStartSPD()); | |
185 | break; | |
186 | } | |
187 | case 3: case 4: { | |
188 | fSDDClustersCharge->Fill(cluster->GetQ()); | |
189 | fSDDNClustersVsModule->Fill(iModule - fGeom->GetStartSDD()); | |
190 | break; | |
191 | } | |
192 | case 5: case 6: { | |
193 | fSSDClustersCharge->Fill(cluster->GetQ()); | |
194 | fSSDNClustersVsModule->Fill(iModule - fGeom->GetStartSSD()); | |
195 | break; | |
196 | } | |
197 | } | |
198 | fNClustersVsLayer->Fill(iLayer); | |
199 | } | |
200 | } | |
201 | ||
202 | delete clusters; | |
203 | itsLoader->UnloadRecPoints(); | |
204 | ||
205 | ||
206 | itsLoader->LoadTracks(); | |
207 | TTree* tracks = itsLoader->TreeT(); | |
208 | if (!tracks) return; | |
209 | AliITStrackV2* track = new AliITStrackV2; | |
210 | tracks->SetBranchAddress("tracks", &track); | |
211 | ||
212 | fNTracks->Fill(tracks->GetEntries()); | |
213 | for (Int_t i = 0; i < tracks->GetEntries(); i++) { | |
214 | tracks->GetEntry(i); | |
215 | fTrackPt->Fill(track->Pt()); | |
216 | fTrackEta->Fill(track->Eta()); | |
217 | fTrackPhi->Fill(track->Phi() * TMath::RadToDeg()); | |
218 | fTrackDEdxVsP->Fill(track->P(), track->GetdEdx()); | |
219 | } | |
220 | ||
221 | AliLoader* tpcLoader = runLoader->GetLoader("TPCLoader"); | |
222 | if (tpcLoader) { | |
223 | tpcLoader->LoadTracks(); | |
224 | TTree* tracksTPC = tpcLoader->TreeT(); | |
225 | if (tracksTPC) { | |
226 | fNTracksITSTPC->Fill(tracks->GetEntries(), tracksTPC->GetEntries()); | |
227 | } | |
228 | tpcLoader->UnloadTracks(); | |
229 | } | |
230 | ||
231 | delete track; | |
232 | itsLoader->UnloadTracks(); | |
233 | } |