]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MONITOR/AliMonitorTPC.cxx
Removing warnings. From now on they are considered as errors
[u/mrichter/AliRoot.git] / MONITOR / AliMonitorTPC.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 TPC          //
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 "AliMonitorDataTPC.h"
31 #include "AliMonitorHisto.h"
32 #include "AliMonitorTPC.h"
33 #include "AliMonitorTrend.h"
34 #include "AliRawReader.h"
35 #include "AliRunLoader.h"
36 #include "AliTPCClustersRow.h"
37 #include "AliTPCParam.h"
38 #include "AliTPCRawStream.h"
39 #include "AliTPCclusterMI.h"
40
41 ClassImp(AliMonitorTPC) 
42
43
44 //_____________________________________________________________________________
45 AliMonitorTPC::AliMonitorTPC(AliTPCParam* param)
46 {
47 // create a TPC monitor object with the given parameters
48
49   fParam = param;
50   fData = new AliMonitorDataTPC(10000);
51 }
52
53 //_____________________________________________________________________________
54 AliMonitorTPC::AliMonitorTPC(const AliMonitorTPC& monitor) :
55   AliMonitor(monitor)
56 {
57   AliFatal("copy constructor not implemented");
58 }
59
60 //_____________________________________________________________________________
61 AliMonitorTPC& AliMonitorTPC::operator = (const AliMonitorTPC& /*monitor*/)
62 {
63   AliFatal("assignment operator not implemented");
64   return *this;
65 }
66
67 //_____________________________________________________________________________
68 AliMonitorTPC::~AliMonitorTPC()
69 {
70   delete fData;
71 }
72
73
74 //_____________________________________________________________________________
75 void AliMonitorTPC::CreateHistos(TFolder* folder)
76 {
77 // create the TPC monitor histograms
78
79   fFolder = folder->AddFolder("TPC", "TPC");
80
81   fPadsCharge = CreateHisto1("PadsCharge", "charge distribution of pads", 
82                              100, 0, 200, "charge", "#Delta N/N",
83                              AliMonitorHisto::kNormEvents);
84
85   fClustersCharge = CreateHisto1("ClustersCharge", 
86                                  "charge distribution of clusters", 
87                                  100, 0, 1000, "charge", "#Delta N/N",
88                                  AliMonitorHisto::kNormEvents);
89
90   Int_t nRows = fParam->GetNRowLow() + fParam->GetNRowUp();
91   fNClustersVsRow = CreateHisto1("NClustersVsRow", 
92                                  "mean number of clusters per pad row", 
93                                  nRows, -0.5, nRows-0.5,
94                                  "pad row", "<N_{clusters}>",
95                                  AliMonitorHisto::kNormEvents);
96
97   Int_t nSector = fParam->GetNInnerSector();
98   fNClustersVsSector = CreateHisto1("NClustersVsSector", 
99                                     "mean number of clusters per sector", 
100                                     nSector, -0.5, nSector-0.5, 
101                                     "sector", "<N_{clusters}>",
102                                     AliMonitorHisto::kNormEvents);
103
104   fNTracks = CreateTrend("NTracks", "number of tracks per event", 
105                          "N_{tracks}");
106
107   fTrackPt = CreateHisto1("TrackPt", "pt distribution of tracks", 
108                           90, 0, 3, "p_{t} [GeV/c]", "#Delta N/N",
109                           AliMonitorHisto::kNormNone);
110
111   fTrackEta = CreateHisto1("TrackEta", "eta distribution of tracks", 
112                            100, -2, 2, "#eta", "#Delta N/N",
113                            AliMonitorHisto::kNormEntries);
114
115   fTrackPhi = CreateHisto1("TrackPhi", "phi distribution of tracks", 
116                            120, -180, 180, "#phi [#circ]", "#Delta N/N",
117                            AliMonitorHisto::kNormEntries);
118   fTrackPhi->SetDescription("The phi distribution should be flat on average.\nIf it is not flat check for dead TPC sectors.");
119
120   fTrackNCl = CreateHisto1("TrackNCl", "Number of clusters per track", 
121                            200, 0, 200, "N_{clusters}", "#Delta N/N",
122                            AliMonitorHisto::kNormNone);
123
124   fTrackDEdxVsP = CreateHisto2("TrackDEdxVsP", "dE/dx of tracks", 
125                                100, 0, 3, 100, 0, 200, 
126                                "p [GeV/c]", "dE/dx", "#Delta N/N",
127                                AliMonitorHisto::kNormEntries);
128
129   fTrackDEdx = CreateHisto1("TrackDEdx", "dE/dx of tracks with 0.4<p<1.0 GeV/c", 
130                                50, 0, 100, 
131                                "dE/dx", "#Delta N/N",
132                                AliMonitorHisto::kNormEntries);
133
134   fTrackEtaVsPhi = CreateHisto2("TrackEtaVsPhi", "#phi vs #eta", 
135                                20, -1, 1, 25, 0, 360, 
136                                "#eta", "#phi", "#Delta N/N",
137                                AliMonitorHisto::kNormNone);
138
139   fPtEtaVsPhi = CreateHisto2("PtEtaVsPhi", "#phi vs #eta", 
140                                20, -1, 1, 25, 0, 360, 
141                                "#eta", "#phi", "#Delta N/N",
142                                AliMonitorHisto::kNormNone);
143
144 }
145
146
147 //_____________________________________________________________________________
148 void AliMonitorTPC::CreateBranches(TTree* tree)
149 {
150 // create a branch with TPC variables
151
152   tree->Branch("TPC", "AliMonitorDataTPC", &fData);
153 }
154
155
156 //_____________________________________________________________________________
157 void AliMonitorTPC::FillHistos(AliRunLoader* runLoader, 
158                                AliRawReader* rawReader, AliESD* esd)
159 {
160 // fill the TPC monitor histogrms
161
162   rawReader->Reset();
163   AliTPCRawStream input(rawReader);
164   while (input.Next()) {
165     fPadsCharge->Fill(input.GetSignal());
166   }
167
168   AliLoader* tpcLoader = runLoader->GetLoader("TPCLoader");
169   if (!tpcLoader) return;
170
171
172   tpcLoader->LoadRecPoints();
173   TTree* clusters = tpcLoader->TreeR();
174   if (!clusters) return;
175   AliTPCClustersRow* clustersRow = new AliTPCClustersRow;
176   clustersRow->SetClass("AliTPCclusterMI");
177   clusters->SetBranchAddress("Segment", &clustersRow);
178
179   for (Int_t i = 0; i < clusters->GetEntries(); i++) {
180     clusters->GetEntry(i);
181     Int_t iSector, iRow;
182     fParam->AdjustSectorRow(clustersRow->GetID(), iSector, iRow);
183     if (iSector >= fParam->GetNInnerSector()) {
184       iRow += fParam->GetNRowLow();
185       iSector -= fParam->GetNInnerSector();
186     }
187
188     TObjArray* array = clustersRow->GetArray();
189     for (Int_t j = 0; j < array->GetEntriesFast(); j++) {
190       AliTPCclusterMI* cluster = (AliTPCclusterMI*) array->At(j);
191       fClustersCharge->Fill(cluster->GetQ());
192       fNClustersVsRow->Fill(iRow);
193       fNClustersVsSector->Fill(iSector);
194     }
195   }
196   fNClustersVsSector->ScaleErrorBy(10.);
197
198   delete clustersRow;
199   tpcLoader->UnloadRecPoints();
200
201
202   Int_t nTracks = 0;
203   for (Int_t i = 0; i < esd->GetNumberOfTracks(); i++) {
204     AliESDtrack* track = esd->GetTrack(i);
205     if (!track || ((track->GetStatus() | AliESDtrack::kTPCin) == 0)) continue;
206     nTracks++;
207
208     Double_t pxyz[3];
209     track->GetInnerPxPyPz(pxyz);
210     TVector3 pTrack(pxyz);
211     Double_t p = pTrack.Mag();
212     Double_t pt = pTrack.Pt();
213     Double_t eta = pTrack.Eta();
214     Double_t phi = pTrack.Phi() * TMath::RadToDeg();
215
216     fTrackPt->Fill(pt);
217     fTrackEta->Fill(eta);
218     fTrackPhi->Fill(phi);
219     if (pt > 3.) {
220       fTrackEtaVsPhi->Fill(eta, phi);
221       fPtEtaVsPhi->Fill(eta, phi, pTrack.Pt());
222     }
223     fTrackNCl->Fill(track->GetTPCclusters(NULL));
224     fTrackDEdxVsP->Fill(p, track->GetTPCsignal());
225     if(p>0.4 && p<1.0)
226       fTrackDEdx->Fill(track->GetTPCsignal());
227
228     fData->SetData(i, pt, eta, phi); 
229   }
230   fNTracks->Fill(nTracks);
231   fData->SetNTracks(nTracks);
232 }