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 "AliITSRawStreamSPD.h"
28 #include "AliITSRawStreamSDD.h"
29 #include "AliITSRawStreamSSD.h"
30 #include "AliITSclusterV2.h"
31 #include "AliITStrackV2.h"
34 ClassImp(AliMonitorITS)
37 //_____________________________________________________________________________
38 AliMonitorITS::AliMonitorITS(AliITSgeom* geom)
40 // create a ITS monitor object with the given geometry
46 //_____________________________________________________________________________
47 void AliMonitorITS::CreateHistos(TFolder* folder)
49 // create the ITS monitor histograms
51 fFolder = folder->AddFolder("ITS", "ITS");
53 fSDDDigitsCharge = CreateHisto1("SDDDigitsCharge",
54 "charge distribution of SDD digits",
55 128, 0, 256, "charge", "#Delta N/N",
56 AliMonitorHisto::kNormEvents);
58 fSSDDigitsCharge = CreateHisto1("SSDDigitsCharge",
59 "charge distribution of SSD digits",
60 100, 0, 200, "charge", "#Delta N/N",
61 AliMonitorHisto::kNormEvents);
63 fSDDClustersCharge = CreateHisto1("SDDClustersCharge",
64 "charge distribution of SDD clusters",
65 100, 0, 200, "charge", "#Delta N/N",
66 AliMonitorHisto::kNormEvents);
68 fSSDClustersCharge = CreateHisto1("SSDClustersCharge",
69 "charge distribution of SSD clusters",
70 100, 0, 400, "charge", "#Delta N/N",
71 AliMonitorHisto::kNormEvents);
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);
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);
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);
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);
101 fNTracks = CreateHisto1("NTracks", "number of tracks per event",
102 300, 0, 30000, "N_{tracks}", "#Delta N/N",
103 AliMonitorHisto::kNormEvents);
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);
111 fTrackPt = CreateHisto1("TrackPt", "pt distribution of tracks",
112 90, 0, 3, "p_{t} [GeV/c]", "#Delta N/N",
113 AliMonitorHisto::kNormEntries);
115 fTrackEta = CreateHisto1("TrackEta", "eta distribution of tracks",
116 100, -2, 2, "#eta", "#Delta N/N",
117 AliMonitorHisto::kNormEntries);
119 fTrackPhi = CreateHisto1("TrackPhi", "phi distribution of tracks",
120 120, 0, 360, "#phi [#circ]", "#Delta N/N",
121 AliMonitorHisto::kNormEntries);
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);
130 //_____________________________________________________________________________
131 void AliMonitorITS::FillHistos(AliRunLoader* runLoader,
132 AliRawReader* rawReader)
134 // fill the ITS monitor histogrms
137 AliITSRawStreamSDD inputSDD(rawReader);
138 while (inputSDD.Next()) {
139 fSDDDigitsCharge->Fill(inputSDD.GetSignal());
143 AliITSRawStreamSSD inputSSD(rawReader);
144 while (inputSSD.Next()) {
145 fSSDDigitsCharge->Fill(inputSSD.GetSignal());
148 AliLoader* itsLoader = runLoader->GetLoader("ITSLoader");
149 if (!itsLoader) return;
152 itsLoader->LoadRecPoints();
153 TTree* clustersTree = itsLoader->TreeR();
154 if (!clustersTree) return;
155 TClonesArray* clusters = new TClonesArray("AliITSclusterV2");
156 clustersTree->SetBranchAddress("Clusters", &clusters);
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);
163 for (Int_t j = 0; j < clusters->GetEntriesFast(); j++) {
164 AliITSclusterV2* cluster = (AliITSclusterV2*) clusters->At(j);
167 fSPDNClustersVsModule->Fill(iModule - fGeom->GetStartSPD());
171 fSDDClustersCharge->Fill(cluster->GetQ());
172 fSDDNClustersVsModule->Fill(iModule - fGeom->GetStartSDD());
176 fSSDClustersCharge->Fill(cluster->GetQ());
177 fSSDNClustersVsModule->Fill(iModule - fGeom->GetStartSSD());
181 fNClustersVsLayer->Fill(iLayer);
186 itsLoader->UnloadRecPoints();
189 itsLoader->LoadTracks();
190 TTree* tracks = itsLoader->TreeT();
192 AliITStrackV2* track = new AliITStrackV2;
193 tracks->SetBranchAddress("tracks", &track);
195 fNTracks->Fill(tracks->GetEntries());
196 for (Int_t i = 0; i < tracks->GetEntries(); 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());
204 AliLoader* tpcLoader = runLoader->GetLoader("TPCLoader");
206 tpcLoader->LoadTracks();
207 TTree* tracksTPC = tpcLoader->TreeT();
209 fNTracksITSTPC->Fill(tracks->GetEntries(), tracksTPC->GetEntries());
211 tpcLoader->UnloadTracks();
215 itsLoader->UnloadTracks();