FairRoot/PandaRoot
PndCAGlobalPerformance.cxx
Go to the documentation of this file.
1 // $Id: PndCAGlobalPerformance.cxx,v 1.11 2010/08/18 20:46:09 ikulakov Exp $
2 // **************************************************************************
3 // This file is property of and copyright by the ALICE HLT Project *
4 // ALICE Experiment at CERN, All rights reserved. *
5 // *
6 // Primary Authors: Sergey Gorbunov <sergey.gorbunov@kip.uni-heidelberg.de> *
7 // Ivan Kisel <kisel@kip.uni-heidelberg.de> *
8 // for The ALICE HLT Project. *
9 // *
10 // Developed by: Igor Kulakov <I.Kulakov@gsi.de> *
11 // Maksym Zyzak <M.Zyzak@gsi.de> *
12 // *
13 // Permission to use, copy, modify and distribute this software and its *
14 // documentation strictly for non-commercial purposes is hereby granted *
15 // without fee, provided that the above copyright notice appears in all *
16 // copies and that both the copyright notice and this permission notice *
17 // appear in the supporting documentation. The authors make no claims *
18 // about the suitability of this software for any purpose. It is *
19 // provided "as is" without express or implied warranty. *
20 // *
21 //***************************************************************************
22 // List of most recent changes: *
23 // *
24 // 29-05-17 Adjustment of initialization for PANDA data (Irina Rostovtseva) *
25 //***************************************************************************
26 
27 #ifdef DO_TPCCATRACKER_EFF_PERFORMANCE
28 
29 #include "PndCACounters.h"
30 
32 #include "PndCAGlobalPerformance.h"
33 
34 
35 #include "PndCAGBHit.h"
36 #include "PndCAMCTrack.h"
37 #ifndef HLTCA_STANDALONE
38 #include "PndCAMCPoint.h"
39 #endif
40 #include "PndCAGBTrack.h"
41 #include "PndCAGBTracker.h"
42 
43 //#include "PndCADisplay.h"
44 
45 #include "TMath.h"
46 #include "TROOT.h"
47 #include "Riostream.h"
48 #include "TFile.h"
49 #include "TH1.h"
50 #include "TH2.h"
51 #include "TProfile.h"
52 #include "TStyle.h"
53 
54 void PndCAGlobalPerformance::SetNewEvent(const PndCAGBTracker * const tracker,
55  vector<PndCAHitLabel> *hitLabels,
56  vector<PndCAMCTrack> *mcTracks,
57  vector<PndCALocalMCPoint> *localMCPoints)
58 {
59  PndCATrackPerformanceBase::SetNewEvent(tracker, hitLabels, mcTracks, localMCPoints);
60 
61 
62  nRecoTracks = fTracker->NTracks();
63 
64 
65 } // void PndCAGlobalPerformance::SetNewEvent
66 
68 {
69  for ( int imc = 0; imc < nMCTracks; imc++ ) (*fMCTracks)[imc].SetNHits( 0 );
70  cout<<"nHits "<<fTracker->NHits()<<" n labels "<<(*fHitLabels).size()<<endl;
71 
72  for ( int ih = 0; ih < fTracker->NHits(); ih++ ) { // TODO: do we need to calculate consequtive hits??
73  const PndCAGBHit &hit = fTracker->Hit( ih );
74  const PndCAHitLabel &l = (*fHitLabels)[hit.ID()];
75  if( hit.IsLeft() ){
76  cout<<" mirrored hit in array, something wrong!! "<<endl;
77  continue;
78  }
79  //cout<<ih<<" "<<hit.ID()<<" "<<l.fLab[0]<<" "<<l.fLab[1]<<" "<<l.fLab[2]<<endl;
80  if ( l.fLab[0] >= 0 ) (*fMCTracks)[l.fLab[0]].SetNHits( (*fMCTracks)[l.fLab[0]].NHits() + 1 );
81  if ( l.fLab[1] >= 0 ) (*fMCTracks)[l.fLab[1]].SetNHits( (*fMCTracks)[l.fLab[1]].NHits() + 1 );
82  if ( l.fLab[2] >= 0 ) (*fMCTracks)[l.fLab[2]].SetNHits( (*fMCTracks)[l.fLab[2]].NHits() + 1 );
83  }
84  mcData.resize(nMCTracks);
85  for ( int imc = 0; imc < nMCTracks; imc++ ) {
86  PndCAMCTrack &mc = (*fMCTracks)[imc];
87  PndCAPerformanceMCTrackData &mcTrackData = mcData[imc];
88  mc.SetSet( 0 );
89  mc.SetNReconstructed( 0 );
90  mc.SetNTurns( 1 );
92  //if ( mc.NHits() == mc.NHitRows() )
93 // if ( mc.NMCPoints() >= PParameters::MinimumMCPointsForMCTrack )
95  mcTrackData.SetAsReconstructable();
96 
97  // sets of tracks 0-OutSet, 1-ExtraSet, 2-RefSet, 3-ExtraSecSet, 4-ExtraPrimSet, 5-RefSecSet, 6-RefPrimSet, 7-LongRefPrimSet
98  if ( mc.P() >= PParameters::ExtraThreshold ) {
99  if ( mc.P() >= PParameters::RefThreshold ) {
100  mc.SetSet( 2 );
101  mcTrackData.SetSet( 2 );
102  if ( mc.MotherId() != -1 ) {
103  mc.SetSet( 5 );
104  mcTrackData.SetSet( 5 );
105 
106  }
107  else {
108  mc.SetSet( 6 );
109  mcTrackData.SetSet( 6 );
110  if ( mc.NHits() == fTracker->NStations() ) {
111  mc.SetSet( 7 );
112  mcTrackData.SetSet( 7 );
113  }
114  }
115  }
116  else{
117  mc.SetSet( 1 );
118  mcTrackData.SetSet( 1 );
119  if ( mc.MotherId() != -1 ) {
120  mc.SetSet( 3 );
121  mcTrackData.SetSet( 3 );
122 
123  }
124  else {
125  mc.SetSet( 4 );
126  mcTrackData.SetSet( 4 );
127  }
128  }
129  }
130  } // for iMC
131 
132 } // void PndCAGlobalPerformance::CheckMCTracks()
133 
134 
136 {
137  cout<<"n reco trcks: "<<nRecoTracks<<" "<<fTracker->NTracks()<<endl;
138  recoData.resize(nRecoTracks);
139  for ( int itr = 0; itr < nRecoTracks; itr++ ) {
140  int traLabels = -1;
141  double traPurity = 0;
142  const PndCAGBTrack &tCA = fTracker->Track( itr );
143  const int nhits = tCA.NHits();
144  int *lb = new int[nhits*3];
145  int nla = 0;
146  for ( int ihit = 0; ihit < nhits; ihit++ ) {
147  const int index = fTracker->TrackHit( tCA.FirstHitRef() + ihit );
148  const PndCAHitLabel &l = (*fHitLabels)[fTracker->Hit( index ).ID()];
149  if ( l.fLab[0] >= 0 ) lb[nla++] = l.fLab[0];
150  if ( l.fLab[1] >= 0 ) lb[nla++] = l.fLab[1];
151  if ( l.fLab[2] >= 0 ) lb[nla++] = l.fLab[2];
152  }
153  std::sort( lb, lb + nla );
154  int labmax = -1, labcur = -1, lmax = 0, lcurr = 0;
155  for ( int i = 0; i < nla; i++ ) {
156  if ( lb[i] != labcur ) {
157  if ( labcur >= 0 && lmax < lcurr ) {
158  lmax = lcurr;
159  labmax = labcur;
160  }
161  labcur = lb[i];
162  lcurr = 0;
163  }
164  lcurr++;
165  }
166  if ( labcur >= 0 && lmax < lcurr ) {
167  lmax = lcurr;
168  labmax = labcur;
169  }
170  lmax = 0;
171  for ( int ihit = 0; ihit < nhits; ihit++ ) {
172  const int index = fTracker->TrackHit( tCA.FirstHitRef() + ihit );
173  const PndCAHitLabel &l = (*fHitLabels)[fTracker->Hit( index ).ID()];
174  if ( l.fLab[0] == labmax || l.fLab[1] == labmax || l.fLab[2] == labmax
175  ) lmax++;
176  }
177  traLabels = labmax;
178  traPurity = ( ( nhits > 0 ) ? double( lmax ) / double( nhits ) : 0 );
179  if ( lb ) delete[] lb;
180 
181  recoData[itr].SetMCTrack(traLabels, traPurity, nhits);
182  if ( recoData[itr].IsReco(PParameters::MinTrackPurity, PParameters::MinimumHitsForRecoTrack) ) mcData[traLabels].AddReconstructed(itr);
183  } // for iReco
184 } // void PndCAGlobalPerformance::MatchTracks()
185 
186 
188 {
189  for ( int iRTr = 0; iRTr < nRecoTracks; iRTr++ ) {
190  if ( recoData[iRTr].IsGhost(PParameters::MinTrackPurity) )
191  fEff.ghosts++;
192  }
193 
194  for ( int iMCTr = 0; iMCTr < nMCTracks; iMCTr++ ) {
195  PndCAPerformanceMCTrackData &mc = mcData[iMCTr];
196  if ( !mc.IsReconstructable() ) continue;
197  const bool reco = mc.IsReconstructed();
198  const int nclones = mc.GetNClones();
199 
200  bool killed = 0; // isn't used
201  double ratio_length = 0;
202  double ratio_fakes = 0;
203  const vector<int>& rTIds = mc.RecoTrackIds();
204  if (reco){
205  double mc_length = (*fMCTracks)[iMCTr].NMCRows();
206  for (unsigned int irt = 0; irt < rTIds.size(); irt++) {
207  PndCAPerformanceRecoTrackData& rd = recoData[rTIds[irt]];
208  ratio_length += static_cast<double>( rd.NHits() )*rd.GetPurity() / mc_length;
209  ratio_fakes += 1 - rd.GetPurity();
210  }
211  }
212 
213  if ( mc.GetSet() == 0){ // rest, out track
214  fEff.Inc(reco, killed, ratio_length, ratio_fakes, nclones, "rest");
215  }
216  else{ // good
217  fEff.Inc(reco, killed, ratio_length, ratio_fakes, nclones, "total");
218  switch( mc.GetSet() ){
219  case 1:
220  fEff.Inc(reco, killed, ratio_length, ratio_fakes, nclones, "extra");
221  break;
222  case 2:
223  fEff.Inc(reco, killed, ratio_length, ratio_fakes, nclones, "ref");
224  break;
225  case 3:
226  fEff.Inc(reco, killed, ratio_length, ratio_fakes, nclones, "extra");
227  fEff.Inc(reco, killed, ratio_length, ratio_fakes, nclones, "extra_sec");
228  break;
229  case 4:
230  fEff.Inc(reco, killed, ratio_length, ratio_fakes, nclones, "extra");
231  fEff.Inc(reco, killed, ratio_length, ratio_fakes, nclones, "extra_prim");
232  break;
233  case 5:
234  fEff.Inc(reco, killed, ratio_length, ratio_fakes, nclones, "ref");
235  fEff.Inc(reco, killed, ratio_length, ratio_fakes, nclones, "ref_sec");
236  break;
237  case 6:
238  fEff.Inc(reco, killed, ratio_length, ratio_fakes, nclones, "ref");
239  fEff.Inc(reco, killed, ratio_length, ratio_fakes, nclones, "ref_prim");
240  break;
241  break;
242  case 7:
243  fEff.Inc(reco, killed, ratio_length, ratio_fakes, nclones, "ref");
244  fEff.Inc(reco, killed, ratio_length, ratio_fakes, nclones, "ref_prim");
245  fEff.Inc(reco, killed, ratio_length, ratio_fakes, nclones, "ref_prim_long");
246  break;
247  }
248  }
249  } // for iMCTr
250  PndCATrackPerformanceBase::EfficiencyPerformance();
251 } // void PndCAGlobalPerformance::EfficiencyPerformance( )
252 
254 {
255  PndCATrackPerformanceBase::FillHistos();
256 
257  //
258 
259  const int NMCTracks = (*fMCTracks).size();
260  vector<int> mcTrackNRecoHits;
261  vector<int> nHitsVsRow;
262  vector<int> nMCPointsVsRow;
263  const int Multiplicity = (*fMCTracks).size();
264 
265  mcTrackNRecoHits.resize(NMCTracks, 0);
266  nHitsVsRow.resize(fTracker->NStations());
267  nMCPointsVsRow.resize(fTracker->NStations());
268  for(int iH=0; iH < fTracker->NHits(); iH++){
269  const PndCAGBHit &hit = fTracker->Hit( iH );
270 
271  nHitsVsRow[hit.IRow()]++;
272 
273  const PndCAHitLabel &l = (*fHitLabels)[hit.ID()];
274  if ( l.fLab[0] >= 0 ) mcTrackNRecoHits[l.fLab[0]]++;
275  if ( l.fLab[1] >= 0 ) mcTrackNRecoHits[l.fLab[1]]++;
276  if ( l.fLab[2] >= 0 ) mcTrackNRecoHits[l.fLab[2]]++;
277  }
278 
279  for(int i=0; i < NMCTracks; i++){
280  PndCAMCTrack &mcT = (*fMCTracks)[i];
281 
282  GetHisto("mcTrackNRecoHits")->Fill( mcTrackNRecoHits[i] );
283  GetHisto("nMCPointsVsMCMom")->Fill( mcT.P(), mcT.NMCPoints() );
284 
285  if ( mcT.NMCPoints() > 0 ) {
286  double mcEx = mcT.Px();
287  double mcEy = mcT.Py();
288  double mcEz = mcT.Pz();
289  double mcEt = TMath::Sqrt( mcEx * mcEx + mcEy * mcEy );
290 
291  const double dZdS = mcEz/mcEt;
292 
293  GetHisto("nHitsOverNMCPointsVsMCMom")->Fill( mcT.P(), float(mcTrackNRecoHits[i])/float(mcT.NMCPoints()) );
294  GetHisto("nHitsOverNMCPointsVsMCDzDs")->Fill( dZdS, float(mcTrackNRecoHits[i])/float(mcT.NMCPoints()) );
295  }
296  }
297 
298  for(unsigned int i=0; i < (*fLocalMCPoints).size(); i++){
299  nMCPointsVsRow[(*fLocalMCPoints)[i].IRow()]++;
300  }
301  for(int i=0; i < fTracker->NStations(); i++){
302  if ( nMCPointsVsRow[i] > 0 ) {
303  GetHisto("nHitsOverNMCPointsVsRow")->Fill( i, float(nHitsVsRow[i])/float(nMCPointsVsRow[i]) );
304  GetHisto("nHitsOverNMCPointsVsNMCTracks")->Fill( Multiplicity, float(nHitsVsRow[i])/float(nMCPointsVsRow[i]) );
305  }
306  }
307 
308  for(int iRTr=0; iRTr < nRecoTracks; iRTr++){ // TODO: make common
309  PndCAPerformanceRecoTrackData &recoD = recoData[iRTr];
310 
311  const PndCAGBTrack &recoTr = fTracker->Track( iRTr ); // TODO: make common
312  PndCAMCTrack &mcTr = (*fMCTracks)[ recoD.GetMCTrackId() ];
313 
314 #ifdef USE_CA_FIT // use 3 iterational fit, which ends on inner station
315 // PndCATrackParam param = recoTr.OuterParam();
316  PndCATrackParam param = recoTr.InnerParam();
317 #else
318  PndCATrackParam param = recoTr.InnerParam();
319 // PndCATrackParam param = recoTr.OuterParam();
320 #endif
321 
322  double RecoPt = 1. / fabs(param.QPt());
323  double RecoMom = RecoPt * sqrt(1. + param.DzDs()*param.DzDs());
324  // fNVsMom->Fill( param.GetY());
325  // fLengthVsMom->Fill( param.GetY(), t.NHits());
326  double prob = param.GetChi2() != -1 ? TMath::Prob( param.GetChi2(), param.GetNDF() ) : -1;
327 
328  GetHisto("purity")->Fill( recoData[iRTr].GetPurity() );
329  if ( recoD.IsGhost(PParameters::MinTrackPurity) ) {
330  GetHisto("ghostsLength")->Fill( recoTr.NHits() );
331  GetHisto("ghostsRMom")->Fill( RecoMom );
332  GetHisto("ghostsMCMom")->Fill( mcTr.P() );
333  GetHisto("ghostsRPt")->Fill( RecoPt );
334  GetHisto("ghostsMCPt")->Fill( mcTr.Pt() );
335  GetHisto("ghostsLengthAndMCMom")->Fill( recoTr.NHits(), mcTr.P() );
336  GetHisto("ghostsLengthAndRMom")->Fill( recoTr.NHits(), RecoMom );
337  GetHisto("ghostsChi2")->Fill( param.GetChi2() );
338  GetHisto("ghostsProb")->Fill( prob );
339  }
340  else {
341  GetHisto("recosLength")->Fill( recoTr.NHits() );
342  GetHisto("recosRMom")->Fill( RecoMom );
343  GetHisto("recosMCMom")->Fill( mcTr.P() );
344  GetHisto("recosRPt")->Fill( RecoPt );
345  GetHisto("recosMCPt")->Fill( mcTr.Pt() );
346  GetHisto("recosLengthAndMCMom")->Fill( recoTr.NHits() , mcTr.P() );
347  GetHisto("recosLengthAndRMom")->Fill( recoTr.NHits() , RecoMom );
348  GetHisto("recosChi2")->Fill( param.GetChi2() );
349  GetHisto("recosProb")->Fill( prob );
350  if ( abs(mcTr.P()) > PParameters::RefThreshold )
351  GetHisto("recosRefProb")->Fill( prob );
352  if( mcTr.P() > 0.5 ) GetHisto("nHitsRecoTOverNHitsMCT")->Fill( float(recoTr.NHits()) / mcTrackNRecoHits[recoD.GetMCTrackId()] );
353  }
354  }
355 } // void PndCAGlobalPerformance::FillHistos()
356 
357 #endif //DO_TPCCATRACKER_EFF_PERFORMANCE
358 
float Pz() const
Definition: PndCAMCTrack.h:45
static const float ExtraThreshold
static const float MinTrackPurity
float QPt() const
float DzDs() const
int NMCRows() const
Definition: PndCAMCTrack.h:60
int IRow() const
Definition: PndCAGBHit.h:57
Int_t i
Definition: run_full.C:25
float P() const
Definition: PndCAMCTrack.h:46
void SetSet(int v)
Definition: PndCAMCTrack.h:75
float Py() const
Definition: PndCAMCTrack.h:44
int FirstHitRef() const
Definition: PndCAGBTrack.h:31
friend F32vec4 sqrt(const F32vec4 &a)
Definition: P4_F32vec4.h:29
static T Sqrt(const T &x)
Definition: PndCAMath.h:37
virtual void EfficiencyPerformance()
int NHits() const
Definition: PndCAMCTrack.h:53
int NHitRows() const
Definition: PndCAMCTrack.h:61
virtual void CheckMCTracks()
Efficiency.
int NMCPoints() const
Definition: PndCAMCTrack.h:54
virtual void MatchTracks()
int MotherId() const
Definition: PndCAMCTrack.h:35
void SetNReconstructed(int v)
Definition: PndCAMCTrack.h:74
Information about reconstruction of MCTrack.
friend F32vec4 fabs(const F32vec4 &a)
Definition: P4_F32vec4.h:47
float Pt() const
Definition: PndCAMCTrack.h:47
bool IsGhost(float minPurity=0)
Information about reconstruction of Reconstructed Track.
float GetChi2() const
int reco()
PndSdsMCPoint * hit
Definition: anasim.C:70
virtual void FillHistos()
Histograms.
void SetNTurns(int v)
Definition: PndCAMCTrack.h:76
const std::vector< int > & RecoTrackIds() const
virtual void SetNewEvent(const PndCAGBTracker *const Tracker, vector< PndCAHitLabel > *hitLabels, vector< PndCAMCTrack > *mcTracks, vector< PndCALocalMCPoint > *localMCPoints)
int NHits() const
Definition: PndCAGBTrack.h:30
const PndCATrackParam & InnerParam() const
Definition: PndCAGBTrack.h:33
bool IsLeft() const
Definition: PndCAGBHit.h:65
float Px() const
Definition: PndCAMCTrack.h:43
int ID() const
Definition: PndCAGBHit.h:58
int GetNDF() const
static const float RefThreshold