FairRoot/PandaRoot
PndGenfitAdapters.cxx
Go to the documentation of this file.
1 #include"PndGenfitAdapters.h"
2 
3 #include <iostream>
4 
5 #include"GFTrack.h"
6 #include"GFAbsTrackRep.h"
7 #include"GFTrackCand.h"
8 #include"PndTrack.h"
9 #include"PndTrackCand.h"
10 #include"GFDetPlane.h"
11 #include"GFException.h"
12 #include"TMatrixT.h"
13 #include"FairTrackParP.h"
14 
15 #include"GeaneTrackRep.h"
16 #include"RKTrackRep.h"
17 #include <cmath>
18 
20  PndTrackCand* retVal = new PndTrackCand();
21  unsigned int nhits = cand->getNHits();
22  unsigned detId,hitId;
23  double rho;
24  for(unsigned int i=0;i<nhits;++i){
25  cand->getHit(i,detId,hitId,rho);
26  retVal->AddHit(detId,hitId,rho);
27  }
28  retVal->setMcTrackId(cand->getMcTrackId());
29  return retVal;
30 }
31 
33  GFTrackCand* retVal = new GFTrackCand();
34  unsigned int nhits = cand->GetNHits();
35  for(unsigned int i=0;i<nhits;++i){
36  PndTrackCandHit candHit = cand->GetSortedHit(i);
37  retVal->addHit(candHit.GetDetId(),candHit.GetHitId(),candHit.GetRho(),i);
38  }
39  retVal->setMcTrackId(cand->getMcTrackId());
40  return retVal;
41 }
42 
44  GFAbsTrackRep* clone = tr->getCardinalRep()->clone();
45  TMatrixT<double> firstState = clone->getFirstState();
46  TMatrixT<double> lastState = clone->getLastState();
47  TMatrixT<double> firstCov = clone->getFirstCov();
48  TMatrixT<double> lastCov = clone->getLastCov();
49  GFDetPlane firstPlane = clone->getFirstPlane();
50  GFDetPlane lastPlane = clone->getLastPlane();
51 
52  GFAbsTrackRep* gtr;
53  if (dynamic_cast<GeaneTrackRep*>(clone)!=NULL)
54  gtr = dynamic_cast<GeaneTrackRep*>(clone);
55  else if (dynamic_cast<RKTrackRep*>(clone)!=NULL)
56  gtr = dynamic_cast<RKTrackRep*>(clone);
57  else {
58  std::cerr << " GenfitGFAbsTrackRep2PndTrack() can currently only handle GeaneTrackRep and RKTrackRep" << std::endl;
59  throw;
60  }
61 
62  //make the FairTrackParP for first and last hit
63  double firstCova[15];
64  int count=0;;
65  for(int i=0; i<5;++i){
66  for(int j=i;j<5;++j){
67  firstCova[count++]=firstCov[i][j];
68  }
69  }
70  double lastCova[15];
71  count=0;;
72  for(int i=0; i<5;++i){
73  for(int j=i;j<5;++j){
74  lastCova[count++]=lastCov[i][j];
75  }
76  }
77 
78  // calculation of spu = sign[p(DJ x DK)]
79  double first_pro(0), last_pro(0), first_spu, last_spu;
80  bool exc(false);
81 
82  try{
83  first_pro = gtr->getMom(firstPlane).Dot(firstPlane.getNormal());
84  last_pro = gtr->getMom(lastPlane).Dot(lastPlane.getNormal());
85  }
86  catch (GFException& e){
87  exc=true;
88  std::cerr<<"*** PndGenfitAdapters::GenfitTrack2PndTrack" << "\t" << "could not convert GenfitTrack to PndTrack"<<std::endl;
89  e.what();
90  }
91 
92  first_spu = (first_pro>=0) ? 1 : -1;
93  last_spu = (last_pro>=0) ? 1 : -1;
94 
95  FairTrackParP first(firstState[3][0],firstState[4][0],firstState[1][0],firstState[2][0],firstState[0][0],firstCova,firstPlane.getO(),firstPlane.getU(),firstPlane.getV(),first_spu);
96  FairTrackParP last(lastState[3][0],lastState[4][0],lastState[1][0],lastState[2][0],lastState[0][0],lastCova,lastPlane.getO(),lastPlane.getU(),lastPlane.getV(),last_spu);
97 
98  //copy the trackCand
99  GFTrackCand genfitCand = tr->getCand();
100  PndTrackCand* pndCand = GenfitTrackCand2PndTrackCand(&genfitCand);
101  PndTrack* retVal = new PndTrack(first,last,*pndCand);
102  retVal->SetChi2(tr->getChiSqu());
103  retVal->SetNDF(tr->getNDF());
104  if (tr->getNDF()==0)
105  {
106  retVal->SetFlag(-15);
107  }
108  else if (exc) {
109  retVal->SetFlag(-1);
110  }
111  else {
112  retVal->SetFlag(1);
113  }
114 
115  delete pndCand;
116  return retVal;
117 
118 }
Base Class for genfit track representations. Defines interface for track parameterizations.
Definition: GFAbsTrackRep.h:80
GFDetPlane getFirstPlane() const
int getMcTrackId() const
Definition: PndTrackCand.h:60
TVector3 getV() const
Definition: GFDetPlane.h:77
Int_t i
Definition: run_full.C:25
Detector plane genfit geometry class.
Definition: GFDetPlane.h:59
Track object for genfit. genfit algorithms work on these objects.
Definition: GFTrack.h:60
unsigned int getNDF() const
Get NDF.
Definition: GFTrack.h:258
double getChiSqu() const
Get chi2.
Definition: GFTrack.h:252
PndTrackCandHit GetSortedHit(UInt_t i)
Definition: PndTrackCand.h:54
void setMcTrackId(int i)
Definition: PndTrackCand.h:72
virtual const char * what() const
standard error message handling for exceptions. use like "std::cerr << e.what();" ...
Definition: GFException.cxx:47
unsigned int getNHits() const
Definition: GFTrackCand.h:113
void SetChi2(Double_t d)
Definition: PndTrack.h:39
virtual TVector3 getMom(const GFDetPlane &pl)=0
TMatrixT< double > getFirstCov() const
void setMcTrackId(int i)
set the MCT track id, for MC simulations
Definition: GFTrackCand.h:149
void AddHit(UInt_t detId, UInt_t hitId, Double_t rho)
void getHit(unsigned int i, unsigned int &detId, unsigned int &hitId) const
Get detector ID and cluster index (hitId) for hit number i.
Definition: GFTrackCand.h:84
Double_t GetRho() const
TMatrixT< double > getLastState() const
TMatrixT< double > getLastCov() const
TMatrixT< double > getFirstState() const
Track candidate – a list of cluster indices.
Definition: GFTrackCand.h:55
GFTrackCand * PndTrackCand2GenfitTrackCand(PndTrackCand *cand)
int getMcTrackId() const
get the MCT track id, for MC simulations - def. value -1
Definition: GFTrackCand.h:129
GFDetPlane getLastPlane() const
GFAbsTrackRep * getCardinalRep() const
Get cardinal track representation.
Definition: GFTrack.h:202
const GFTrackCand & getCand() const
Definition: GFTrack.h:142
Exception class for error handling in GENFIT (provides storage for diagnostic information) ...
Definition: GFException.h:45
PndTrack * GenfitTrack2PndTrack(const GFTrack *tr)
void addHit(unsigned int detId, unsigned int hitId, double rho=0., unsigned int planeId=0)
Definition: GFTrackCand.cxx:44
void SetFlag(Int_t i)
Definition: PndTrack.h:38
void SetNDF(Int_t i)
Definition: PndTrack.h:40
UInt_t GetNHits() const
Definition: PndTrackCand.h:59
TVector3 getU() const
Definition: GFDetPlane.h:76
int count
Int_t GetHitId() const
TVector3 getNormal() const
Definition: GFDetPlane.cxx:138
Int_t GetDetId() const
virtual GFAbsTrackRep * clone() const =0
TVector3 getO() const
Definition: GFDetPlane.h:75
PndTrackCand * GenfitTrackCand2PndTrackCand(const GFTrackCand *cand)