]>
Commit | Line | Data |
---|---|---|
04fa961a | 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 "AliITSRawStreamSPD.h" | |
28 | #include "AliITSRawStreamSDD.h" | |
29 | #include "AliITSRawStreamSSD.h" | |
30 | #include "AliITSclusterV2.h" | |
31 | #include "AliITStrackV2.h" | |
32 | ||
33 | ||
34 | ClassImp(AliMonitorITS) | |
35 | ||
36 | ||
37 | //_____________________________________________________________________________ | |
38 | AliMonitorITS::AliMonitorITS(AliITSgeom* geom) | |
39 | { | |
40 | // create a ITS monitor object with the given geometry | |
41 | ||
42 | fGeom = geom; | |
43 | } | |
44 | ||
45 | ||
46 | //_____________________________________________________________________________ | |
47 | void AliMonitorITS::CreateHistos(TFolder* folder) | |
48 | { | |
49 | // create the ITS monitor histograms | |
50 | ||
51 | fFolder = folder->AddFolder("ITS", "ITS"); | |
52 | ||
53 | fSDDDigitsCharge = CreateHisto1("SDDDigitsCharge", | |
54 | "charge distribution of SDD digits", | |
55 | 128, 0, 256, "charge", "#Delta N/N", | |
56 | AliMonitorHisto::kNormEvents); | |
57 | ||
58 | fSSDDigitsCharge = CreateHisto1("SSDDigitsCharge", | |
59 | "charge distribution of SSD digits", | |
60 | 100, 0, 200, "charge", "#Delta N/N", | |
61 | AliMonitorHisto::kNormEvents); | |
62 | ||
63 | fSDDClustersCharge = CreateHisto1("SDDClustersCharge", | |
64 | "charge distribution of SDD clusters", | |
65 | 100, 0, 200, "charge", "#Delta N/N", | |
66 | AliMonitorHisto::kNormEvents); | |
67 | ||
68 | fSSDClustersCharge = CreateHisto1("SSDClustersCharge", | |
69 | "charge distribution of SSD clusters", | |
70 | 100, 0, 400, "charge", "#Delta N/N", | |
71 | AliMonitorHisto::kNormEvents); | |
72 | ||
73 | Int_t nModulesSPD = fGeom->GetLastSPD() - fGeom->GetStartSPD() + 1; | |
74 | fSPDNClustersVsModule = CreateHisto1("SPDNClustersVsModule", | |
75 | "mean number of clusters per SPD module", | |
76 | nModulesSPD, -0.5, nModulesSPD-0.5, | |
77 | "module", "<N_{clusters}>", | |
78 | AliMonitorHisto::kNormEvents); | |
79 | ||
80 | Int_t nModulesSDD = fGeom->GetLastSDD() - fGeom->GetStartSDD() + 1; | |
81 | fSDDNClustersVsModule = CreateHisto1("SDDNClustersVsModule", | |
82 | "mean number of clusters per SDD module", | |
83 | nModulesSDD, -0.5, nModulesSDD-0.5, | |
84 | "module", "<N_{clusters}>", | |
85 | AliMonitorHisto::kNormEvents); | |
86 | ||
87 | Int_t nModulesSSD = fGeom->GetLastSSD() - fGeom->GetStartSSD() + 1; | |
88 | fSSDNClustersVsModule = CreateHisto1("SSDNClustersVsModule", | |
89 | "mean number of clusters per SSD module", | |
90 | nModulesSSD, -0.5, nModulesSSD-0.5, | |
91 | "module", "<N_{clusters}>", | |
92 | AliMonitorHisto::kNormEvents); | |
93 | ||
94 | Int_t nLayer = fGeom->GetNlayers(); | |
95 | fNClustersVsLayer = CreateHisto1("NClustersVsLayer", | |
96 | "mean number of clusters per layer", | |
97 | nLayer, 0.5, nLayer+0.5, | |
98 | "layer", "<N_{clusters}>", | |
99 | AliMonitorHisto::kNormEvents); | |
100 | ||
101 | fNTracks = CreateHisto1("NTracks", "number of tracks per event", | |
102 | 300, 0, 30000, "N_{tracks}", "#Delta N/N", | |
103 | AliMonitorHisto::kNormEvents); | |
104 | ||
105 | fNTracksITSTPC = CreateHisto2("NTracksTPCITS", | |
106 | "correlation of number of ITS and TPC tracks per event", | |
107 | 30, 0, 30000, 30, 0, 30000, | |
108 | "N_{tracks}(ITS)", "N_{tracks}(TPC)", | |
109 | "#Delta N/N", AliMonitorHisto::kNormEvents); | |
110 | ||
111 | fTrackPt = CreateHisto1("TrackPt", "pt distribution of tracks", | |
112 | 90, 0, 3, "p_{t} [GeV/c]", "#Delta N/N", | |
113 | AliMonitorHisto::kNormEntries); | |
114 | ||
115 | fTrackEta = CreateHisto1("TrackEta", "eta distribution of tracks", | |
116 | 100, -2, 2, "#eta", "#Delta N/N", | |
117 | AliMonitorHisto::kNormEntries); | |
118 | ||
119 | fTrackPhi = CreateHisto1("TrackPhi", "phi distribution of tracks", | |
120 | 120, 0, 360, "#phi [#circ]", "#Delta N/N", | |
121 | AliMonitorHisto::kNormEntries); | |
122 | ||
123 | fTrackDEdxVsP = CreateHisto2("TrackDEdxVsP", "dE/dx of tracks", | |
124 | 100, 0, 3, 100, 0, 150, | |
125 | "p [GeV/c]", "dE/dx", "#Delta N/N", | |
126 | AliMonitorHisto::kNormEntries); | |
127 | } | |
128 | ||
129 | ||
130 | //_____________________________________________________________________________ | |
131 | void AliMonitorITS::FillHistos(AliRunLoader* runLoader, | |
132 | AliRawReader* rawReader) | |
133 | { | |
134 | // fill the ITS monitor histogrms | |
135 | ||
136 | rawReader->Reset(); | |
137 | AliITSRawStreamSDD inputSDD(rawReader); | |
138 | while (inputSDD.Next()) { | |
139 | fSDDDigitsCharge->Fill(inputSDD.GetSignal()); | |
140 | } | |
141 | ||
142 | rawReader->Reset(); | |
143 | AliITSRawStreamSSD inputSSD(rawReader); | |
144 | while (inputSSD.Next()) { | |
145 | fSSDDigitsCharge->Fill(inputSSD.GetSignal()); | |
146 | } | |
147 | ||
148 | AliLoader* itsLoader = runLoader->GetLoader("ITSLoader"); | |
149 | if (!itsLoader) return; | |
150 | ||
151 | ||
152 | itsLoader->LoadRecPoints(); | |
153 | TTree* clustersTree = itsLoader->TreeR(); | |
154 | if (!clustersTree) return; | |
155 | TClonesArray* clusters = new TClonesArray("AliITSclusterV2"); | |
156 | clustersTree->SetBranchAddress("Clusters", &clusters); | |
157 | ||
158 | for (Int_t iModule = 0; iModule < clustersTree->GetEntries(); iModule++) { | |
159 | clustersTree->GetEntry(iModule); | |
160 | Int_t iLayer, iLadder, iDetector; | |
161 | fGeom->GetModuleId(iModule, iLayer, iLadder, iDetector); | |
162 | ||
163 | for (Int_t j = 0; j < clusters->GetEntriesFast(); j++) { | |
164 | AliITSclusterV2* cluster = (AliITSclusterV2*) clusters->At(j); | |
165 | switch (iLayer) { | |
166 | case 1: case 2: { | |
167 | fSPDNClustersVsModule->Fill(iModule - fGeom->GetStartSPD()); | |
168 | break; | |
169 | } | |
170 | case 3: case 4: { | |
171 | fSDDClustersCharge->Fill(cluster->GetQ()); | |
172 | fSDDNClustersVsModule->Fill(iModule - fGeom->GetStartSDD()); | |
173 | break; | |
174 | } | |
175 | case 5: case 6: { | |
176 | fSSDClustersCharge->Fill(cluster->GetQ()); | |
177 | fSSDNClustersVsModule->Fill(iModule - fGeom->GetStartSSD()); | |
178 | break; | |
179 | } | |
180 | } | |
181 | fNClustersVsLayer->Fill(iLayer); | |
182 | } | |
183 | } | |
184 | ||
185 | delete clusters; | |
186 | itsLoader->UnloadRecPoints(); | |
187 | ||
188 | ||
189 | itsLoader->LoadTracks(); | |
190 | TTree* tracks = itsLoader->TreeT(); | |
191 | if (!tracks) return; | |
192 | AliITStrackV2* track = new AliITStrackV2; | |
193 | tracks->SetBranchAddress("tracks", &track); | |
194 | ||
195 | fNTracks->Fill(tracks->GetEntries()); | |
196 | for (Int_t i = 0; i < tracks->GetEntries(); i++) { | |
197 | tracks->GetEntry(i); | |
198 | fTrackPt->Fill(track->Pt()); | |
199 | fTrackEta->Fill(track->Eta()); | |
200 | fTrackPhi->Fill(track->Phi() * TMath::RadToDeg()); | |
201 | fTrackDEdxVsP->Fill(track->P(), track->GetdEdx()); | |
202 | } | |
203 | ||
204 | AliLoader* tpcLoader = runLoader->GetLoader("TPCLoader"); | |
205 | if (tpcLoader) { | |
206 | tpcLoader->LoadTracks(); | |
207 | TTree* tracksTPC = tpcLoader->TreeT(); | |
208 | if (tracksTPC) { | |
209 | fNTracksITSTPC->Fill(tracks->GetEntries(), tracksTPC->GetEntries()); | |
210 | } | |
211 | tpcLoader->UnloadTracks(); | |
212 | } | |
213 | ||
214 | delete track; | |
215 | itsLoader->UnloadTracks(); | |
216 | } |