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