LArSoft  v07_13_02
Liquid Argon Software toolkit - http://larsoft.org/
Track3DKalman_module.cc
Go to the documentation of this file.
1 //
3 // \file Track3DKalman.cxx
4 //
5 // \author echurch@fnal.gov
6 //
7 // This algorithm is designed to reconstruct 3D tracks through
8 // GENFIT's Kalman filter.
10 
11 // C++ includes
12 #include <math.h>
13 #include <algorithm>
14 #include <iostream>
15 #include <fstream>
16 
18 
19 
20 // Framework includes
23 #include "fhiclcpp/ParameterSet.h"
31 
32 // art extensions
34 
35 // LArSoft includes
36 
45 
46 
47 // ROOT includes
48 #include "TVectorD.h"
49 #include "TFile.h"
50 #include "TGeoManager.h"
51 #include "TF1.h"
52 #include "TGraph.h"
53 #include "TMath.h"
54 
55 // GENFIT includes
62 #include "larreco/Genfit/GFTrack.h"
64 #include "larreco/Genfit/GFDaf.h"
65 
67 #include <TTree.h>
68 #include <TMatrixT.h>
69 
71 
72 #include "CLHEP/Random/RandFlat.h"
73 #include "CLHEP/Random/RandGaussQ.h"
74 
75 #include <vector>
76 #include <string>
77 
78 //#include "RecoBase/SpacePoint.h"
79 
80 namespace trkf {
81 
82  class Track3DKalman : public art::EDProducer {
83 
84  public:
85 
86  explicit Track3DKalman(fhicl::ParameterSet const& pset);
88 
90  void produce(art::Event& evt);
91  void beginJob();
92  void endJob();
93  void reconfigure(fhicl::ParameterSet const& p);
94 
95  private:
96 
97  std::string fSpacePtsModuleLabel;// label for input collection
98  std::string fGenieGenModuleLabel;// label for input MC single particle generator
99  std::string fG4ModuleLabel;// label for input MC single particle generator
101 
102  // TFile *fileGENFIT;
103  TTree *tree;
104 
105  TMatrixT<Double_t> *stMCT;
106  TMatrixT<Double_t> *covMCT;
107  TMatrixT<Double_t> *stREC;
108  TMatrixT<Double_t> *covREC;
109  Float_t chi2;
110  Float_t chi2ndf;
111 
112  Float_t *fpRECt3D;
113  Float_t *fpRECL;
114  Float_t *fpREC;
115  Float_t *fpMCT;
116  int nfail;
117  int ndf;
118  unsigned int evtt;
119  unsigned int nTrks;
120  unsigned int fptsNo;
121  Float_t *fshx;
122  Float_t *fshy;
123  Float_t *fshz;
124  unsigned int fDimSize; // if necessary will get this from pset in constructor.
125 
126  std::vector<double> fPosErr;
127  std::vector<double> fMomErr;
128  std::vector<double> fMomStart;
131 
132 
133  }; // class Track3DKalman
134 
135 } // end namespace
136 
138 {
139  const double* xyz1 = h1->XYZ();
140  const double* xyz2 = h2->XYZ();
141  return xyz1[2] < xyz2[2];
142 }
143 
144 namespace trkf {
145 
146 //-------------------------------------------------
148 {
149 
150  this->reconfigure(pset);
151 
152  // create a default random engine; obtain the random seed from NuRandomService,
153  // unless overridden in configuration with key "Seed"
155 
156  produces< std::vector<recob::Track> >();
157  produces< std::vector<recob::SpacePoint> >();
158  produces< art::Assns<recob::Track, recob::Cluster> >();
159  produces< art::Assns<recob::Track, recob::SpacePoint> >();
160  produces< art::Assns<recob::Track, recob::Hit> >();
161 }
162 
163 //-------------------------------------------------
165 {
166 
167  fSpacePtsModuleLabel = pset.get< std::string >("SpacePtsModuleLabel");
168  fGenieGenModuleLabel = pset.get< std::string >("GenieGenModuleLabel");
169  fG4ModuleLabel = pset.get< std::string >("G4ModuleLabel");
170 
171  fPosErr = pset.get< std::vector < double > >("PosErr3"); // resolution. cm
172  fMomErr = pset.get< std::vector < double > >("MomErr3"); // GeV
173  fMomStart = pset.get< std::vector < double > >("MomStart3"); // Will be unit norm'd.
174  fGenfPRINT = pset.get< bool >("GenfPRINT");
175 
176 }
177 
178 //-------------------------------------------------
180 {
181 
182  /*
183  delete stMCT ;
184  delete covMCT;
185  delete stREC;
186  delete covREC;
187 
188  delete fpMCT;
189  delete fpREC;
190  delete fpRECL;
191  delete fpRECt3D;
192  */
193 }
194 
195 //-------------------------------------------------
197 {
198 
199 
201 
202  stMCT = new TMatrixT<Double_t>(5,1);
203  covMCT = new TMatrixT<Double_t>(5,5);
204  stREC = new TMatrixT<Double_t>(5,1);
205  covREC = new TMatrixT<Double_t>(5,5);
206 
207  fpMCT = new Float_t[4];
208  fpREC = new Float_t[4];
209  fpRECL = new Float_t[4];
210  fpRECt3D = new Float_t[4];
211  fDimSize = 400; // if necessary will get this from pset in constructor.
212 
213  fshx = new Float_t[fDimSize];
214  fshy = new Float_t[fDimSize];
215  fshz = new Float_t[fDimSize];
216 
217 // //TFile fileGENFIT("GENFITout.root","RECREATE");
218 
219 
220  tree = tfs->make<TTree>("GENFITttree","GENFITttree");
221  //tree->Branch("stMCT",&stMCT,"stMCT[5]/F"); // "TMatrixT<Double_t>"
222 
223  tree->Branch("stMCT","TMatrixD",&stMCT,64000,0);
224  //tree->Branch("covMCT",&covMCT,"covMCT[25]/F");
225  tree->Branch("covMCT","TMatrixD",&covMCT,64000,0);
226  //tree->Branch("stREC",&stREC,"stREC[5]/F");
227  tree->Branch("stREC","TMatrixD",&stREC,64000,0);
228  //tree->Branch("covREC",&covREC,"covREC[25]/F");
229  tree->Branch("covREC","TMatrixD",&covREC,64000,0);
230 
231 
232  tree->Branch("chi2",&chi2,"chi2/F");
233  tree->Branch("nfail",&nfail,"nfail/I");
234  tree->Branch("ndf",&ndf,"ndf/I");
235  tree->Branch("evtNo",&evtt,"evtNo/I");
236  tree->Branch("chi2ndf",&chi2ndf,"chi2ndf/F");
237 
238  tree->Branch("trkNo",&nTrks,"trkNo/I");
239  tree->Branch("ptsNo",&fptsNo,"ptsNo/I");
240  tree->Branch("shx",fshx,"shx[ptsNo]/F");
241  tree->Branch("shy",fshy,"shy[ptsNo]/F");
242  tree->Branch("shz",fshz,"shz[ptsNo]/F");
243 
244  tree->Branch("pMCT",fpMCT,"pMCT[4]/F");
245  tree->Branch("pRECKalF",fpREC,"pRECKalF[4]/F");
246  tree->Branch("pRECKalL",fpRECL,"pRECKalL[4]/F");
247  tree->Branch("pRECt3D",fpRECt3D,"pRECt3D[4]/F");
248 
249 
250  //TGeoManager* geomGENFIT = new TGeoManager("Geometry", "Geane geometry");
251  //TGeoManager::Import("config/genfitGeom.root");
252  // gROOT->Macro("config/Geane.C");
253 
254 }
255 
256 //-------------------------------------------------
258 {
259  if (!rep) delete rep;
260  if (!repMC) delete repMC;
261 }
262 
263 
264 //------------------------------------------------------------------------------------//
266 {
267 
268  rep=0;
269  repMC=0;
270 
271  // get services
273  // get the random number generator service and make some CLHEP generators
275  CLHEP::HepRandomEngine &engine = rng->getEngine();
276  CLHEP::RandGaussQ gauss(engine);
277 
279  // Make a std::unique_ptr<> for the thing you want to put into the event
280  // because that handles the memory management for you
282  std::unique_ptr<std::vector<recob::Track> > tcol(new std::vector<recob::Track>);
283  std::unique_ptr< std::vector<recob::SpacePoint> > spcol (new std::vector<recob::SpacePoint>);
284  std::unique_ptr< art::Assns<recob::Track, recob::SpacePoint> > tspassn(new art::Assns<recob::Track, recob::SpacePoint>);
285  std::unique_ptr< art::Assns<recob::Track, recob::Cluster> > tcassn (new art::Assns<recob::Track, recob::Cluster>);
286  std::unique_ptr< art::Assns<recob::Track, recob::Hit> > thassn (new art::Assns<recob::Track, recob::Hit>);
287 
288  // define TPC parameters
289  TString tpcName = geom->GetLArTPCVolumeName();
290 
291 
292  // get input Hit object(s).
293  art::Handle< std::vector<recob::Track> > trackListHandle;
294  evt.getByLabel(fSpacePtsModuleLabel,trackListHandle);
295 
297 
298  if (!evt.isRealData())
299  {
300 
301  // std::cout << "Track3DKalman: This is MC." << std::endl;
302  // std::cout<<"Run "<<evt.run()<<" Event "<<evt.id().event()<<std::endl;
303 
304  art::Handle< std::vector<simb::MCTruth> > mctruthListHandle;
305  evt.getByLabel(fGenieGenModuleLabel,mctruthListHandle);
306 
307  for (unsigned int ii = 0; ii < mctruthListHandle->size(); ++ii)
308  {
309  art::Ptr<simb::MCTruth> mctparticle(mctruthListHandle,ii);
310  mclist.push_back(mctparticle);
311  }
312 
313  }
314 
315  //create collection of spacepoints that will be used when creating the Track object
316  std::vector< art::Ptr<recob::SpacePoint> > spacepoints;
318  // std::cout<<"Run "<<evt.run()<<" Event "<<evt.id().event()<<std::endl;
319  mf::LogInfo("Track3DKalman: ") << "There are " << trackListHandle->size() << " Track3Dreco/SpacePt tracks/groups (whichever) in this event.";
320 
321  art::FindManyP<recob::SpacePoint> fmsp(trackListHandle, evt, fSpacePtsModuleLabel);
322  art::FindManyP<recob::Cluster> fmc (trackListHandle, evt, fSpacePtsModuleLabel);
323  art::FindManyP<recob::Hit> fmh (trackListHandle, evt, fSpacePtsModuleLabel);
324 
325  for(unsigned int ii = 0; ii < trackListHandle->size(); ++ii)
326  {
327  art::Ptr<recob::Track> track(trackListHandle, ii);
328  trackIn.push_back(track);
329 
330  }
331 
332  TVector3 MCOrigin;
333  TVector3 MCMomentum;
334  // TVector3 posErr(.05,.05,.05); // resolution. 0.5mm
335  // TVector3 momErr(.1,.1,0.2); // GeV
336  TVector3 posErr(fPosErr[0],fPosErr[1],fPosErr[2]); // resolution. 0.5mm
337  TVector3 momErr(fMomErr[0],fMomErr[1],fMomErr[2]); // GeV
338 
339  // This is strictly for MC
341  if (!evt.isRealData())
342  {
343  // Below breaks are stupid, I realize. But rather than keep all the MC
344  // particles I just take the first primary, e.g., muon and only keep its
345  // info in the branches of the Ttree. I could generalize later, ...
346  for( unsigned int ii = 0; ii < mclist.size(); ++ii )
347  {
348  //art::Ptr<const simb::MCTruth> mc(mctruthListHandle,i);
349  art::Ptr<simb::MCTruth> mc(mclist[ii]);
350  for(int jj = 0; jj < mc->NParticles(); ++jj)
351  {
353  mf::LogInfo("Track3DKalman: ") << "FROM MC TRUTH, the particle's pdg code is: "<<part.PdgCode()<< " with energy = "<<part.E() <<", with energy = "<<part.E()<< " and vtx and momentum in Global (not volTPC) coords are " ;
354  MCOrigin.SetXYZ(part.Vx(),part.Vy(),part.Vz()); // V for Vertex
355  MCMomentum.SetXYZ(part.Px(),part.Py(),part.Pz());
356  MCOrigin.Print();
357  MCMomentum.Print();
358  repMC = new genf::RKTrackRep(MCOrigin,
359  MCMomentum,
360  posErr,
361  momErr,
362  part.PdgCode());
363  break;
364  }
365  break;
366  }
367  //for saving of MC truth
368  stMCT->ResizeTo(repMC->getState());
369  *stMCT = repMC->getState();
370  covMCT-> ResizeTo(repMC->getCov());
371  *covMCT = repMC->getCov();
372  mf::LogInfo("Track3DKalman: ") <<" repMC, covMC are ... " ;
373  repMC->getState().Print();
374  repMC->getCov().Print();
375 
376  } // !isRealData
377 
379 
380  nTrks=0;
381  while(trackIter!=trackIn.end()) {
382  spacepoints.clear();
383  spacepoints = fmsp.at(nTrks);
384 
385  mf::LogInfo("Track3DKalman: ") << "found "<< spacepoints.size() <<" 3D spacepoint(s).";
386 
387  // Add the 3D track to the vector of the reconstructed tracks
388  if(spacepoints.size()>0){
389  // Insert the GENFIT/Kalman stuff here then make the tracks. Units are cm, GeV.
390  const double resolution = 0.5; // dunno, 5 mm
391  const int numIT = 3; // 3->1, EC, 6-Jan-2011. Back, 7-Jan-2011.
392 
393 
394  //TVector3 mom(0.0,0.0,2.0);
395  TVector3 mom(fMomStart[0],fMomStart[1],fMomStart[2]);
396  //mom.SetMag(1.);
397  TVector3 momM(mom);
398  momM.SetX(gauss.fire(momM.X(),momErr.X()/* *momM.X() */));
399  momM.SetY(gauss.fire(momM.Y(),momErr.Y()/* *momM.Y() */));
400  momM.SetZ(gauss.fire(momM.Z(),momErr.Z()/* *momM.Z() */));
401  //std::cout << "Track3DKalman: sort spacepoints by z
402 
403  std::sort(spacepoints.begin(), spacepoints.end(), sp_sort_3dz);
404 
405  //std::sort(spacepoints.begin(), spacepoints.end(), sp_sort_3dx); // Reverse sort!
406 
408  genf::GFDetPlane planeG((TVector3)(spacepoints[0]->XYZ()),momM);
409 
410 
411  // std::cout<<"Track3DKalman about to do GAbsTrackRep."<<std::endl;
412  // Initialize with 1st spacepoint location and a guess at the momentum.
413  rep = new genf::RKTrackRep(//posM-.5/momM.Mag()*momM,
414  (TVector3)(spacepoints[0]->XYZ()),
415  momM,
416  posErr,
417  momErr,
418  13); // mu- hypothesis
419  // std::cout<<"Track3DKalman: about to do GFTrack. repDim is " << rep->getDim() <<std::endl;
420 
421 
422  genf::GFTrack fitTrack(rep);//initialized with smeared rep
423  // Gonna sort in x cuz I want to essentially transform here to volTPC coords.
424  // volTPC coords, cuz that's what the Geant3/Geane stepper wants, as that's its understanding
425  // from the Geant4 geometry, which it'll use. EC, 7-Jan-2011.
426  int ihit = 0;
427  fptsNo = 0;
428  for (unsigned int point=0;point<spacepoints.size();++point)
429  {
430 
431  TVector3 spt3(spacepoints[point]->XYZ());
432  if (point > 0 )
433  {
434  double eps (0.1);
435  TVector3 magNew(spt3[0],spt3[1],spt3[2]);
436  TVector3 magLast(spacepoints.back()->XYZ()[0],
437  spacepoints.back()->XYZ()[1],
438  spacepoints.back()->XYZ()[2]
439  );
440  if (!(magNew.Mag()>=magLast.Mag()+eps ||
441  magNew.Mag()<=magLast.Mag()-eps)
442  )
443  continue;
444  }
445 
446  if (point%20) // Jump out of loop except on every 20th pt.
447  {
448  //continue;
449  // Icarus paper suggests we may wanna decimate our data in order to give
450  // trackfitter a better idea of multiple-scattering. EC, 7-Jan-2011.
451  //if (std::abs(spt3[0]-spacepoints.at(point-1).XYZ()[0]) < 2) continue;
452  }
453  if (fptsNo<fDimSize)
454  {
455  fshx[fptsNo] = spt3[0];
456  fshy[fptsNo] = spt3[1];
457  fshz[fptsNo] = spt3[2];
458  }
459  fptsNo++;
460 
461 
462  LOG_DEBUG("Track3DKalman: ") << "ihit xyz..." << spt3[0]<<","<< spt3[1]<<","<< spt3[2];
463  fitTrack.addHit(new genf::PointHit(spt3,resolution),
464  1,//dummy detector id
465  ihit++
466  );
467  }
468 
469  // std::cout<<"Track3DKalman about to do GFKalman."<<std::endl;
470  genf::GFKalman k;
471  //k.setBlowUpFactor(500); // Instead of 500 out of box. EC, 6-Jan-2011.
472  //k.setInitialDirection(+1); // Instead of 1 out of box. EC, 6-Jan-2011.
473  k.setNumIterations(numIT);
474  // std::cout<<"Track3DKalman back from setNumIterations."<<std::endl;
475  try{
476  // std::cout<<"Track3DKalman about to processTrack."<<std::endl;
477  k.processTrack(&fitTrack);
478  //std::cout<<"Track3DKalman back from processTrack."<<std::endl;
479  }
480  catch(GFException& e){
481  mf::LogError("Track3DKalman: ") << "just caught a GFException."<<std::endl;
482  e.what();
483  mf::LogError("Track3DKalman: ") << "Exceptions won't be further handled ->exit(1) "<<__LINE__;
484 
485  // exit(1);
486  }
487 
488  if(rep->getStatusFlag()==0) // 0 is successful completion
489  {
490  LOG_DEBUG("Track3DKalman: ") << __FILE__ << " " << __LINE__ ;
491  LOG_DEBUG("Track3DKalman: ") << "Track3DKalman.cxx: Original plane:";
492 
493  if(fGenfPRINT) planeG.Print();
494  LOG_DEBUG("Track3DKalman: ") << "Current (fit) reference Plane:";
496 
497  LOG_DEBUG("Track3DKalman: ") << "Track3DKalman.cxx: Last reference Plane:";
498  if(fGenfPRINT) rep->getLastPlane().Print();
499 
500  if(fGenfPRINT)
501  {
502  if(planeG!=rep->getReferencePlane())
503  LOG_DEBUG("Track3DKalman: ") <<"Track3DKalman: Original hit plane (not surprisingly) not current reference Plane!"<<std::endl;
504  }
505 
506  stREC->ResizeTo(rep->getState());
507  *stREC = rep->getState();
508  covREC->ResizeTo(rep->getCov());
509  *covREC = rep->getCov();
510  if(fGenfPRINT)
511  {
512  LOG_DEBUG("Track3DKalman: ") << " Final State and Cov:";
513  stREC->Print();
514  covREC->Print();
515  }
516  chi2 = rep->getChiSqu();
517  ndf = rep->getNDF();
518  nfail = fitTrack.getFailedHits();
519  chi2ndf = chi2/ndf;
520  TVector3 dircoss = (*trackIter)->VertexDirection<TVector3>();
521 
522  for (int ii=0;ii<3;++ii)
523  {
524  fpMCT[ii] = MCMomentum[ii]/MCMomentum.Mag();
525  fpREC[ii] = rep->getReferencePlane().getNormal()[ii];
526  fpRECL[ii] = rep->getLastPlane().getNormal()[ii];
527  fpRECt3D[ii] = dircoss[ii];
528  }
529  fpMCT[3] = MCMomentum.Mag();
530  fpREC[3] = -1.0/(*stREC)[0][0];
531 
532  evtt = (unsigned int) evt.id().event();
533  mf::LogInfo("Track3DKalman: ") << "Track3DKalman about to do tree->Fill(). Chi2/ndf is " << chi2/ndf << ". All in volTPC coords .... pMCT[0-3] is " << fpMCT[0] << ", " << fpMCT[1] << ", " << fpMCT[2] << ", " << fpMCT[3] << ". pREC[0-3] is " << fpREC[0] << ", "<< fpREC[1] << ", " << fpREC[2] << ", " << fpREC[3] << ".";
534 
535  tree->Fill();
536 
537  // Get the clusters associated to each track in induction and collection view
538  std::vector< art::Ptr<recob::Cluster> > clusters = fmc.at(nTrks);
539  std::vector< art::Ptr<recob::Hit> > hits = fmh.at(nTrks);
540 
541  // Use new Track constructor... EC, 21-Apr-2011.
542  std::vector<TVector3> dircos(spacepoints.size());
543  dircos[0] = TVector3(fpREC);
544  dircos.back() = TVector3(fpRECL);
545 
546  // fill a vector of TVector3 with the space points used to make this track
547  std::vector<TVector3> xyz(spacepoints.size());
548  size_t spStart = spcol->size();
549  for(size_t tv = 0; tv < spacepoints.size(); ++tv){
550  xyz[tv] = TVector3(spacepoints[tv]->XYZ());
551  spcol->push_back(*(spacepoints[tv].get()));
552  }
553  size_t spEnd = spcol->size();
554 
557  recob::Track::Flags_t(xyz.size()), false),
558  0, -1., 0, recob::tracking::SMatrixSym55(), recob::tracking::SMatrixSym55(), tcol->size()-1));
559 
560  // associate the track with its clusters and tracks
561  util::CreateAssn(*this, evt, *tcol, clusters, *tcassn);
562  util::CreateAssn(*this, evt, *tcol, hits, *thassn);
563 
564  // associate the track to the space points
565  util::CreateAssn(*this, evt, *tcol, *spcol, *tspassn, spStart, spEnd);
566 
567  } // getStatusFlag
568  } // spacepoints.size()>0
569 
570  //
571  //std::cout<<"Track3DKalman found "<< tcol->size() <<" 3D track(s)"<<std::endl;
572  if(trackIter!=trackIn.end()) trackIter++;
573 
574  nTrks++;
575 
576  } // end loop over Track3Dreco/SpacePt tracks/groups (whichever) we brought into this event.
577 
578  evt.put(std::move(tcol));
579  evt.put(std::move(spcol));
580  evt.put(std::move(tcassn));
581  evt.put(std::move(thassn));
582  evt.put(std::move(tspassn));
583 
584 } // end method
585 
587 
588 
589 
590 } // end namespace trkf
unsigned int getNDF() const
static bool sp_sort_3dz(const art::Ptr< recob::SpacePoint > &h1, const art::Ptr< recob::SpacePoint > &h2)
const GFDetPlane & getReferencePlane() const
GFDetPlane getLastPlane() const
ROOT::Math::SMatrix< Double32_t, 5, 5, ROOT::Math::MatRepSym< Double32_t, 5 > > SMatrixSym55
Definition: TrackingTypes.h:85
MaybeLogger_< ELseverityLevel::ELsev_info, false > LogInfo
std::vector< Vector_t > convertCollToVector(std::vector< Vector > const &coll)
Definition: TrackingTypes.h:58
iterator begin()
Definition: PtrVector.h:223
Declaration of signal hit object.
std::vector< double > fMomStart
virtual const char * what() const
standard error message handling for exceptions. use like "std::cerr << e.what();" ...
Definition: GFException.cxx:46
TrackTrajectory::Flags_t Flags_t
Definition: Track.h:71
int NParticles() const
Definition: MCTruth.h:72
void Print(std::ostream &out=std::cout) const
Definition: GFDetPlane.cxx:242
const TMatrixT< Double_t > & getState() const
MaybeLogger_< ELseverityLevel::ELsev_error, false > LogError
bool isRealData() const
Definition: Event.h:83
void produce(art::Event &evt)
std::vector< double > fMomErr
genf::GFAbsTrackRep * rep
Base Class for genfit track representations. Defines interface for track parameterizations.
Definition: GFAbsTrackRep.h:85
ProductID put(std::unique_ptr< PROD > &&product)
Definition: Event.h:102
TVector3 getNormal() const
Definition: GFDetPlane.cxx:140
base_engine_t & getEngine() const
void hits()
Definition: readHits.C:15
base_engine_t & createEngine(seed_t seed)
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:42
Provides recob::Track data product.
void push_back(Ptr< U > const &p)
Definition: PtrVector.h:441
A trajectory in space reconstructed from hits.
T get(std::string const &key) const
Definition: ParameterSet.h:231
TString part[npart]
Definition: Style.C:32
iterator end()
Definition: PtrVector.h:237
void setNumIterations(Int_t i)
Set number of iterations for Kalman Filter.
Definition: GFKalman.h:90
int getFailedHits(int repId=-1)
return the number of failed Hits in track fit repId == -1 will use cardinal rep
Definition: GFTrack.h:148
void processTrack(GFTrack *)
Performs fit on a GFTrack.
Definition: GFKalman.cxx:48
static GFFieldManager * getInstance()
genf::GFAbsTrackRep * repMC
bool CreateAssn(PRODUCER const &prod, art::Event &evt, std::vector< T > const &a, art::Ptr< U > const &b, art::Assns< U, T > &assn, std::string a_instance, size_t indx=UINT_MAX)
Creates a single one-to-one association.
const Double32_t * XYZ() const
Definition: SpacePoint.h:65
std::vector< Point_t > convertCollToPoint(std::vector< Point > const &coll)
Definition: TrackingTypes.h:55
Declaration of cluster object.
const TMatrixT< Double_t > & getCov() const
Track3DKalman(fhicl::ParameterSet const &pset)
An art service to assist in the distribution of guaranteed unique seeds to all engines within an art ...
size_type size() const
Definition: PtrVector.h:308
TH1F * h2
Definition: plot.C:46
data_t::const_iterator const_iterator
Definition: PtrVector.h:61
const simb::MCParticle & GetParticle(int i) const
Definition: MCTruth.h:73
void reconfigure(fhicl::ParameterSet const &p)
T * make(ARGS...args) const
static art::ServiceHandle< art::RandomNumberGenerator > & rng()
Utility object to perform functions of association.
Exception class for error handling in GENFIT (provides storage for diagnostic information) ...
Definition: GFException.h:50
TMatrixT< Double_t > * stMCT
bool getByLabel(std::string const &label, std::string const &productInstanceName, Handle< PROD > &result) const
Definition: DataViewImpl.h:344
TMatrixT< Double_t > * stREC
void init(GFAbsBField *b)
set the magntic field here. Magnetic field classes must be derived from GFAbsBField ...
#define LOG_DEBUG(id)
TH1F * h1
Definition: plot.C:43
EventNumber_t event() const
Definition: EventID.h:117
double getChiSqu() const
void addHit(GFAbsRecoHit *theHit)
deprecated!
Definition: GFTrack.h:299
TCEvent evt
Definition: DataStructs.cxx:5
Float_t e
Definition: plot.C:34
Float_t track
Definition: plot.C:34
Tools and modules for checking out the basics of the Monte Carlo.
std::vector< double > fPosErr
EventID id() const
Definition: Event.h:56
art framework interface to geometry description
Track from a non-cascading particle.A recob::Track consists of a recob::TrackTrajectory, plus additional members relevant for a "fitted" track:
Definition: Track.h:52
TMatrixT< Double_t > * covMCT
std::string GetLArTPCVolumeName(geo::TPCID const &tpcid) const
Return the name of specified LAr TPC volume.
TMatrixT< Double_t > * covREC