]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MONITOR/AliMonitorHLTHough.cxx
- fixing compilation warnings
[u/mrichter/AliRoot.git] / MONITOR / AliMonitorHLTHough.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
17 ///////////////////////////////////////////////////////////////////////////////
18 //                                                                           //
19 // This class creates and fills monitor histograms for HLT Hough transform   //
20 //                                                                           //
21 ///////////////////////////////////////////////////////////////////////////////
22
23
24 #include "AliMonitorHLTHough.h"
25 #include "AliMonitorHisto.h"
26 #include "AliMonitorTrend.h"
27 #include "AliTPCParam.h"
28 #include "AliLog.h"
29 #include <TFolder.h>
30 #include <stdlib.h>
31 #include <AliHLTMemHandler.h>
32 #include <AliHLTTrackArray.h>
33 #include <AliHLTSpacePointData.h>
34 #include <AliHLTHoughTrack.h>
35 #include <AliHLTTransform.h>
36
37 //_____________________________________________________________________________
38 AliMonitorHLTHough::AliMonitorHLTHough(AliTPCParam* param):
39   AliMonitor(),
40   fParam(param),
41   fClustersCharge(NULL),
42   fNClustersVsRow(NULL),
43   fNClustersVsSector(NULL),
44   fNTracks(NULL),
45   fTrackPt(NULL),
46   fTrackEta(NULL),
47   fTrackPhi(NULL),
48   fTrackNHits(NULL),
49   fTrackDEdxVsP(NULL),
50   fTrackDEdx(NULL),
51   fTrackEtaVsPhi(NULL),
52   fPtEtaVsPhi(NULL)
53 {
54 // create a HLT monitor object with the given parameters
55
56 }
57
58 //_____________________________________________________________________________
59 void AliMonitorHLTHough::CreateHistos(TFolder* folder)
60 {
61 // create the HLT Hough transform monitor histograms
62
63   fFolder = folder->AddFolder("HLTHOUGH", "HLTHOUGH");
64
65   fClustersCharge = CreateHisto1("ClustersCharge", 
66                                  "charge distribution of clusters", 
67                                  100, 0, 1000, "charge", "#Delta N/N",
68                                  AliMonitorHisto::kNormEvents);
69
70   Int_t nRows = fParam->GetNRowLow() + fParam->GetNRowUp();
71   fNClustersVsRow = CreateHisto1("NClustersVsRow", 
72                                  "mean number of clusters per pad row", 
73                                  nRows, -0.5, nRows-0.5,
74                                  "pad row", "<N_{clusters}>",
75                                  AliMonitorHisto::kNormEvents);
76
77   Int_t nSector = fParam->GetNInnerSector();
78   fNClustersVsSector = CreateHisto1("NClustersVsSector", 
79                                     "mean number of clusters per sector", 
80                                     nSector, -0.5, nSector-0.5, 
81                                     "sector", "<N_{clusters}>",
82                                     AliMonitorHisto::kNormEvents);
83
84   fNTracks = CreateTrend("NTracks", "number of tracks per event", 
85                          "N_{tracks}");
86
87   fTrackPt = CreateHisto1("TrackPt", "pt distribution of tracks", 
88                           90, 0, 3, "p_{t} [GeV/c]", "#Delta N/N",
89                           AliMonitorHisto::kNormNone);
90
91   fTrackEta = CreateHisto1("TrackEta", "eta distribution of tracks", 
92                            100, -2, 2, "#eta", "#Delta N/N",
93                            AliMonitorHisto::kNormEntries);
94
95   fTrackPhi = CreateHisto1("TrackPhi", "phi distribution of tracks", 
96                            120, 0, 360, "#phi [#circ]", "#Delta N/N",
97                            AliMonitorHisto::kNormEntries);
98
99   fTrackNHits = CreateHisto1("TrackNHits", "Number of hits per track", 
100                            200, 0, 200, "N_{hits}", "#Delta N/N",
101                            AliMonitorHisto::kNormNone);
102
103   fTrackDEdxVsP = CreateHisto2("TrackDEdxVsP", "dE/dx of tracks", 
104                                100, 0, 3, 100, 0, 1000, 
105                                "p [GeV/c]", "dE/dx", "#Delta N/N",
106                                AliMonitorHisto::kNormEntries);
107
108   fTrackDEdx = CreateHisto1("TrackDEdx", "dE/dx for tracks with 0.4<p<1.0 GeV/c", 
109                                50, 0, 300, 
110                                "dE/dx", "#Delta N/N",
111                                AliMonitorHisto::kNormEntries);
112
113   fTrackEtaVsPhi = CreateHisto2("TrackEtaVsPhi", "#phi vs #eta", 
114                                20, -1, 1, 25, 0, 360, 
115                                "#eta", "#phi", "#Delta N/N",
116                                AliMonitorHisto::kNormNone);
117
118   fPtEtaVsPhi = CreateHisto2("PtEtaVsPhi", "#phi vs #eta", 
119                                20, -1, 1, 25, 0, 360, 
120                                "#eta", "#phi", "#Delta N/N",
121                                AliMonitorHisto::kNormNone);
122
123 }
124
125
126 //_____________________________________________________________________________
127 void AliMonitorHLTHough::FillHistos(AliRunLoader* /*runLoader*/, 
128                                     AliRawReader* /*rawReader*/, 
129                                     AliESDEvent* /*esd*/)
130 {
131 // fill the HLT Hough transform monitor histograms
132
133   AliHLTMemHandler clusterHandler[36][6];
134   AliHLTSpacePointData *clusters[36][6];
135   for (Int_t iSector = 0; iSector < fParam->GetNInnerSector(); iSector++) {
136     for (Int_t iPatch = 0; iPatch < 6; iPatch++) {
137       char fileName[256];
138       sprintf(fileName, "hlt/fitter/points_%d_%d.raw", iSector,iPatch);
139       if (!clusterHandler[iSector][iPatch].SetBinaryInput(fileName)) {
140         AliWarning(Form("could not open file %s", fileName));
141         continue;
142       }
143       clusters[iSector][iPatch] = (AliHLTSpacePointData*) clusterHandler[iSector][iPatch].Allocate();
144       UInt_t nClusters = 0;
145       clusterHandler[iSector][iPatch].Binary2Memory(nClusters, clusters[iSector][iPatch]);
146
147       for (UInt_t iCluster = 0; iCluster < nClusters; iCluster++) {
148         AliHLTSpacePointData& cluster = clusters[iSector][iPatch][iCluster];
149         fClustersCharge->Fill(cluster.fCharge);
150         fNClustersVsRow->Fill(cluster.fPadRow);
151         fNClustersVsSector->Fill(iSector);
152       }
153
154       clusterHandler[iSector][iPatch].CloseBinaryInput();
155     }
156   }
157
158   fNClustersVsSector->ScaleErrorBy(10.);
159
160
161   AliHLTMemHandler memHandler;
162   Int_t nHoughTracks = 0;
163
164   char fileName[256];
165   sprintf(fileName, "hlt/fitter/tracks.raw");
166   if (!memHandler.SetBinaryInput(fileName)) {
167     AliWarning("could not open file hlt/fitter/tracks.raw");
168     return;
169   }
170   AliHLTTrackArray* tracks = new AliHLTTrackArray;
171   memHandler.Binary2TrackArray(tracks);
172
173   nHoughTracks += tracks->GetNTracks();
174   for (Int_t iTrack = 0; iTrack < tracks->GetNTracks(); iTrack++) {
175     AliHLTHoughTrack *track = (AliHLTHoughTrack*)tracks->GetCheckedTrack(iTrack);
176     if(!track) continue;
177
178     track->CalculateHelix();
179     //  track->Rotate(iSector);
180
181     fTrackPt->Fill(track->GetPt());
182     fTrackEta->Fill(track->GetPseudoRapidity());
183     fTrackPhi->Fill(track->GetPsi() * TMath::RadToDeg());
184     if(track->GetPt()>3.) {
185       fTrackEtaVsPhi->Fill(track->GetPseudoRapidity(),track->GetPsi() * TMath::RadToDeg());
186       fPtEtaVsPhi->Fill(track->GetPseudoRapidity(),track->GetPsi() * TMath::RadToDeg(),track->GetPt());
187     }
188     fTrackNHits->Fill(track->GetNHits());
189
190     // Track dEdx
191     Int_t nc=track->GetNHits();
192     UInt_t *hits = track->GetHitNumbers();
193     Float_t sampleDEdx[159];
194     for (Int_t iHit = 0; iHit < nc; iHit++) {
195       UInt_t hitID = hits[iHit];
196       Int_t iSector = (hitID>>25) & 0x7f;
197       Int_t iPatch = (hitID>>22) & 0x7;
198       UInt_t position = hitID&0x3fffff;
199       UChar_t padrow = clusters[iSector][iPatch][position].fPadRow;
200       Float_t pWidth = AliHLTTransform::GetPadPitchWidthLow();
201       if (padrow>63)
202         pWidth = AliHLTTransform::GetPadPitchWidthUp(); 
203       Float_t corr=1.; if (padrow>63) corr=0.67;
204       sampleDEdx[iHit] = clusters[iSector][iPatch][position].fCharge/pWidth*corr;
205       Double_t crossingangle = track->GetCrossingAngle(padrow,iSector);
206       Double_t s = sin(crossingangle);
207       Double_t t = track->GetTgl();
208       sampleDEdx[iHit] *= TMath::Sqrt((1-s*s)/(1+t*t));
209     }
210
211     /* Cook dEdx */
212     Int_t i;
213     Int_t swap;//stupid sorting
214     do {
215       swap=0;
216       for (i=0; i<nc-1; i++) {
217         if (sampleDEdx[i]<=sampleDEdx[i+1]) continue;
218         Float_t tmp=sampleDEdx[i];
219         sampleDEdx[i]=sampleDEdx[i+1]; sampleDEdx[i+1]=tmp;
220         swap++;
221       }
222     } while (swap);
223
224     Double_t low=0.05; Double_t up=0.7;
225     Int_t nl=Int_t(low*nc), nu=Int_t(up*nc);
226     Float_t trackDEdx=0;
227     for (i=nl; i<=nu; i++) trackDEdx += sampleDEdx[i];
228     trackDEdx /= (nu-nl+1);
229
230     fTrackDEdxVsP->Fill(track->GetP(),trackDEdx);
231     if(track->GetP()>0.4 && track->GetP()<1.0)
232       fTrackDEdx->Fill(trackDEdx);
233   }
234   fNTracks->Fill(nHoughTracks);
235
236   delete tracks;
237   memHandler.CloseBinaryInput();
238 }