FairRoot/PandaRoot
PndSdsRecoHit.cxx
Go to the documentation of this file.
1 //-----------------------------------------------------------
2 // File and Version Information:
3 // $Id$
4 //
5 // Description:
6 // Implementation of class PndSdsRecoHit
7 // see PndSdsRecoHit.h for details
8 //
9 // Environment:
10 // Software developed for the PANDA Detector at FAIR.
11 //
12 // Author List:
13 // Sebastian Neubert TUM (original author)
14 // Ralf Kliemt TUD (application to PndSds)
15 //
16 //-----------------------------------------------------------
17 
18 // C/C++ Headers ----------------------
19 // root Headers ----------------------
20 #include "TMatrixT.h"
21 #include "TMath.h"
22 // Collaborating Class Headers --------
23 #include "FairMCPoint.h"
24 #include "GeaneTrackRep.h"
25 #include "RKTrackRep.h"
26 #include "GFDetPlane.h"
27 // This Class' Header ------------------
28 #include "PndSdsRecoHit.h"
29 #include "PndSdsMCPoint.h"
30 #include "PndSdsHit.h"
31 #include "PndGeoHandling.h"
32 #include "TGeoManager.h"
33 #include "FairRootManager.h"
34 
35 // Class Member definitions -----------
36 
38 
39 
41 {
42 }
43 
45 : GFRecoHitIfc<GFPlanarHitPolicy>(fNparHitRep)
46 {
47 }
48 
49 
51 : GFRecoHitIfc<GFPlanarHitPolicy>(fNparHitRep)
52 {
53  std::cout<<" -I- PndSdsRecoHit::PndSdsRecoHit(PndSdsMCPoint*) called."<<std::endl;
54 
55  fHitCoord[0][0] = point->GetX();
56  fHitCoord[1][0] = point->GetY();
57 
58  // we set the covariances to (50mu)^2 by hand.
59  fHitCov[0][0] = 0.0050 * 0.0050;//cm //TODO: cm is rigt?
60  fHitCov[1][1] = 0.0050 * 0.0050;//cm //TODO: cm is rigt?
61 
62  TVector3 o(0.,0.,point->GetZ()),
63  u(1.,0.,0.),
64  v(0.,1.,0.);
65 
67 
68 }
69 
71 : GFRecoHitIfc<GFPlanarHitPolicy>(fNparHitRep)
72 {
73 
74  // std::cout<<" -I- PndSdsRecoHit::PndSdsRecoHit(PndSdsHit*) called."<<std::endl;
75  // std::cout<<*hit<<std::endl;
76 
77  Int_t id = hit->GetSensorID();
78 
79  // FairRootManager* ioman = FairRootManager::Instance();
80  // TString fGeoFile = ioman->GetInFile()->GetName();
82  TString path = fGeoH->GetPath(id);
83  // std::cout<<"Detector path: "<<path.Data()<<std::endl;
84  TVector3 oo, uu, vv;
85  fGeoH->GetOUVShortId(id, oo,uu,vv);
86 
87  TVector3 position = hit->GetPosition();
88  TVector3 localpos = fGeoH->MasterToLocalShortId(position, id);
89 
90  fHitCoord[0][0] = localpos.X();
91  fHitCoord[1][0] = localpos.Y();
92 
93  TMatrixD cova = fGeoH->MasterToLocalErrorsShortId(hit->GetCov(), id);
94  // project only the 2 dimensions of cov.
95  fHitCov[0][0] = cova[0][0];
96  fHitCov[0][1] = cova[0][1];
97  fHitCov[1][0] = cova[1][0];
98  fHitCov[1][1] = cova[1][1];
99 
100 // std::cout<<" -I- PndSdsRecoHit::PndSdsRecoHit: Wrote a hit with"
101 // <<"\n(x,y) = ("<<localpos.X()<<","<<localpos.Y()<<")."
102 // <<"\nCovariance Matrix is";
103 // fHitCov.Print();
104 // std::cout<<"From 3D hit matrix";
105 // cova.Print();
106 
107  fPolicy.setDetPlane(GFDetPlane(oo,uu,vv));
108 }
109 //============================================================================
110 
111 
112 
113 TMatrixT<double>
115 {
116 
117  if (dynamic_cast<const RKTrackRep*>(stateVector) != NULL) {
118  // Uses TrackParP (q/p,v',w',v,w)
119  // coordinates are defined by detplane!
120  TMatrixT<double> HMatrix(2,5);
121 
122  HMatrix[0][0] = 0.;
123  HMatrix[0][1] = 0.;
124  HMatrix[0][2] = 0.;
125  HMatrix[0][3] = 1.;
126  HMatrix[0][4] = 0.;
127 
128  HMatrix[1][0] = 0.;
129  HMatrix[1][1] = 0.;
130  HMatrix[1][2] = 0.;
131  HMatrix[1][3] = 0.;
132  HMatrix[1][4] = 1.;
133 
134  return HMatrix;
135  }
136 
137  // !! TODO I copied this from the DemoRecoHit - check validity!!!
138  if (dynamic_cast<const GeaneTrackRep*>(stateVector) != NULL) {
139  // Uses TrackParP (q/p,v',w',v,w)
140  // coordinates are defined by detplane!
141  TMatrixT<double> HMatrix(fNparHitRep,5);
142 
143  HMatrix[0][0] = 0.;
144  HMatrix[0][1] = 0.;
145  HMatrix[0][2] = 0.;
146  HMatrix[0][3] = 1.;
147  HMatrix[0][4] = 0.;
148 
149  HMatrix[1][0] = 0.;
150  HMatrix[1][1] = 0.;
151  HMatrix[1][2] = 0.;
152  HMatrix[1][3] = 0.;
153  HMatrix[1][4] = 1.;
154  return HMatrix;
155  }
156  else {
157  std::cerr << "PndSdsRecoHit can only handle state"
158  << " vectors of type GeaneTrackRep -> abort"
159  << std::endl;
160  throw;
161  }
162 
163 }
164 
165 Double_t
167  const TMatrixT<Double_t>& ) // stateVector state// [R.K.03/2017] unused variable(s)
168 {
169  throw;
170 }
171 
Base Class for genfit track representations. Defines interface for track parameterizations.
Definition: GFAbsTrackRep.h:80
PndGeoHandling * fGeoH
Definition: PndSdsRecoHit.h:64
virtual TMatrixT< double > getHMatrix(const GFAbsTrackRep *stateVector)
Get transformation matrix. Transformation between hit coordinates and track representation coordinate...
void GetOUVShortId(Int_t shortId, TVector3 &o, TVector3 &u, TVector3 &v)
RecoHit interface template class. Provides comfortable interface to create RecoHits.
Definition: GFRecoHitIfc.h:60
TVector3 GetPosition() const
Definition: PndSdsHit.h:93
Detector plane genfit geometry class.
Definition: GFDetPlane.h:59
TMatrixT< double > fHitCoord
Vector of raw coordinates of hit.
Definition: GFAbsRecoHit.h:76
virtual Double_t residualScalar(GFAbsTrackRep *stateVector, const TMatrixT< Double_t > &state)
static const Int_t fNparHitRep
Definition: PndSdsRecoHit.h:63
void setDetPlane(const GFDetPlane &p)
Set physical detector plane. Needs to be called before hit can be used.
virtual ~PndSdsRecoHit()
TMatrixD MasterToLocalErrorsShortId(const TMatrixD &master, const Int_t &shortId)
TString GetPath(Int_t shortID)
for a given shortID the path is returned
__m128 v
Definition: P4_F32vec4.h:4
TMatrixD GetCov() const
Definition: PndSdsHit.h:98
Double_t
TClonesArray * point
Definition: anaLmdDigi.C:29
Policy class implementing a planar hit geometry.
static PndGeoHandling * Instance()
TVector3 MasterToLocalShortId(const TVector3 &master, const Int_t &shortId)
int hit(Int_t nEvents=0, TString inFile="sim.root", TString parFile="par.root", TString inDigi="digi.root", TString outFile="hit.root", Int_t timeBased=0)
Definition: hit.C:1
ClassImp(PndSdsRecoHit)
Int_t GetSensorID() const
Definition: PndSdsHit.h:90
TMatrixT< double > TMatrixD
Definition: PndLmdDim.h:52
TMatrixT< double > fHitCov
Covariance of raw hit coordinates.
Definition: GFAbsRecoHit.h:79