FairRoot/PandaRoot
PndRiemannTrack.h
Go to the documentation of this file.
1 //-----------------------------------------------------------
2 // File and Version Information:
3 // $Id$
4 //
5 // Description:
6 // Track on Riemann Sphere
7 // Circle parameters can be calculated from plane parameters
8 // plane(c,nx,ny,nz);
9 //
10 //
11 // Environment:
12 // Software developed for the PANDA Detector at FAIR.
13 //
14 // Author List:
15 // Sebastian Neubert TUM (original author)
16 //
17 //
18 //-----------------------------------------------------------
19 
20 #ifndef PNDRIEMANNTRACK_HH
21 #define PNDRIEMANNTRACK_HH
22 
23 // Base Class Headers ----------------
24 #include "TObject.h"
25 
26 // Collaborating Class Headers -------
27 #include <vector>
28 #include "TVectorD.h"
29 #include "TVector3.h"
30 #include "TMatrixD.h"
31 #include "TMath.h"
32 
33 // Collaborating Class Declarations --
34 #include "PndRiemannHit.h"
35 #include "PndTrack.h"
36 #include "PndTrackCand.h"
37 
38 #include "PndSttHit.h"
39 #include "PndSttTube.h"
40 
41 #include <iostream>
42 #include <iomanip>
43 #include <algorithm>
44 
45 
46 class PndRiemannTrack : public TObject{
47 public:
48 
49  // Constructors/Destructors ---------
51  PndRiemannTrack(PndTrackCand* trackCand);
53 
54 
55  // Accessors -----------------------
56  const TVectorD& n() const{return fn;}
57  double c() const {return fc;}
58  const TVectorD& av() const{return fav;}
59 
60  TVectorD orig() const;
61  double dX(){return TMath::Sqrt(fcovRXY[1][1]);}
62  double dY(){return TMath::Sqrt(fcovRXY[2][2]);}
63  double r() const;
64  double dR();
65  double dip();
66  double dipangle(); //< dipangle theta
67  double dDip();
68  double Pt(double B); //< transvers momentum Pt calculated by the magnetic field B [tesla]
69  double Pl(double B);
70  double P(double B);
71  double sign() const;
72  double getSZm() const {return fm;}
73  double getSZt() const {return ft;}
74  unsigned int getNumHits() {return fHits.size();}
75  PndRiemannHit* getHit(unsigned int i) {PndRiemannHit* myHit = &(fHits[i]); return myHit;}
77  std::vector<PndRiemannHit> getHits() const{return fHits;};
78  TVector3 getPforHit(int i, double B);
80  FairTrackParP getTrackParPForHit(Int_t i, Double_t B);
81  Int_t getCharge(Double_t B);
82 
83  double calcZPosByS(double s);
84  TVector3 calcPosByS(double s);
85  void calcStartStopAlpha();
86  double calcAlpha(PndRiemannHit* myHit);
87 
88  void calcSForHits();
89 
90  int calcIntersection(PndRiemannTrack& track, TVector3& p1, TVector3& p2);
91 
92  void sortHits(){ std::sort(fHits.begin(), fHits.end());}
93 
94  double weight() const {return fweight;};
95  TMatrixD covPlane() const {return fcovPlane;};
96  TMatrixD jacRXY() const {return fjacRXY;};
97  TMatrixD covRXY() const {return fcovRXY;};
98  double m() const {return fm;}
99  double mError() const {return fmError;}
100  double t() const {return ft;}
101  double tError() const {return ftError;}
102 
103  // Modifiers -----------------------
104  void addHit(PndRiemannHit& hit);
105  void addPndTrackCand(PndTrackCand* trackCand);
106  void init(double x0, double y0, double R,
107  double dip, double z0);
108 
109  // Operations ----------------------
110  void refit(bool withErrorCalc = true);
111  double dist(PndRiemannHit* hit);
112  double distError(PndRiemannHit* hit);
113  double distCircle(PndRiemannHit* hit);
114  double ChiSquareDistCircle();
115  void szFit(bool withErrorCalc = true);
116  double calcChi2Plane();
117  double calcSZChi2(PndRiemannHit* hit); //calculates the chi2 of the track plus the additional hit
118  double szDist(PndRiemannHit* hit);
119  double szError(PndRiemannHit* hit);
120  double szChi2() const{return fChi2;};
121 
122  void SetVerbose(int i){fVerbose = i;}
123  void SetVertexCut(double cut){fVertexCut = cut;}
124  TVector3 calcErrorPosByS(Double_t s, Double_t dS);
125 
126 
127  void correctSttHits();
130  void PrintHits();
131 
132  virtual void Print(std::ostream& out = std::cout){
133 
134  out << std::setprecision(6) << "Riemann Track: Radius " << r() << " +/- " << dR() << " Origin: " << orig()[0] << " +/- " << dX() << " / " << orig()[1] << " +/- " << dY() << std::endl;
135  out << "RiemannTrack: Normal: " << n()[0] << "/" << n()[1] << "/" << n()[2] << " c: " << c() << std::endl;
136  out << "Dip: " << dip() << " +/- " << dDip() << " StartAlpha: " << fStartAlpha << " StopAlpha: " << fStopAlpha << std::endl;
137  out << "Chi2Plane: " << calcChi2Plane() << std::endl;
138  PrintHits();
139  }
140 
142  track.Print(out);
143  return out;
144  }
145 
146 private:
147 
148  // Private Data Members ------------
149  TVectorD fn;
150  TVectorD fav;
151  double fc;
152 
153  double fm;
154  double ft;
155  double fmError;
156  double ftError;
157  double fChi2;
161 
162  int fVerbose;
163 
164  bool fFitDone;
167  double fweight;
168  bool ftrefit;
169  double fVertexCut;
170 
171  std::vector<PndRiemannHit> fHits;
174 
175  std::map<TString, Int_t> fBranchNameMap;
178 
179 
180  void calcJacRXY();
181 
184  TVectorD calcErrorXY1XY2(TVector3& line, TVector3& dLine, TVector3& offset, TVector3& dOffset);
185  Double_t calcErrorS(TVector2& XY, TVector2& dXY, PndRiemannTrack* track);
186 
187  Int_t GetBranchId(TString branchName);
188 
189 
190 
191 
192  // Private Methods -----------------
193 
194 
195 public:
196  ClassDef(PndRiemannTrack,3)
197 
198 };
199 
200 #endif
201 
202 //--------------------------------------------------------------
203 // $Log$
204 //--------------------------------------------------------------
Double_t z0
Definition: checkhelixhit.C:62
Double_t x0
Definition: checkhelixhit.C:70
PndRiemannHit * getLastHit()
unsigned int getNumHits()
double fweight
sum over all weights (1/(sigmaXY*sigmaXY))
double calcSZChi2(PndRiemannHit *hit)
double distError(PndRiemannHit *hit)
PndRiemannHit correctSttSkewedHit(PndSttHit *mySttHit, PndSttTube *myTube)
void init(double x0, double y0, double R, double dip, double z0)
double calcZPosByS(double s)
Int_t i
Definition: run_full.C:25
PndRiemannHit correctSttHit(PndSttHit *mySttHit)
double szChi2() const
PndRiemannTrack track
Definition: RiemannTest.C:33
std::vector< PndRiemannHit > fHits
static T Sqrt(const T &x)
Definition: PndCAMath.h:37
double distCircle(PndRiemannHit *hit)
TLorentzVector s
Definition: Pnd2DStar.C:50
double szDist(PndRiemannHit *hit)
friend std::ostream & operator<<(std::ostream &out, PndRiemannTrack &track)
double P(double B)
double dist(PndRiemannHit *hit)
TVector3 offset(2, 0, 0)
Int_t GetBranchId(TString branchName)
double fmError
Error of fit.
std::map< TString, Int_t > fBranchNameMap
double fChi2
Chisquare of sz fit.
void calcJacRXY()
calcualtes fjacRXY
double getSZt() const
TVector3 calcErrorPosByS(Double_t s, Double_t dS)
double sign() const
TVector3 calcErrorLineOffset(PndRiemannTrack &track)
const TVectorD & n() const
double cut[MAX]
Definition: autocutx.C:36
TVector3 calcErrorLineNorm(PndRiemannTrack &track)
FairTrackParP getTrackParPForHit(Int_t i, Double_t B)
double m() const
double t() const
double tError() const
Int_t getCharge(Double_t B)
Double_t
void addPndTrackCand(PndTrackCand *trackCand)
void SetVertexCut(double cut)
TVectorD fn
normal vector to plane;
Double_t y0
Definition: checkhelixhit.C:71
TVectorD calcErrorXY1XY2(TVector3 &line, TVector3 &dLine, TVector3 &offset, TVector3 &dOffset)
void refit(bool withErrorCalc=true)
TMatrixD covRXY() const
TFile * out
Definition: reco_muo.C:20
void SetVerbose(int i)
double Pt(double B)
double getSZm() const
TPad * p2
Definition: hist-t7.C:117
double ChiSquareDistCircle()
double r() const
double weight() const
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
int calcIntersection(PndRiemannTrack &track, TVector3 &p1, TVector3 &p2)
TMatrixD covPlane() const
double mError() const
double c() const
virtual void Print(std::ostream &out=std::cout)
TVectorD orig() const
TPad * p1
Definition: hist-t7.C:116
TMatrixD fjacRXY
jacobian matrix to transform from c,n1,n2,n2 to r,x,y
void szFit(bool withErrorCalc=true)
void addHit(PndRiemannHit &hit)
Double_t calcErrorS(TVector2 &XY, TVector2 &dXY, PndRiemannTrack *track)
double fc
distance of plane to origin
TMatrixD jacRXY() const
TMatrixD fcovPlane
full covarince matrix of the plane;
TVector3 getPforHit(int i, double B)
std::vector< PndRiemannHit > getHits() const
double fm
parameters of sz-fit
const TVectorD & av() const
TVectorD fav
average over all hits
Double_t R
Definition: checkhelixhit.C:61
PndTrack getPndTrack(Double_t B)
double calcAlpha(PndRiemannHit *myHit)
TVector3 calcPosByS(double s)
PndRiemannHit * getHit(unsigned int i)
TMatrixT< double > TMatrixD
Definition: PndLmdDim.h:52
double ftError
Error of fit.
double szError(PndRiemannHit *hit)
double Pl(double B)