LArSoft  v06_85_00
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
63 #include "larreco/Genfit/GFTrack.h"
65 #include "larreco/Genfit/GFDaf.h"
66 
68 #include <TTree.h>
69 #include <TMatrixT.h>
70 
72 
73 #include "CLHEP/Random/RandFlat.h"
74 #include "CLHEP/Random/RandGaussQ.h"
75 
76 #include <vector>
77 #include <string>
78 
79 //#include "RecoBase/SpacePoint.h"
80 
81 namespace trkf {
82 
83  class Track3DKalman : public art::EDProducer {
84 
85  public:
86 
87  explicit Track3DKalman(fhicl::ParameterSet const& pset);
89 
91  void produce(art::Event& evt);
92  void beginJob();
93  void endJob();
94  void reconfigure(fhicl::ParameterSet const& p);
95 
96  private:
97 
98  std::string fSpacePtsModuleLabel;// label for input collection
99  std::string fGenieGenModuleLabel;// label for input MC single particle generator
100  std::string fG4ModuleLabel;// label for input MC single particle generator
102 
103  // TFile *fileGENFIT;
104  TTree *tree;
105 
106  TMatrixT<Double_t> *stMCT;
107  TMatrixT<Double_t> *covMCT;
108  TMatrixT<Double_t> *stREC;
109  TMatrixT<Double_t> *covREC;
110  Float_t chi2;
111  Float_t chi2ndf;
112 
113  Float_t *fpRECt3D;
114  Float_t *fpRECL;
115  Float_t *fpREC;
116  Float_t *fpMCT;
117  int nfail;
118  int ndf;
119  unsigned int evtt;
120  unsigned int nTrks;
121  unsigned int fptsNo;
122  Float_t *fshx;
123  Float_t *fshy;
124  Float_t *fshz;
125  unsigned int fDimSize; // if necessary will get this from pset in constructor.
126 
127  std::vector<double> fPosErr;
128  std::vector<double> fMomErr;
129  std::vector<double> fMomStart;
132 
133 
134  }; // class Track3DKalman
135 
136 } // end namespace
137 
139 {
140  const double* xyz1 = h1->XYZ();
141  const double* xyz2 = h2->XYZ();
142  return xyz1[2] < xyz2[2];
143 }
144 
145 namespace trkf {
146 
147 //-------------------------------------------------
149 {
150 
151  this->reconfigure(pset);
152 
153  // create a default random engine; obtain the random seed from NuRandomService,
154  // unless overridden in configuration with key "Seed"
156 
157  produces< std::vector<recob::Track> >();
158  produces< std::vector<recob::SpacePoint> >();
159  produces< art::Assns<recob::Track, recob::Cluster> >();
160  produces< art::Assns<recob::Track, recob::SpacePoint> >();
161  produces< art::Assns<recob::Track, recob::Hit> >();
162 }
163 
164 //-------------------------------------------------
166 {
167 
168  fSpacePtsModuleLabel = pset.get< std::string >("SpacePtsModuleLabel");
169  fGenieGenModuleLabel = pset.get< std::string >("GenieGenModuleLabel");
170  fG4ModuleLabel = pset.get< std::string >("G4ModuleLabel");
171 
172  fPosErr = pset.get< std::vector < double > >("PosErr3"); // resolution. cm
173  fMomErr = pset.get< std::vector < double > >("MomErr3"); // GeV
174  fMomStart = pset.get< std::vector < double > >("MomStart3"); // Will be unit norm'd.
175  fGenfPRINT = pset.get< bool >("GenfPRINT");
176 
177 }
178 
179 //-------------------------------------------------
181 {
182 
183  /*
184  delete stMCT ;
185  delete covMCT;
186  delete stREC;
187  delete covREC;
188 
189  delete fpMCT;
190  delete fpREC;
191  delete fpRECL;
192  delete fpRECt3D;
193  */
194 }
195 
196 //-------------------------------------------------
198 {
199 
200 
202 
203  stMCT = new TMatrixT<Double_t>(5,1);
204  covMCT = new TMatrixT<Double_t>(5,5);
205  stREC = new TMatrixT<Double_t>(5,1);
206  covREC = new TMatrixT<Double_t>(5,5);
207 
208  fpMCT = new Float_t[4];
209  fpREC = new Float_t[4];
210  fpRECL = new Float_t[4];
211  fpRECt3D = new Float_t[4];
212  fDimSize = 400; // if necessary will get this from pset in constructor.
213 
214  fshx = new Float_t[fDimSize];
215  fshy = new Float_t[fDimSize];
216  fshz = new Float_t[fDimSize];
217 
218 // //TFile fileGENFIT("GENFITout.root","RECREATE");
219 
220 
221  tree = tfs->make<TTree>("GENFITttree","GENFITttree");
222  //tree->Branch("stMCT",&stMCT,"stMCT[5]/F"); // "TMatrixT<Double_t>"
223 
224  tree->Branch("stMCT","TMatrixD",&stMCT,64000,0);
225  //tree->Branch("covMCT",&covMCT,"covMCT[25]/F");
226  tree->Branch("covMCT","TMatrixD",&covMCT,64000,0);
227  //tree->Branch("stREC",&stREC,"stREC[5]/F");
228  tree->Branch("stREC","TMatrixD",&stREC,64000,0);
229  //tree->Branch("covREC",&covREC,"covREC[25]/F");
230  tree->Branch("covREC","TMatrixD",&covREC,64000,0);
231 
232 
233  tree->Branch("chi2",&chi2,"chi2/F");
234  tree->Branch("nfail",&nfail,"nfail/I");
235  tree->Branch("ndf",&ndf,"ndf/I");
236  tree->Branch("evtNo",&evtt,"evtNo/I");
237  tree->Branch("chi2ndf",&chi2ndf,"chi2ndf/F");
238 
239  tree->Branch("trkNo",&nTrks,"trkNo/I");
240  tree->Branch("ptsNo",&fptsNo,"ptsNo/I");
241  tree->Branch("shx",fshx,"shx[ptsNo]/F");
242  tree->Branch("shy",fshy,"shy[ptsNo]/F");
243  tree->Branch("shz",fshz,"shz[ptsNo]/F");
244 
245  tree->Branch("pMCT",fpMCT,"pMCT[4]/F");
246  tree->Branch("pRECKalF",fpREC,"pRECKalF[4]/F");
247  tree->Branch("pRECKalL",fpRECL,"pRECKalL[4]/F");
248  tree->Branch("pRECt3D",fpRECt3D,"pRECt3D[4]/F");
249 
250 
251  //TGeoManager* geomGENFIT = new TGeoManager("Geometry", "Geane geometry");
252  //TGeoManager::Import("config/genfitGeom.root");
253  // gROOT->Macro("config/Geane.C");
254 
255 }
256 
257 //-------------------------------------------------
259 {
260  if (!rep) delete rep;
261  if (!repMC) delete repMC;
262 }
263 
264 
265 //------------------------------------------------------------------------------------//
267 {
268 
269  rep=0;
270  repMC=0;
271 
272  // get services
274  // get the random number generator service and make some CLHEP generators
276  CLHEP::HepRandomEngine &engine = rng->getEngine();
277  CLHEP::RandGaussQ gauss(engine);
278 
280  // Make a std::unique_ptr<> for the thing you want to put into the event
281  // because that handles the memory management for you
283  std::unique_ptr<std::vector<recob::Track> > tcol(new std::vector<recob::Track>);
284  std::unique_ptr< std::vector<recob::SpacePoint> > spcol (new std::vector<recob::SpacePoint>);
285  std::unique_ptr< art::Assns<recob::Track, recob::SpacePoint> > tspassn(new art::Assns<recob::Track, recob::SpacePoint>);
286  std::unique_ptr< art::Assns<recob::Track, recob::Cluster> > tcassn (new art::Assns<recob::Track, recob::Cluster>);
287  std::unique_ptr< art::Assns<recob::Track, recob::Hit> > thassn (new art::Assns<recob::Track, recob::Hit>);
288 
289  // define TPC parameters
290  TString tpcName = geom->GetLArTPCVolumeName();
291 
292 
293  // get input Hit object(s).
294  art::Handle< std::vector<recob::Track> > trackListHandle;
295  evt.getByLabel(fSpacePtsModuleLabel,trackListHandle);
296 
298 
299  if (!evt.isRealData())
300  {
301 
302  // std::cout << "Track3DKalman: This is MC." << std::endl;
303  // std::cout<<"Run "<<evt.run()<<" Event "<<evt.id().event()<<std::endl;
304 
305  art::Handle< std::vector<simb::MCTruth> > mctruthListHandle;
306  evt.getByLabel(fGenieGenModuleLabel,mctruthListHandle);
307 
308  for (unsigned int ii = 0; ii < mctruthListHandle->size(); ++ii)
309  {
310  art::Ptr<simb::MCTruth> mctparticle(mctruthListHandle,ii);
311  mclist.push_back(mctparticle);
312  }
313 
314  }
315 
316  //create collection of spacepoints that will be used when creating the Track object
317  std::vector< art::Ptr<recob::SpacePoint> > spacepoints;
319  // std::cout<<"Run "<<evt.run()<<" Event "<<evt.id().event()<<std::endl;
320  mf::LogInfo("Track3DKalman: ") << "There are " << trackListHandle->size() << " Track3Dreco/SpacePt tracks/groups (whichever) in this event.";
321 
322  art::FindManyP<recob::SpacePoint> fmsp(trackListHandle, evt, fSpacePtsModuleLabel);
323  art::FindManyP<recob::Cluster> fmc (trackListHandle, evt, fSpacePtsModuleLabel);
324  art::FindManyP<recob::Hit> fmh (trackListHandle, evt, fSpacePtsModuleLabel);
325 
326  for(unsigned int ii = 0; ii < trackListHandle->size(); ++ii)
327  {
328  art::Ptr<recob::Track> track(trackListHandle, ii);
329  trackIn.push_back(track);
330 
331  }
332 
333  TVector3 MCOrigin;
334  TVector3 MCMomentum;
335  // TVector3 posErr(.05,.05,.05); // resolution. 0.5mm
336  // TVector3 momErr(.1,.1,0.2); // GeV
337  TVector3 posErr(fPosErr[0],fPosErr[1],fPosErr[2]); // resolution. 0.5mm
338  TVector3 momErr(fMomErr[0],fMomErr[1],fMomErr[2]); // GeV
339 
340  // This is strictly for MC
342  if (!evt.isRealData())
343  {
344  // Below breaks are stupid, I realize. But rather than keep all the MC
345  // particles I just take the first primary, e.g., muon and only keep its
346  // info in the branches of the Ttree. I could generalize later, ...
347  for( unsigned int ii = 0; ii < mclist.size(); ++ii )
348  {
349  //art::Ptr<const simb::MCTruth> mc(mctruthListHandle,i);
350  art::Ptr<simb::MCTruth> mc(mclist[ii]);
351  for(int jj = 0; jj < mc->NParticles(); ++jj)
352  {
354  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 " ;
355  MCOrigin.SetXYZ(part.Vx(),part.Vy(),part.Vz()); // V for Vertex
356  MCMomentum.SetXYZ(part.Px(),part.Py(),part.Pz());
357  MCOrigin.Print();
358  MCMomentum.Print();
359  repMC = new genf::RKTrackRep(MCOrigin,
360  MCMomentum,
361  posErr,
362  momErr,
363  part.PdgCode());
364  break;
365  }
366  break;
367  }
368  //for saving of MC truth
369  stMCT->ResizeTo(repMC->getState());
370  *stMCT = repMC->getState();
371  covMCT-> ResizeTo(repMC->getCov());
372  *covMCT = repMC->getCov();
373  mf::LogInfo("Track3DKalman: ") <<" repMC, covMC are ... " ;
374  repMC->getState().Print();
375  repMC->getCov().Print();
376 
377  } // !isRealData
378 
380 
381  nTrks=0;
382  while(trackIter!=trackIn.end()) {
383  spacepoints.clear();
384  spacepoints = fmsp.at(nTrks);
385 
386  mf::LogInfo("Track3DKalman: ") << "found "<< spacepoints.size() <<" 3D spacepoint(s).";
387 
388  // Add the 3D track to the vector of the reconstructed tracks
389  if(spacepoints.size()>0){
390  // Insert the GENFIT/Kalman stuff here then make the tracks. Units are cm, GeV.
391  const double resolution = 0.5; // dunno, 5 mm
392  const int numIT = 3; // 3->1, EC, 6-Jan-2011. Back, 7-Jan-2011.
393 
394 
395  //TVector3 mom(0.0,0.0,2.0);
396  TVector3 mom(fMomStart[0],fMomStart[1],fMomStart[2]);
397  //mom.SetMag(1.);
398  TVector3 momM(mom);
399  momM.SetX(gauss.fire(momM.X(),momErr.X()/* *momM.X() */));
400  momM.SetY(gauss.fire(momM.Y(),momErr.Y()/* *momM.Y() */));
401  momM.SetZ(gauss.fire(momM.Z(),momErr.Z()/* *momM.Z() */));
402  //std::cout << "Track3DKalman: sort spacepoints by z
403 
404  std::sort(spacepoints.begin(), spacepoints.end(), sp_sort_3dz);
405 
406  //std::sort(spacepoints.begin(), spacepoints.end(), sp_sort_3dx); // Reverse sort!
407 
409  genf::GFDetPlane planeG((TVector3)(spacepoints[0]->XYZ()),momM);
410 
411 
412  // std::cout<<"Track3DKalman about to do GAbsTrackRep."<<std::endl;
413  // Initialize with 1st spacepoint location and a guess at the momentum.
414  rep = new genf::RKTrackRep(//posM-.5/momM.Mag()*momM,
415  (TVector3)(spacepoints[0]->XYZ()),
416  momM,
417  posErr,
418  momErr,
419  13); // mu- hypothesis
420  // std::cout<<"Track3DKalman: about to do GFTrack. repDim is " << rep->getDim() <<std::endl;
421 
422 
423  genf::GFTrack fitTrack(rep);//initialized with smeared rep
424  // Gonna sort in x cuz I want to essentially transform here to volTPC coords.
425  // volTPC coords, cuz that's what the Geant3/Geane stepper wants, as that's its understanding
426  // from the Geant4 geometry, which it'll use. EC, 7-Jan-2011.
427  int ihit = 0;
428  fptsNo = 0;
429  for (unsigned int point=0;point<spacepoints.size();++point)
430  {
431 
432  TVector3 spt3(spacepoints[point]->XYZ());
433  if (point > 0 )
434  {
435  double eps (0.1);
436  TVector3 magNew(spt3[0],spt3[1],spt3[2]);
437  TVector3 magLast(spacepoints.back()->XYZ()[0],
438  spacepoints.back()->XYZ()[1],
439  spacepoints.back()->XYZ()[2]
440  );
441  if (!(magNew.Mag()>=magLast.Mag()+eps ||
442  magNew.Mag()<=magLast.Mag()-eps)
443  )
444  continue;
445  }
446 
447  if (point%20) // Jump out of loop except on every 20th pt.
448  {
449  //continue;
450  // Icarus paper suggests we may wanna decimate our data in order to give
451  // trackfitter a better idea of multiple-scattering. EC, 7-Jan-2011.
452  //if (std::abs(spt3[0]-spacepoints.at(point-1).XYZ()[0]) < 2) continue;
453  }
454  if (fptsNo<fDimSize)
455  {
456  fshx[fptsNo] = spt3[0];
457  fshy[fptsNo] = spt3[1];
458  fshz[fptsNo] = spt3[2];
459  }
460  fptsNo++;
461 
462 
463  LOG_DEBUG("Track3DKalman: ") << "ihit xyz..." << spt3[0]<<","<< spt3[1]<<","<< spt3[2];
464  fitTrack.addHit(new genf::PointHit(spt3,resolution),
465  1,//dummy detector id
466  ihit++
467  );
468  }
469 
470  // std::cout<<"Track3DKalman about to do GFKalman."<<std::endl;
471  genf::GFKalman k;
472  //k.setBlowUpFactor(500); // Instead of 500 out of box. EC, 6-Jan-2011.
473  //k.setInitialDirection(+1); // Instead of 1 out of box. EC, 6-Jan-2011.
474  k.setNumIterations(numIT);
475  // std::cout<<"Track3DKalman back from setNumIterations."<<std::endl;
476  try{
477  // std::cout<<"Track3DKalman about to processTrack."<<std::endl;
478  k.processTrack(&fitTrack);
479  //std::cout<<"Track3DKalman back from processTrack."<<std::endl;
480  }
481  catch(GFException& e){
482  mf::LogError("Track3DKalman: ") << "just caught a GFException."<<std::endl;
483  e.what();
484  mf::LogError("Track3DKalman: ") << "Exceptions won't be further handled ->exit(1) "<<__LINE__;
485 
486  // exit(1);
487  }
488 
489  if(rep->getStatusFlag()==0) // 0 is successful completion
490  {
491  LOG_DEBUG("Track3DKalman: ") << __FILE__ << " " << __LINE__ ;
492  LOG_DEBUG("Track3DKalman: ") << "Track3DKalman.cxx: Original plane:";
493 
494  if(fGenfPRINT) planeG.Print();
495  LOG_DEBUG("Track3DKalman: ") << "Current (fit) reference Plane:";
497 
498  LOG_DEBUG("Track3DKalman: ") << "Track3DKalman.cxx: Last reference Plane:";
499  if(fGenfPRINT) rep->getLastPlane().Print();
500 
501  if(fGenfPRINT)
502  {
503  if(planeG!=rep->getReferencePlane())
504  LOG_DEBUG("Track3DKalman: ") <<"Track3DKalman: Original hit plane (not surprisingly) not current reference Plane!"<<std::endl;
505  }
506 
507  stREC->ResizeTo(rep->getState());
508  *stREC = rep->getState();
509  covREC->ResizeTo(rep->getCov());
510  *covREC = rep->getCov();
511  if(fGenfPRINT)
512  {
513  LOG_DEBUG("Track3DKalman: ") << " Final State and Cov:";
514  stREC->Print();
515  covREC->Print();
516  }
517  chi2 = rep->getChiSqu();
518  ndf = rep->getNDF();
519  nfail = fitTrack.getFailedHits();
520  chi2ndf = chi2/ndf;
521  double dircoss[3],dircose[3];
522  (*trackIter)->Direction(dircoss,dircose);
523 
524  for (int ii=0;ii<3;++ii)
525  {
526  fpMCT[ii] = MCMomentum[ii]/MCMomentum.Mag();
527  fpREC[ii] = rep->getReferencePlane().getNormal()[ii];
528  fpRECL[ii] = rep->getLastPlane().getNormal()[ii];
529  fpRECt3D[ii] = dircoss[ii];
530  }
531  fpMCT[3] = MCMomentum.Mag();
532  fpREC[3] = -1.0/(*stREC)[0][0];
533 
534  evtt = (unsigned int) evt.id().event();
535  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] << ".";
536 
537  tree->Fill();
538 
539  // Get the clusters associated to each track in induction and collection view
540  std::vector< art::Ptr<recob::Cluster> > clusters = fmc.at(nTrks);
541  std::vector< art::Ptr<recob::Hit> > hits = fmh.at(nTrks);
542 
543  // Use new Track constructor... EC, 21-Apr-2011.
544  std::vector<TVector3> dircos(spacepoints.size());
545  dircos[0] = TVector3(fpREC);
546  dircos.back() = TVector3(fpRECL);
547 
548  // fill a vector of TVector3 with the space points used to make this track
549  std::vector<TVector3> xyz(spacepoints.size());
550  size_t spStart = spcol->size();
551  for(size_t tv = 0; tv < spacepoints.size(); ++tv){
552  xyz[tv] = TVector3(spacepoints[tv]->XYZ());
553  spcol->push_back(*(spacepoints[tv].get()));
554  }
555  size_t spEnd = spcol->size();
556 
557  tcol->push_back(recob::Track(xyz, dircos, std::vector< std::vector<double> >(0),
558  std::vector<double>(2, util::kBogusD), 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
MaybeLogger_< ELseverityLevel::ELsev_info, false > LogInfo
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
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
auto vector(Vector const &v)
Returns a manipulator which will print the specified array.
Definition: DumpUtils.h:265
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
void push_back(Ptr< U > const &p)
Definition: PtrVector.h:441
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.
Declaration of cluster object.
const TMatrixT< Double_t > & getCov() const
Track3DKalman(fhicl::ParameterSet const &pset)
Provides recob::Track data product.
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()
const double * XYZ() const
Definition: SpacePoint.h:64
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
constexpr double kBogusD
obviously bogus double value
double getChiSqu() const
void addHit(GFAbsRecoHit *theHit)
deprecated!
Definition: GFTrack.h:299
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:51
TMatrixT< Double_t > * covMCT
std::string GetLArTPCVolumeName(geo::TPCID const &tpcid) const
Return the name of specified LAr TPC volume.
TMatrixT< Double_t > * covREC