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