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