FairRoot/PandaRoot
PndRecoKalmanFit2.cxx
Go to the documentation of this file.
1 //-----------------------------------------------------------
2 // File and Version Information:
3 // $Id$
4 //
5 // Description:
6 // Implementation of class PndRecoKalmanFit
7 // see PndRecoKalmanFit.h for details
8 //
9 // Environment:
10 // Software developed for the PANDA Detector at FAIR.
11 //
12 // Author List:
13 // Stefano Spataro, UNI Torino
14 // modified by Elisabetta Prencipe, 19/05/2014
15 //-----------------------------------------------------------
16 
17 // Panda Headers ----------------------
18 
19 // This Class' Header ------------------
20 #include "PndRecoKalmanFit2.h"
21 
22 // C/C++ Headers ----------------------
23 #include <algorithm>
24 #include <iostream>
25 #include <assert.h>
26 #include <cmath>
27 
28 // Collaborating Class Headers --------
29 #include "FairRootManager.h"
30 #include "FairRuntimeDb.h"
31 #include "FairRunAna.h"
32 #include "TClonesArray.h"
33 
34 #include "Track.h"
35 #include "ConstField.h"
36 
37 #include "PndSdsRecoHit2.h"
38 #include "PndGemRecoHit2.h"
39 #include "PndSttRecoHit2.h"
40 #include "PndFtsRecoHit2.h"
41 #include "PndMdtRecoHit2.h"
42 #include "PndSttRecoHitProducer2.h"
43 #include "PndFtsRecoHitProducer2.h"
44 #include "PndGeoSttPar.h"
45 #include "PndGeoFtsPar.h"
46 #include "PndSttMapCreator.h"
47 #include "PndFtsMapCreator.h"
48 #include "PndGenfitAdapters2.h"
49 #include "PndTrack.h"
50 #include "PndTrackCand.h"
51 #include "PndDetectorList.h"
52 #include "PndGeoHandling.h"
53 #include "MeasurementFactory.h"
54 #include "KalmanFitter.h"
55 #include "Exception.h"
56 #include "TLorentzVector.h"
57 
58 #include "FairTrackParH.h"
59 
60 #include <MaterialEffects.h>
61 #include <TGeoMaterialInterface.h>
62 #include "RKTrackRep.h"
63 #include "FieldManager.h"
64 #include "PndGenfitField2.h"
65 #include "FairGeanePro.h"
66 
67 // Class Member definitions -----------
68 
69 
70 PndRecoKalmanFit2::PndRecoKalmanFit2(): TNamed("Genfit", "Fit Tracks"),
71  fMvdBranchName(""), fCentralTrackerBranchName(""),
72  fUseGeane(kTRUE), fPropagateToIP(kTRUE), fPropagateDistance(-1.f), fPerpPlane(kFALSE), fNumIt(1), fVerbose(0)
73 {
75 }
76 
78 {
79  //Get ROOT Manager
80  FairRootManager* ioman= FairRootManager::Instance();
81  if(ioman==0)
82  {
83  Error("PndRecoKalmanFit2::Init","RootManager not instantiated!");
84  return kFALSE;
85  }
86 
87  // STT map loading
88  FairRuntimeDb* rtdb = FairRunAna::Instance()->GetRuntimeDb();
89  PndGeoSttPar *sttParameters = (PndGeoSttPar*) rtdb->getContainer("PndGeoSttPar");
90  TClonesArray *tubeArray = NULL;
91  if(sttParameters->GetGeometryType() == 1) {
92  PndSttMapCreator *mapper = new PndSttMapCreator(sttParameters);
93  tubeArray = mapper->FillTubeArray();
94  }
95 
96  // FTS map loading
97  PndGeoFtsPar *ftsParameters = (PndGeoFtsPar*) rtdb->getContainer("PndGeoFtsPar");
98  TClonesArray *ftsTubeArray = NULL;
99  if(ftsParameters->GetGeometryType() != -1) {
100  PndFtsMapCreator *ftsMapper = new PndFtsMapCreator(ftsParameters);
101  ftsTubeArray = ftsMapper->FillTubeArray();
102  }
103 
104  // Build hit factory -----------------------------
106  if (fVerbose<2) genfit::Exception::quiet(true);
107 
108  TClonesArray* stripar; TClonesArray* pixelar;
109  if (fMvdBranchName == "Mix")
110  {
111  stripar=(TClonesArray*) ioman->GetObject("MVDHitsStripMix");
112  if(stripar!=0)
113  {
114  fTheRecoHitFactory->addProducer(FairRootManager::Instance()->GetBranchId("MVDHitsStripMix"),new genfit::MeasurementProducer<PndSdsHit,PndSdsRecoHit2>(stripar));
115  std::cout << "*** PndRecoKalmanFit2::Init" << "\t" << "MVDHitsStripMix array found" << std::endl;
116  }
117 
118  pixelar=(TClonesArray*) ioman->GetObject("MVDHitsPixelMix");
119  if(pixelar!=0)
120  {
121  fTheRecoHitFactory->addProducer(FairRootManager::Instance()->GetBranchId("MVDHitsPixelMix"),new genfit::MeasurementProducer<PndSdsHit,PndSdsRecoHit2>(pixelar));
122  std::cout << "*** PndRecoKalmanFit2::Init" << "\t" << "MVDHitsPixelMix array found" << std::endl;
123  }
124  }
125  else
126  {
127  stripar=(TClonesArray*) ioman->GetObject("MVDHitsStrip");
128  if(stripar!=0)
129  {
130  fTheRecoHitFactory->addProducer(FairRootManager::Instance()->GetBranchId("MVDHitsStrip"),new genfit::MeasurementProducer<PndSdsHit,PndSdsRecoHit2>(stripar));
131  std::cout << "*** PndRecoKalmanFit2::Init" << "\t" << "MVDHitsStrip array found" << std::endl;
132  }
133 
134  pixelar=(TClonesArray*) ioman->GetObject("MVDHitsPixel");
135  if(pixelar!=0)
136  {
137  fTheRecoHitFactory->addProducer(FairRootManager::Instance()->GetBranchId("MVDHitsPixel"),new genfit::MeasurementProducer<PndSdsHit,PndSdsRecoHit2>(pixelar));
138  std::cout << "*** PndRecoKalmanFit2::Init" << "\t" << "MVDHitsPixel array found" << std::endl;
139  }
140  }
141 
142  TClonesArray* sthit; //TClonesArray *sttr; //[R.K. 01/2017] unused variable?
143 
144  if (fCentralTrackerBranchName == "Mix")
145  {
146  sthit=(TClonesArray*) ioman->GetObject("STTHitMix");
147  if(sthit!=0)
148  {
149  fTheRecoHitFactory->addProducer(FairRootManager::Instance()->GetBranchId("STTHitMix"),new PndSttRecoHitProducer2<PndSttHit,PndSttRecoHit2>(sthit, tubeArray));
150  std::cout << "*** PndRecoKalmanFit2::Init" << "\t" << "STTHitMix array found" << std::endl;
151  }
152  }
153  else
154  {
155  sthit=(TClonesArray*) ioman->GetObject("STTHit");
156  if(sthit!=0)
157  {
158  fTheRecoHitFactory->addProducer(FairRootManager::Instance()->GetBranchId("STTHit"),new PndSttRecoHitProducer2<PndSttHit,PndSttRecoHit2>(sthit, tubeArray));
159  std::cout << "*** PndRecoKalmanFit2::Init" << "\t" << "SttHit array found" << std::endl;
160  }
161  }
162 
163  TClonesArray* gemar=(TClonesArray*) ioman->GetObject("GEMHit");
164  if(gemar!=0)
165  {
166  fTheRecoHitFactory->addProducer(FairRootManager::Instance()->GetBranchId("GEMHit"),new genfit::MeasurementProducer<PndGemHit,PndGemRecoHit2>(gemar));
167  std::cout << "*** PndRecoKalmanFit2::Init" << "\t" << "GEMHit array found" << std::endl;
168  }
169 
170  TClonesArray* mdtar=(TClonesArray*) ioman->GetObject("MdtHit");
171  if(mdtar!=0)
172  {
173  fTheRecoHitFactory->addProducer(FairRootManager::Instance()->GetBranchId("MdtHit"),new genfit::MeasurementProducer<PndMdtHit,PndMdtRecoHit2>(mdtar));
174  std::cout << "*** PndRecoKalmanFit2::Init" << "\t" << "MdtHit array found" << std::endl;
175  }
176 
177  TClonesArray* ftsar=(TClonesArray*) ioman->GetObject("FTSHit");
178  if(ftsar!=0)
179  {
180  fTheRecoHitFactory->addProducer(FairRootManager::Instance()->GetBranchId("FTSHit"),new PndFtsRecoHitProducer2<PndFtsHit,PndFtsRecoHit2>(ftsar, ftsTubeArray));
181  std::cout << "*** PndRecoKalmanFit2::Init" << "\t" << "FtsHit array found" << std::endl;
182  }
183 
184  if (fUseGeane)
185  {
186  fPro = new FairGeanePro();
187  if (fVerbose==0) fPro->SetPrintErrors(kFALSE);
188  }
189  else
190  {
191  Error("PndRecoKalmanFit2::Init","Only GEANE Propagation available!!!");
192  return kFALSE;
193  }
194 
196 
199 
200  std::cout << "===PndRecoKalmanFit2::Init() finished ===================================================" << std::endl;
201 
202  return kTRUE;
203 }
204 
205 
207 
208 PndTrack* PndRecoKalmanFit2::Fit(PndTrack *tBefore, Int_t PDG) {
209  PndTrack* tAfter = NULL;
210  if (fVerbose > 0)
211  std::cout << "PndRecoKalmanFit2::Fit" << std::endl;
212  if (fabs(tBefore->GetParamFirst().GetPz()) < 1e-9) {
213  tAfter = tBefore;
214  tAfter->SetFlag(-10);
215  std::cout << "tAfter = tBefore " << tBefore->GetParamFirst().GetPz()
216  << std::endl;
217  return tAfter; // flag -10 : pz==0
218  }
219 
220  //Int_t fCharge = tBefore->GetParamFirst().GetQ(); //[R.K. 01/2017] unused variable?
221  Int_t PDGCode = PDG;
222  TVector3 StartPos(tBefore->GetParamFirst().GetX(),
223  tBefore->GetParamFirst().GetY(), tBefore->GetParamFirst().GetZ());
224  TVector3 StartMom(tBefore->GetParamFirst().GetPx(),
225  tBefore->GetParamFirst().GetPy(), tBefore->GetParamFirst().GetPz());
226 
227  TMatrixDSym covSeed(6);
228  covSeed(0, 0) = tBefore->GetParamFirst().GetDX()
229  * tBefore->GetParamFirst().GetDX();
230  covSeed(1, 1) = tBefore->GetParamFirst().GetDY()
231  * tBefore->GetParamFirst().GetDY();
232  covSeed(2, 2) = tBefore->GetParamFirst().GetDZ()
233  * tBefore->GetParamFirst().GetDZ();
234 
235  covSeed(3, 3) = tBefore->GetParamFirst().GetDPx()
236  * tBefore->GetParamFirst().GetDPx();
237  covSeed(4, 4) = tBefore->GetParamFirst().GetDPy()
238  * tBefore->GetParamFirst().GetDPy();
239  covSeed(5, 5) = tBefore->GetParamFirst().GetDPz()
240  * tBefore->GetParamFirst().GetDPz();
241 
242  FairTrackParP par = tBefore->GetParamFirst();
243  Int_t ierr = 0;
244  FairTrackParH helix(&par, ierr);
245  FairGeanePro fPro0;
246  if (fVerbose == 0)
247  fPro0.SetPrintErrors(kFALSE);
248  FairTrackParH fRes;
249 
250 
251 
252  if (fPropagateToIP) {
253  // Calculating params at PCA to Origin
254 
255  fPro0.SetPoint(TVector3(0, 0, 0));
256  fPro0.PropagateToPCA(1, -1);
257  Bool_t rc = fPro0.Propagate(&helix, &fRes, PDGCode);
258  if (rc) {
259  StartPos.SetXYZ(fRes.GetX(), fRes.GetY(), fRes.GetZ());
260  StartMom.SetXYZ(fRes.GetPx(), fRes.GetPy(), fRes.GetPz());
261 
262  covSeed(0, 0) = fRes.GetDX() * fRes.GetDX();
263  covSeed(1, 1) = fRes.GetDY() * fRes.GetDY();
264  covSeed(2, 2) = fRes.GetDZ() * fRes.GetDZ();
265 
266  covSeed(3, 3) = fRes.GetDPx() * fRes.GetDPx();
267  covSeed(4, 4) = fRes.GetDPy() * fRes.GetDPy();
268  covSeed(5, 5) = fRes.GetDPz() * fRes.GetDPz();
269  }
270  } else if (fPropagateDistance > 0.f) {
271  // Calculating params at fPropagateDistance cm before the first hit
272 // FairTrackParP par = tBefore->GetParamFirst();
273 // Int_t ierr = 0;
274 // FairTrackParH *helix = new FairTrackParH(&par, ierr);
275 // FairGeanePro *fPro0 = new FairGeanePro();
276 // if (fVerbose == 0)
277 // fPro0->SetPrintErrors(kFALSE);
278 // FairTrackParH *fRes = new FairTrackParH();
279  fPro0.PropagateToLength(-fPropagateDistance);
280  std::cout << "PndRekoKalmanFit2::Fit helix: ";
281  helix.GetMomentum().Print();
282  helix.GetPosition().Print();
283  std::cout << std::endl;
284  Bool_t rc = fPro0.Propagate(&helix, &fRes, PDGCode);
285  if (rc) {
286  StartPos.SetXYZ(fRes.GetX(), fRes.GetY(), fRes.GetZ());
287  StartMom.SetXYZ(fRes.GetPx(), fRes.GetPy(), fRes.GetPz());
288 
289  covSeed(0, 0) = fRes.GetDX() * fRes.GetDX();
290  covSeed(1, 1) = fRes.GetDY() * fRes.GetDY();
291  covSeed(2, 2) = fRes.GetDZ() * fRes.GetDZ();
292 
293  covSeed(3, 3) = fRes.GetDPx() * fRes.GetDPx();
294  covSeed(4, 4) = fRes.GetDPy() * fRes.GetDPy();
295  covSeed(5, 5) = fRes.GetDPz() * fRes.GetDPz();
296  }
297  }
298 
299  TVector3 plane_v1, plane_v2;
300  if (fPerpPlane) {
301  plane_v1 = StartMom.Orthogonal();
302  plane_v2 = StartPos.Cross(plane_v1);
303  } else {
304  plane_v1.SetXYZ(1., 0., 0.);
305  plane_v2.SetXYZ(0., 1., 0.);
306  }
307 
308  if (StartMom.Mag2() == 0) {
309  std::cout
310  << "*** PndRecoKalmanFit2::Fit\tMomentum seed is ZERO. Cannot fit. ***"
311  << std::endl;
312  return tBefore;
313  }
314 
315  PndTrackCand trackCand = tBefore->GetTrackCand();
316 
317  genfit::AbsTrackRep *rep = new genfit::RKTrackRep(PDGCode);
318  // PndTrackCand does not store seed, then PndTrackCand2Genfit2TrackCand cannot convert the seed.
319  // You need to set the seed afterwards, taking it from PndTrack (setCovSeed/setPosMomSeedAndPdgCode)
320  genfit::TrackCand* gfCand = PndTrackCand2Genfit2TrackCand(&trackCand); // TODO: link TrackCand to track
321  gfCand->setCovSeed(covSeed);
322  gfCand->setPosMomSeedAndPdgCode(StartPos, StartMom, PDGCode);
323 
324  genfit::Track* trk = new genfit::Track(*gfCand, *fTheRecoHitFactory, rep);
325  delete (gfCand);
326 
327  // Start Fitter
328  try {
330  } catch (genfit::Exception& e) {
331 // delete(trk);
332  std::cout << "*** PndRecoKalmanFit2::Fit" << "\t"
333  << "FITTER EXCEPTION ***" << std::endl;
334  std::cerr << e.what();
335  }
336  if (fVerbose > 0) {
337  std::cout << "*** PndRecoKalmanFit2::Fit" << "\t" << "SUCCESSFULL FIT!"
338  << std::endl;
339  if (fVerbose > 2)
340  trk->getFitStatus()->Print();
341  }
342 
343  try {
344  tAfter = (PndTrack*) Genfit2Track2PndTrack(trk);
345  } catch (genfit::Exception& e) {
346  std::cout << "*** PndRecoKalmanFit2::Fit" << "\t"
347  << "CONVERSION EXCEPTION ***" << std::endl;
348  std::cerr << e.what();
349  tAfter = tBefore;
350  tAfter->SetFlag(-2); // flag -2: conversion failed
351 // delete(trk);
352  }
353 // tAfter = tBefore;
354  delete(trk);
355 
356  if (fVerbose > 0)
357  std::cout << "*** PndRecoKalmanFit2::Fit" << "\t" << "Fitting done"
358  << std::endl;
359 
360  return tAfter;
361 }
362 
static void quiet(bool b=true)
&quot;std::cerr &lt;&lt; e.what();&quot; will not write anything.
Definition: Exception.h:78
Track candidate – seed values and indices.
Definition: TrackCand.h:69
genfit::MeasurementFactory< genfit::AbsMeasurement > * fTheRecoHitFactory
void addProducer(int detID, AbsMeasurementProducer< measurement_T > *hitProd)
Register a producer module to the factory.
int fVerbose
Definition: poormantracks.C:24
FairGeanePro * fPro
static FieldManager * getInstance()
Get singleton instance.
Definition: FieldManager.h:112
TString fCentralTrackerBranchName
Name of the TCA for MVD.
Int_t fNumIt
Flag to use as initial plane the one perpendicular to the track.
Int_t fVerbose
Number of iterations.
Magnetic field.
Int_t GetGeometryType()
Definition: PndGeoFtsPar.h:37
Collection of TrackPoint objects, AbsTrackRep objects and FitStatus objects.
Definition: Track.h:71
TString fMvdBranchName
Geane Propagator.
Abstract base class for a track representation.
Definition: AbsTrackRep.h:66
Double_t par[3]
virtual void setMaxIterations(unsigned int n)
Set the maximum number of iterations.
Float_t fPropagateDistance
Flag to propagate to the interaction point.
void processTrack(Track *, bool resortHits=false)
Int_t GetGeometryType()
Definition: PndGeoSttPar.h:28
virtual void Print(const Option_t *="") const
genfit::TrackCand * PndTrackCand2Genfit2TrackCand(PndTrackCand *cand)
Bool_t fPerpPlane
Distance in [cm] to back-propagate the parameters, negative number means no backpropagation.
void init(AbsBField *b)
set the magnetic field here. Magnetic field classes must be derived from AbsBField.
Definition: FieldManager.h:78
AbsMaterialInterface implementation for use with ROOT&#39;s TGeoManager.
FairMCTracks * Track
Definition: drawEveTracks.C:8
static MaterialEffects * getInstance()
PndTrackCand GetTrackCand()
Definition: PndTrack.h:47
FitStatus * getFitStatus(const AbsTrackRep *rep=NULL) const
Get FitStatus for a AbsTrackRep. Per default, return FitStatus for cardinalRep.
Definition: Track.h:150
PndTrack * Fit(PndTrack *tBefore, Int_t PDG)
FairRuntimeDb * rtdb
Definition: hit_dirc.C:66
TFile * f
Definition: bump_analys.C:12
TClonesArray * FillTubeArray()
void init(AbsMaterialInterface *matIfc)
set the material interface here. Material interface classes must be derived from AbsMaterialInterface...
friend F32vec4 fabs(const F32vec4 &a)
Definition: P4_F32vec4.h:47
static PndGeoHandling * Instance()
GFTrack * trk
Definition: checkgenfit.C:13
Bool_t fUseGeane
Name of the TCA for central tracker.
Exception class for error handling in GENFIT (provides storage for diagnostic information) ...
Definition: Exception.h:48
void SetFlag(Int_t i)
Definition: PndTrack.h:38
void setPosMomSeedAndPdgCode(const TVector3 &pos, const TVector3 &mom, const int pdgCode)
This function works the same as setPosMomSeed but instead of a charge hypothesis you can set a pdg co...
ClassImp(PndAnaContFact)
TClonesArray * FillTubeArray()
this function will be used in PndFtsHitProducesRealFast
genfit::KalmanFitter fGenFitter
PndTrack * Genfit2Track2PndTrack(const genfit::Track *tr)
Bool_t fPropagateToIP
Flag to use Geane.
AbsTrackRep with 5D track parameterization in plane coordinates: (q/p, u&#39;, v&#39;, u, v) ...
Template class for a measurement producer module.
virtual const char * what() const
Standard error message handling for exceptions. use like &quot;std::cerr &lt;&lt; e.what();&quot;.
void setCovSeed(const TMatrixDSym &cov6D)
set the covariance matrix seed (6D).
Definition: TrackCand.h:175
FairTrackParP GetParamFirst()
Definition: PndTrack.h:49