PandaRoot
PndSttMvdGemTracking.h
Go to the documentation of this file.
1 //****************************************************************************
2 //* This file is part of PandaRoot. *
3 //* *
4 //* PandaRoot is distributed under the terms of the *
5 //* GNU General Public License (GPL) version 3, *
6 //* copied verbatim in the file "LICENSE". *
7 //* *
8 //* Copyright (C) 2006 - 2024 FAIR GmbH and copyright holders of PandaRoot *
9 //* The copyright holders are listed in the file "COPYRIGHTHOLDERS". *
10 //* The authors are listed in the file "AUTHORS". *
11 //****************************************************************************
12 
13 #ifndef PNDSTTMVDGEMTRACKING_H
14 #define PNDSTTMVDGEMTRACKING_H 1
15 
16 #include "PndGeoSttPar.h"
17 #include "PndGemDigiPar.h"
18 #include "PndGemHit.h"
19 #include "PndSttHit.h"
20 #include "PndTrackCand.h"
21 #include "PndTrack.h"
22 #include "PndPersistencyTask.h"
23 #include "FairGeanePro.h"
24 
25 #include "TCanvas.h"
26 #include "TH2F.h"
27 
28 class TClonesArray;
29 
31 
32  public:
35 
36  PndSttMvdGemTracking(Int_t verbose);
37 
40 
42  virtual InitStatus Init();
43 
45  virtual void Exec(Option_t *opt);
46 
48  void SetPersistenc(Bool_t persistence) { SetPersistency(persistence); }
49 
50  virtual void SetParContainers();
51 
52  // -------- MVD + STT ---> + GEM --------------
53  void Copy(PndTrackCand *completeCand, PndTrack *completeTrack, PndTrackCand *sttmvdCand, PndTrack *sttmvd);
54  void SetTimes(Int_t times) { fTimes = times; }
55  void SetMaxDistance(Double_t maxdistance) { fMaxDistance = maxdistance; }
56  void SwitchOnDisplay() { fDisplayOn = true; }
57 
58  void SetupGEMPlanes();
59  void Reset(Int_t nhits, Int_t ntracks);
60  void OrderGemHits(Int_t nhits);
61 
62  Int_t GetTrackIndex(Int_t i);
63  Int_t GetHitIndex(Int_t i);
64  Int_t GetPosIndex(PndGemHit *hit);
65 
66  Int_t CountTracks();
67  Int_t CountHitsInTrack(Int_t itrk);
68  void AddHitToTrack(Int_t ihit, Int_t itrk);
69  void DeleteHitFromTrack(Int_t ihit, Int_t itrk);
70  std::vector<int> GetTracksAssociatedToHit(Int_t ihit);
71  std::vector<int> GetHitsAssociatedToTrack(Int_t itrk);
72  std::vector<int> GetHitsAssociatedToTrackOnPlane(Int_t itrk, Int_t ipos);
73  void AddRemainingHits(Int_t ntracks);
74 
75  void CheckCombinatorial(Int_t nhits, Int_t ntracks);
76  void ForbidMultiAssignedHits(Int_t nhits, Int_t ntracks);
77  void OnlyOneHitToEachTrack(Int_t nhits, Int_t ntracks);
78  void Retrack();
79 
80  Bool_t PropagateToGemPlane(FairTrackParP *tmppar, FairTrackParP *gempar, Int_t ipos);
81  Bool_t PropagateToGemPlaneAsHelix(PndTrack *sttmvd, FairTrackParP *gempar, Int_t ipos);
82 
83  Double_t IsAssignable(FairTrackParP *gempar, PndGemHit *gemhit);
84  std::vector<int> AssignHits(Int_t itrk, FairTrackParP *gempar, Int_t ipos);
85  void EvaluatePerformances(Int_t nhits, Int_t ntracks);
86  void FillTrueDistances();
87  // -----------------------------------------------
88  Int_t SelectPdgCode(PndTrack *sttmvd);
89 
90  void Kalman(TMatrixT<double> extrap, TMatrixT<double> measurement, TMatrixT<double> H, TMatrixT<double> extrap_cov, TMatrixT<double> measurement_cov, TMatrixT<double> &kalman,
91  TMatrixT<double> &kalman_cov);
92  FairTrackParP SetStartParameters(PndTrack *sttmvd, PndTrackCand *sttmvdCand);
93 
94  Int_t GetClosestOnFirst(FairTrackParP *gempar, Int_t ipos, Double_t &closestdistance);
95 
96  // PREFIT -----------------------------
97  Bool_t Prefit(PndTrack *sttmvdTrack, PndTrackCand *sttmvdCand, TVector3 &lastpos, TVector3 &lastmom);
98  Bool_t IntersectionFinder(Double_t xc, Double_t yc, Double_t radius, PndSttHit *stthit, TVector3 &xyz, TVector3 &dxyz);
99  Bool_t Fit(TMatrixT<double> points, Double_t &outxc, Double_t &outyc, Double_t &outradius);
100  Bool_t ZFit(TMatrixT<double> points, Int_t charge, Double_t xc, Double_t yc, Double_t radius, Double_t &fitm, Double_t &fip);
101  Bool_t GetInitialParams(PndTrack *sttmvd, Double_t &xc, Double_t &yc, Double_t &radius, Double_t &fitm, Double_t &fitp);
102  Double_t CalculatePhi(TVector2 v, TVector2 p, double alpha, double Phi0, int charge);
103  Double_t CompareToPreviousPhi(Double_t Fi, Double_t Fi_pre, int charge);
104  Bool_t ZFind(Int_t nhits, TMatrixT<double> points, Double_t xc, Double_t yc, Double_t radius);
105 
106  // COMBINATORIAL EFFECS ----------------
107  void ConsiderCombinatorialEffect(Int_t nhits);
108  void SetCombinatorialDistance(Double_t combidistance) { fCombiDistance = combidistance; }
109  void UpdateMCTrackId(PndTrackCand *completeCand);
110 
111  // CHECK delete this when everything is ok
112  void UseMonteCarlo() { fUseMC = kTRUE; }
113  void WriteHistograms();
114 
115  void SetEvaluateFlag(Bool_t flag) { fEvaluate = flag; }
116 
117  // PDG stuff
118  Int_t GetPdgFromMC(int trackid);
119  Int_t GetChargeCorrectedPdgFromMC(int trackid, int charge);
120  void SetPdgFromMC() { fPdgFromMC = kTRUE; }
121  void SetPdgFromMC(TString trackid)
122  {
123  fStartTrackIDBranchName = trackid;
124  fPdgFromMC = kTRUE;
125  }
126  void SetDefaultPdg(int pdg) { fDefaultPdgCode = pdg; }
127 
128  void SetBranchNames(TString mvdpixel, TString mvdstrip, TString stt, TString gem);
129  void SetTrackBranchNames(TString startingtrack, TString startingtrackcand);
130 
131  private:
133  TClonesArray *fGemHitArray;
135  TClonesArray *fGemPointArray;
136 
138  TClonesArray *fMvdPixelHitArray;
140  TClonesArray *fMvdStripHitArray;
142  TClonesArray *fSttHitArray;
143 
145  TClonesArray *fMvdPointArray;
147  TClonesArray *fSttPointArray;
148 
150  TClonesArray *fMCTrackArray;
152  TClonesArray *fTrackArray;
154  TClonesArray *fTrackIDArray;
156  TClonesArray *fTrackCandArray;
157 
159  TClonesArray *fCompleteTrackCandArray;
161  TClonesArray *fCompleteTrackArray;
162 
164  TClonesArray *fTubeArray;
165 
166  Bool_t fPdgFromMC;
167 
168  PndGeoSttPar *fSttParameters;
169  PndGemDigiPar *fGemParameters;
170 
172  FairGeanePro *fPro;
173 
174  Int_t countgood, countbad, countdoubt, countreconstructablehit, countprinotassigned, countsecnotassigned, countreconstructablepoint;
175  Int_t countplane[8], countprop[2][8], countsttmvd, countsttmvdusable, countnohitonplane[8], evt;
176 
177  std::vector<int> usabletracks;
178 
180  std::vector<int> fOrdering;
181  std::vector<int>::iterator fOrderingIterator;
182 
184  TClonesArray *fSensPositions;
185 
187  Int_t fNPositions;
188 
190  TVector3 fDj, fDk;
191 
193  Int_t fPdgCode, fDefaultPdgCode;
194 
195  // CHECK added -------------------------
196  std::vector<std::pair<int, int>> trackvector;
197  std::vector<int> trackindexes;
198  std::vector<int> notassignedhits;
199  std::vector<int> notassignedtracks;
200 
202  std::map<int, bool> towhichplane;
203 
205  std::map<int, bool> fProTracks;
206 
212  TMatrixT<float> hitmap;
213 
217  TMatrixT<float> hitcounter;
218 
220  TMatrixT<double> distancemap;
221  // ----------------------------------------
222 
223  Double_t fMaxDistance;
224  Int_t fTimes;
225  Int_t fTurn;
226 
227  // DISPLAY
228  Bool_t fDisplayOn;
229  TCanvas *display;
230  TH2F *h[8], *hnotskewed, *hskewed;
231  TH1F *hdist[8], *hdist2[8], *hsigma[8], *hsigma2[8], *hchosen[8], *hchosen2[8], *hmcdist[8], *hmcx[8], *hmcy[8];
232 
233  // MonteCarlo
234  // CHECK delete this when everything is ok
235  Bool_t fUseMC;
236 
237  TString fMvdPixelBranchName, fMvdStripBranchName, fSttBranchName, fGemBranchName;
238  TString fStartTrackBranchName, fStartTrackCandBranchName, fStartTrackIDBranchName;
239 
241  std::map<int, int> fCombiMap;
242  Double_t fCombiDistance;
243 
244  // evaluate performances? yes/no
245  Bool_t fEvaluate;
246 
247  ClassDef(PndSttMvdGemTracking, 1);
248 };
249 
250 #endif
Double_t CompareToPreviousPhi(Double_t Fi, Double_t Fi_pre, int charge)
void OrderGemHits(Int_t nhits)
Int_t GetPosIndex(PndGemHit *hit)
void Copy(PndTrackCand *completeCand, PndTrack *completeTrack, PndTrackCand *sttmvdCand, PndTrack *sttmvd)
void SetCombinatorialDistance(Double_t combidistance)
std::vector< int > GetTracksAssociatedToHit(Int_t ihit)
void CheckCombinatorial(Int_t nhits, Int_t ntracks)
virtual void SetParContainers()
void SetBranchNames(TString mvdpixel, TString mvdstrip, TString stt, TString gem)
Bool_t Fit(TMatrixT< double > points, Double_t &outxc, Double_t &outyc, Double_t &outradius)
Double_t IsAssignable(FairTrackParP *gempar, PndGemHit *gemhit)
void SetPersistency(Bool_t val=kTRUE)
Double_t CalculatePhi(TVector2 v, TVector2 p, double alpha, double Phi0, int charge)
Int_t CountHitsInTrack(Int_t itrk)
Digitization Parameter Class for GEM part.
Definition: PndGemDigiPar.h:42
Bool_t IntersectionFinder(Double_t xc, Double_t yc, Double_t radius, PndSttHit *stthit, TVector3 &xyz, TVector3 &dxyz)
Int_t GetChargeCorrectedPdgFromMC(int trackid, int charge)
Bool_t ZFit(TMatrixT< double > points, Int_t charge, Double_t xc, Double_t yc, Double_t radius, Double_t &fitm, Double_t &fip)
__m128 v
Definition: P4_F32vec4.h:15
unsigned int i
Definition: P4_F32vec4.h:33
Int_t GetClosestOnFirst(FairTrackParP *gempar, Int_t ipos, Double_t &closestdistance)
void EvaluatePerformances(Int_t nhits, Int_t ntracks)
void AddRemainingHits(Int_t ntracks)
Bool_t PropagateToGemPlane(FairTrackParP *tmppar, FairTrackParP *gempar, Int_t ipos)
Int_t GetPdgFromMC(int trackid)
virtual void Exec(Option_t *opt)
void SetMaxDistance(Double_t maxdistance)
Int_t GetHitIndex(Int_t i)
Bool_t Prefit(PndTrack *sttmvdTrack, PndTrackCand *sttmvdCand, TVector3 &lastpos, TVector3 &lastmom)
void SetPdgFromMC(TString trackid)
void Kalman(TMatrixT< double > extrap, TMatrixT< double > measurement, TMatrixT< double > H, TMatrixT< double > extrap_cov, TMatrixT< double > measurement_cov, TMatrixT< double > &kalman, TMatrixT< double > &kalman_cov)
void DeleteHitFromTrack(Int_t ihit, Int_t itrk)
void ForbidMultiAssignedHits(Int_t nhits, Int_t ntracks)
void AddHitToTrack(Int_t ihit, Int_t itrk)
Bool_t ZFind(Int_t nhits, TMatrixT< double > points, Double_t xc, Double_t yc, Double_t radius)
Bool_t PropagateToGemPlaneAsHelix(PndTrack *sttmvd, FairTrackParP *gempar, Int_t ipos)
void SetTrackBranchNames(TString startingtrack, TString startingtrackcand)
virtual InitStatus Init()
std::vector< int > GetHitsAssociatedToTrack(Int_t itrk)
std::vector< int > AssignHits(Int_t itrk, FairTrackParP *gempar, Int_t ipos)
FairTrackParP SetStartParameters(PndTrack *sttmvd, PndTrackCand *sttmvdCand)
std::vector< int > GetHitsAssociatedToTrackOnPlane(Int_t itrk, Int_t ipos)
Int_t SelectPdgCode(PndTrack *sttmvd)
void OnlyOneHitToEachTrack(Int_t nhits, Int_t ntracks)
double alpha
Definition: f_Init.h:31
void UpdateMCTrackId(PndTrackCand *completeCand)
void SetPersistenc(Bool_t persistence)
Bool_t GetInitialParams(PndTrack *sttmvd, Double_t &xc, Double_t &yc, Double_t &radius, Double_t &fitm, Double_t &fitp)
void SetTimes(Int_t times)
void SetEvaluateFlag(Bool_t flag)
void ConsiderCombinatorialEffect(Int_t nhits)
Int_t GetTrackIndex(Int_t i)
void Reset(Int_t nhits, Int_t ntracks)