1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 ///////////////////////////////////////////////////////////////////////////////
20 // This class creates and fills the monitor histograms for the ITS //
22 ///////////////////////////////////////////////////////////////////////////////
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"
37 ClassImp(AliMonitorITS)
40 //_____________________________________________________________________________
41 AliMonitorITS::AliMonitorITS(AliITSgeom* geom)
43 // create a ITS monitor object with the given geometry
49 //_____________________________________________________________________________
50 AliMonitorITS::AliMonitorITS(const AliMonitorITS& monitor) :
53 Fatal("AliMonitorITS", "copy constructor not implemented");
56 //_____________________________________________________________________________
57 AliMonitorITS& AliMonitorITS::operator = (const AliMonitorITS& /*monitor*/)
59 Fatal("operator =", "assignment operator not implemented");
63 //_____________________________________________________________________________
64 void AliMonitorITS::CreateHistos(TFolder* folder)
66 // create the ITS monitor histograms
68 fFolder = folder->AddFolder("ITS", "ITS");
70 fSDDDigitsCharge = CreateHisto1("SDDDigitsCharge",
71 "charge distribution of SDD digits",
72 128, 0, 256, "charge", "#Delta N/N",
73 AliMonitorHisto::kNormEvents);
75 fSSDDigitsCharge = CreateHisto1("SSDDigitsCharge",
76 "charge distribution of SSD digits",
77 100, 0, 200, "charge", "#Delta N/N",
78 AliMonitorHisto::kNormEvents);
80 fSDDClustersCharge = CreateHisto1("SDDClustersCharge",
81 "charge distribution of SDD clusters",
82 100, 0, 200, "charge", "#Delta N/N",
83 AliMonitorHisto::kNormEvents);
85 fSSDClustersCharge = CreateHisto1("SSDClustersCharge",
86 "charge distribution of SSD clusters",
87 100, 0, 400, "charge", "#Delta N/N",
88 AliMonitorHisto::kNormEvents);
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);
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);
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);
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);
118 fNTracks = CreateHisto1("NTracks", "number of tracks per event",
119 300, 0, 30000, "N_{tracks}", "#Delta N/N",
120 AliMonitorHisto::kNormEvents);
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);
128 fTrackPt = CreateHisto1("TrackPt", "pt distribution of tracks",
129 90, 0, 3, "p_{t} [GeV/c]", "#Delta N/N",
130 AliMonitorHisto::kNormEntries);
132 fTrackEta = CreateHisto1("TrackEta", "eta distribution of tracks",
133 100, -2, 2, "#eta", "#Delta N/N",
134 AliMonitorHisto::kNormEntries);
136 fTrackPhi = CreateHisto1("TrackPhi", "phi distribution of tracks",
137 120, 0, 360, "#phi [#circ]", "#Delta N/N",
138 AliMonitorHisto::kNormEntries);
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);
147 //_____________________________________________________________________________
148 void AliMonitorITS::FillHistos(AliRunLoader* runLoader,
149 AliRawReader* rawReader)
151 // fill the ITS monitor histogrms
154 AliITSRawStreamSDD inputSDD(rawReader);
155 while (inputSDD.Next()) {
156 fSDDDigitsCharge->Fill(inputSDD.GetSignal());
160 AliITSRawStreamSSD inputSSD(rawReader);
161 while (inputSSD.Next()) {
162 fSSDDigitsCharge->Fill(inputSSD.GetSignal());
165 AliLoader* itsLoader = runLoader->GetLoader("ITSLoader");
166 if (!itsLoader) return;
169 itsLoader->LoadRecPoints();
170 TTree* clustersTree = itsLoader->TreeR();
171 if (!clustersTree) return;
172 TClonesArray* clusters = new TClonesArray("AliITSclusterV2");
173 clustersTree->SetBranchAddress("Clusters", &clusters);
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);
180 for (Int_t j = 0; j < clusters->GetEntriesFast(); j++) {
181 AliITSclusterV2* cluster = (AliITSclusterV2*) clusters->At(j);
184 fSPDNClustersVsModule->Fill(iModule - fGeom->GetStartSPD());
188 fSDDClustersCharge->Fill(cluster->GetQ());
189 fSDDNClustersVsModule->Fill(iModule - fGeom->GetStartSDD());
193 fSSDClustersCharge->Fill(cluster->GetQ());
194 fSSDNClustersVsModule->Fill(iModule - fGeom->GetStartSSD());
198 fNClustersVsLayer->Fill(iLayer);
203 itsLoader->UnloadRecPoints();
206 itsLoader->LoadTracks();
207 TTree* tracks = itsLoader->TreeT();
209 AliITStrackV2* track = new AliITStrackV2;
210 tracks->SetBranchAddress("tracks", &track);
212 fNTracks->Fill(tracks->GetEntries());
213 for (Int_t i = 0; i < tracks->GetEntries(); 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());
221 AliLoader* tpcLoader = runLoader->GetLoader("TPCLoader");
223 tpcLoader->LoadTracks();
224 TTree* tracksTPC = tpcLoader->TreeT();
226 fNTracksITSTPC->Fill(tracks->GetEntries(), tracksTPC->GetEntries());
228 tpcLoader->UnloadTracks();
232 itsLoader->UnloadTracks();