Public Member Functions | Private Member Functions | Private Attributes | List of all members
trkf::Track3DKalman Class Reference
Inheritance diagram for trkf::Track3DKalman:
art::EDProducer art::detail::Producer art::detail::LegacyModule art::Modifier art::ModuleBase art::ProductRegistryHelper

Public Member Functions

 Track3DKalman (fhicl::ParameterSet const &pset)
 
- Public Member Functions inherited from art::EDProducer
 EDProducer (fhicl::ParameterSet const &pset)
 
template<typename Config >
 EDProducer (Table< Config > const &config)
 
std::string workerType () const
 
- Public Member Functions inherited from art::detail::Producer
virtual ~Producer () noexcept
 
 Producer (fhicl::ParameterSet const &)
 
 Producer (Producer const &)=delete
 
 Producer (Producer &&)=delete
 
Produceroperator= (Producer const &)=delete
 
Produceroperator= (Producer &&)=delete
 
void doBeginJob (SharedResources const &resources)
 
void doEndJob ()
 
void doRespondToOpenInputFile (FileBlock const &fb)
 
void doRespondToCloseInputFile (FileBlock const &fb)
 
void doRespondToOpenOutputFiles (FileBlock const &fb)
 
void doRespondToCloseOutputFiles (FileBlock const &fb)
 
bool doBeginRun (RunPrincipal &rp, ModuleContext const &mc)
 
bool doEndRun (RunPrincipal &rp, ModuleContext const &mc)
 
bool doBeginSubRun (SubRunPrincipal &srp, ModuleContext const &mc)
 
bool doEndSubRun (SubRunPrincipal &srp, ModuleContext const &mc)
 
bool doEvent (EventPrincipal &ep, ModuleContext const &mc, std::atomic< std::size_t > &counts_run, std::atomic< std::size_t > &counts_passed, std::atomic< std::size_t > &counts_failed)
 
- Public Member Functions inherited from art::Modifier
 ~Modifier () noexcept
 
 Modifier ()
 
 Modifier (Modifier const &)=delete
 
 Modifier (Modifier &&)=delete
 
Modifieroperator= (Modifier const &)=delete
 
Modifieroperator= (Modifier &&)=delete
 
- Public Member Functions inherited from art::ModuleBase
virtual ~ModuleBase () noexcept
 
 ModuleBase ()
 
ModuleDescription const & moduleDescription () const
 
void setModuleDescription (ModuleDescription const &)
 
std::array< std::vector< ProductInfo >, NumBranchTypes > const & getConsumables () const
 
void sortConsumables (std::string const &current_process_name)
 
template<typename T , BranchType BT>
ViewToken< T > consumesView (InputTag const &tag)
 
template<typename T , BranchType BT>
ViewToken< T > mayConsumeView (InputTag const &tag)
 

Private Member Functions

void produce (art::Event &evt) override
 
void beginJob () override
 
void endJob () override
 

Private Attributes

std::string fSpacePtsModuleLabel
 
std::string fGenieGenModuleLabel
 
std::string fG4ModuleLabel
 
bool fGenfPRINT
 
TTree * tree
 
TMatrixT< Double_t > * stMCT
 
TMatrixT< Double_t > * covMCT
 
TMatrixT< Double_t > * stREC
 
TMatrixT< Double_t > * covREC
 
Float_t chi2
 
Float_t chi2ndf
 
Float_t * fpRECt3D
 
Float_t * fpRECL
 
Float_t * fpREC
 
Float_t * fpMCT
 
int nfail
 
int ndf
 
unsigned int evtt
 
unsigned int nTrks
 
unsigned int fptsNo
 
Float_t * fshx
 
Float_t * fshy
 
Float_t * fshz
 
unsigned int fDimSize
 
std::vector< double > fPosErr
 
std::vector< double > fMomErr
 
std::vector< double > fMomStart
 
genf::GFAbsTrackReprepMC
 
genf::GFAbsTrackReprep
 
CLHEP::HepRandomEngine & fEngine
 

Additional Inherited Members

- Public Types inherited from art::EDProducer
using ModuleType = EDProducer
 
using WorkerType = WorkerT< EDProducer >
 
- Public Types inherited from art::detail::Producer
template<typename UserConfig , typename KeysToIgnore = void>
using Table = Modifier::Table< UserConfig, KeysToIgnore >
 
- Public Types inherited from art::Modifier
template<typename UserConfig , typename UserKeysToIgnore = void>
using Table = ProducerTable< UserConfig, detail::ModuleConfig, UserKeysToIgnore >
 
- Static Public Member Functions inherited from art::EDProducer
static void commitEvent (EventPrincipal &ep, Event &e)
 
- Protected Member Functions inherited from art::ModuleBase
ConsumesCollectorconsumesCollector ()
 
template<typename T , BranchType = InEvent>
ProductToken< T > consumes (InputTag const &)
 
template<typename Element , BranchType = InEvent>
ViewToken< Element > consumesView (InputTag const &)
 
template<typename T , BranchType = InEvent>
void consumesMany ()
 
template<typename T , BranchType = InEvent>
ProductToken< T > mayConsume (InputTag const &)
 
template<typename Element , BranchType = InEvent>
ViewToken< Element > mayConsumeView (InputTag const &)
 
template<typename T , BranchType = InEvent>
void mayConsumeMany ()
 

Detailed Description

Definition at line 62 of file Track3DKalman_module.cc.

Constructor & Destructor Documentation

trkf::Track3DKalman::Track3DKalman ( fhicl::ParameterSet const &  pset)
explicit

Definition at line 121 of file Track3DKalman_module.cc.

122  : EDProducer{pset}
123  , fSpacePtsModuleLabel{pset.get< std::string >("SpacePtsModuleLabel")}
124  , fGenieGenModuleLabel{pset.get< std::string >("GenieGenModuleLabel")}
125  , fG4ModuleLabel{pset.get< std::string >("G4ModuleLabel")}
126  , fGenfPRINT{pset.get< bool >("GenfPRINT")}
127  , fPosErr{pset.get< std::vector < double > >("PosErr3")} // resolution. cm
128  , fMomErr{pset.get< std::vector < double > >("MomErr3")} // GeV
129  , fMomStart{pset.get< std::vector < double > >("MomStart3")} // Will be unit norm'd.
130  // create a default random engine; obtain the random seed from NuRandomService,
131  // unless overridden in configuration with key "Seed"
132  , fEngine(art::ServiceHandle<rndm::NuRandomService>()->createEngine(*this, pset, "Seed"))
133 {
134  produces< std::vector<recob::Track> >();
135  produces< std::vector<recob::SpacePoint> >();
136  produces< art::Assns<recob::Track, recob::Cluster> >();
137  produces< art::Assns<recob::Track, recob::SpacePoint> >();
138  produces< art::Assns<recob::Track, recob::Hit> >();
139 }
std::string string
Definition: nybbler.cc:12
EDProducer(fhicl::ParameterSet const &pset)
Definition: EDProducer.h:20
std::vector< double > fMomStart
std::vector< double > fMomErr
CLHEP::HepRandomEngine & fEngine
std::vector< double > fPosErr

Member Function Documentation

void trkf::Track3DKalman::beginJob ( )
overrideprivatevirtual

Reimplemented from art::EDProducer.

Definition at line 142 of file Track3DKalman_module.cc.

143 {
145 
146  stMCT = new TMatrixT<Double_t>(5,1);
147  covMCT = new TMatrixT<Double_t>(5,5);
148  stREC = new TMatrixT<Double_t>(5,1);
149  covREC = new TMatrixT<Double_t>(5,5);
150 
151  fpMCT = new Float_t[4];
152  fpREC = new Float_t[4];
153  fpRECL = new Float_t[4];
154  fpRECt3D = new Float_t[4];
155  fDimSize = 400; // if necessary will get this from pset in constructor.
156 
157  fshx = new Float_t[fDimSize];
158  fshy = new Float_t[fDimSize];
159  fshz = new Float_t[fDimSize];
160 
161  tree = tfs->make<TTree>("GENFITttree","GENFITttree");
162  //tree->Branch("stMCT",&stMCT,"stMCT[5]/F"); // "TMatrixT<Double_t>"
163 
164  tree->Branch("stMCT","TMatrixD",&stMCT,64000,0);
165  //tree->Branch("covMCT",&covMCT,"covMCT[25]/F");
166  tree->Branch("covMCT","TMatrixD",&covMCT,64000,0);
167  //tree->Branch("stREC",&stREC,"stREC[5]/F");
168  tree->Branch("stREC","TMatrixD",&stREC,64000,0);
169  //tree->Branch("covREC",&covREC,"covREC[25]/F");
170  tree->Branch("covREC","TMatrixD",&covREC,64000,0);
171 
172 
173  tree->Branch("chi2",&chi2,"chi2/F");
174  tree->Branch("nfail",&nfail,"nfail/I");
175  tree->Branch("ndf",&ndf,"ndf/I");
176  tree->Branch("evtNo",&evtt,"evtNo/I");
177  tree->Branch("chi2ndf",&chi2ndf,"chi2ndf/F");
178 
179  tree->Branch("trkNo",&nTrks,"trkNo/I");
180  tree->Branch("ptsNo",&fptsNo,"ptsNo/I");
181  tree->Branch("shx",fshx,"shx[ptsNo]/F");
182  tree->Branch("shy",fshy,"shy[ptsNo]/F");
183  tree->Branch("shz",fshz,"shz[ptsNo]/F");
184 
185  tree->Branch("pMCT",fpMCT,"pMCT[4]/F");
186  tree->Branch("pRECKalF",fpREC,"pRECKalF[4]/F");
187  tree->Branch("pRECKalL",fpRECL,"pRECKalL[4]/F");
188  tree->Branch("pRECt3D",fpRECt3D,"pRECt3D[4]/F");
189 }
TMatrixT< Double_t > * stMCT
TMatrixT< Double_t > * stREC
TMatrixT< Double_t > * covMCT
TMatrixT< Double_t > * covREC
void trkf::Track3DKalman::endJob ( )
overrideprivatevirtual

Reimplemented from art::EDProducer.

Definition at line 192 of file Track3DKalman_module.cc.

193 {
194  if (!rep) delete rep;
195  if (!repMC) delete repMC;
196 }
genf::GFAbsTrackRep * rep
genf::GFAbsTrackRep * repMC
void trkf::Track3DKalman::produce ( art::Event evt)
overrideprivatevirtual
Todo:
: remove this statement, there is no place for checks on isRealData in reconstruction code

Implements art::EDProducer.

Definition at line 200 of file Track3DKalman_module.cc.

201 {
202  rep=0;
203  repMC=0;
204 
205  // get services
207  CLHEP::RandGaussQ gauss(fEngine);
208 
209  //////////////////////////////////////////////////////
210  // Make a std::unique_ptr<> for the thing you want to put into the event
211  // because that handles the memory management for you
212  //////////////////////////////////////////////////////
213  std::unique_ptr<std::vector<recob::Track> > tcol(new std::vector<recob::Track>);
214  std::unique_ptr< std::vector<recob::SpacePoint> > spcol (new std::vector<recob::SpacePoint>);
215  std::unique_ptr< art::Assns<recob::Track, recob::SpacePoint> > tspassn(new art::Assns<recob::Track, recob::SpacePoint>);
216  std::unique_ptr< art::Assns<recob::Track, recob::Cluster> > tcassn (new art::Assns<recob::Track, recob::Cluster>);
217  std::unique_ptr< art::Assns<recob::Track, recob::Hit> > thassn (new art::Assns<recob::Track, recob::Hit>);
218 
219  // define TPC parameters
220  TString tpcName = geom->GetLArTPCVolumeName();
221 
222 
223  // get input Hit object(s).
224  art::Handle< std::vector<recob::Track> > trackListHandle;
225  evt.getByLabel(fSpacePtsModuleLabel,trackListHandle);
226 
228 
229  if (!evt.isRealData())
230  {
231  art::Handle< std::vector<simb::MCTruth> > mctruthListHandle;
232  evt.getByLabel(fGenieGenModuleLabel,mctruthListHandle);
233 
234  for (unsigned int ii = 0; ii < mctruthListHandle->size(); ++ii)
235  {
236  art::Ptr<simb::MCTruth> mctparticle(mctruthListHandle,ii);
237  mclist.push_back(mctparticle);
238  }
239  }
240 
241  //create collection of spacepoints that will be used when creating the Track object
242  std::vector< art::Ptr<recob::SpacePoint> > spacepoints;
244  // std::cout<<"Run "<<evt.run()<<" Event "<<evt.id().event()<<std::endl;
245  mf::LogInfo("Track3DKalman: ") << "There are " << trackListHandle->size() << " Track3Dreco/SpacePt tracks/groups (whichever) in this event.";
246 
247  art::FindManyP<recob::SpacePoint> fmsp(trackListHandle, evt, fSpacePtsModuleLabel);
248  art::FindManyP<recob::Cluster> fmc (trackListHandle, evt, fSpacePtsModuleLabel);
249  art::FindManyP<recob::Hit> fmh (trackListHandle, evt, fSpacePtsModuleLabel);
250 
251  for(unsigned int ii = 0; ii < trackListHandle->size(); ++ii)
252  {
253  art::Ptr<recob::Track> track(trackListHandle, ii);
254  trackIn.push_back(track);
255 
256  }
257 
258  TVector3 MCOrigin;
259  TVector3 MCMomentum;
260  // TVector3 posErr(.05,.05,.05); // resolution. 0.5mm
261  // TVector3 momErr(.1,.1,0.2); // GeV
262  TVector3 posErr(fPosErr[0],fPosErr[1],fPosErr[2]); // resolution. 0.5mm
263  TVector3 momErr(fMomErr[0],fMomErr[1],fMomErr[2]); // GeV
264 
265  // This is strictly for MC
266  ///\todo: remove this statement, there is no place for checks on isRealData in reconstruction code
267  if (!evt.isRealData())
268  {
269  // Below breaks are stupid, I realize. But rather than keep all the MC
270  // particles I just take the first primary, e.g., muon and only keep its
271  // info in the branches of the Ttree. I could generalize later, ...
272  for( unsigned int ii = 0; ii < mclist.size(); ++ii )
273  {
274  //art::Ptr<const simb::MCTruth> mc(mctruthListHandle,i);
275  art::Ptr<simb::MCTruth> mc(mclist[ii]);
276  for(int jj = 0; jj < mc->NParticles(); ++jj)
277  {
278  simb::MCParticle part(mc->GetParticle(jj));
279  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 " ;
280  MCOrigin.SetXYZ(part.Vx(),part.Vy(),part.Vz()); // V for Vertex
281  MCMomentum.SetXYZ(part.Px(),part.Py(),part.Pz());
282  MCOrigin.Print();
283  MCMomentum.Print();
284  repMC = new genf::RKTrackRep(MCOrigin,
285  MCMomentum,
286  posErr,
287  momErr,
288  part.PdgCode());
289  break;
290  }
291  break;
292  }
293  //for saving of MC truth
294  stMCT->ResizeTo(repMC->getState());
295  *stMCT = repMC->getState();
296  covMCT-> ResizeTo(repMC->getCov());
297  *covMCT = repMC->getCov();
298  mf::LogInfo("Track3DKalman: ") <<" repMC, covMC are ... " ;
299  repMC->getState().Print();
300  repMC->getCov().Print();
301 
302  } // !isRealData
303 
305 
306  nTrks=0;
307  while(trackIter!=trackIn.end()) {
308  spacepoints.clear();
309  spacepoints = fmsp.at(nTrks);
310 
311  mf::LogInfo("Track3DKalman: ") << "found "<< spacepoints.size() <<" 3D spacepoint(s).";
312 
313  // Add the 3D track to the vector of the reconstructed tracks
314  if(spacepoints.size()>0){
315  // Insert the GENFIT/Kalman stuff here then make the tracks. Units are cm, GeV.
316  const double resolution = 0.5; // dunno, 5 mm
317  const int numIT = 3; // 3->1, EC, 6-Jan-2011. Back, 7-Jan-2011.
318 
319 
320  //TVector3 mom(0.0,0.0,2.0);
321  TVector3 mom(fMomStart[0],fMomStart[1],fMomStart[2]);
322  //mom.SetMag(1.);
323  TVector3 momM(mom);
324  momM.SetX(gauss.fire(momM.X(),momErr.X()/* *momM.X() */));
325  momM.SetY(gauss.fire(momM.Y(),momErr.Y()/* *momM.Y() */));
326  momM.SetZ(gauss.fire(momM.Z(),momErr.Z()/* *momM.Z() */));
327  //std::cout << "Track3DKalman: sort spacepoints by z
328 
329  std::sort(spacepoints.begin(), spacepoints.end(), sp_sort_3dz);
330 
331  //std::sort(spacepoints.begin(), spacepoints.end(), sp_sort_3dx); // Reverse sort!
332 
334  genf::GFDetPlane planeG((TVector3)(spacepoints[0]->XYZ()),momM);
335 
336 
337  // std::cout<<"Track3DKalman about to do GAbsTrackRep."<<std::endl;
338  // Initialize with 1st spacepoint location and a guess at the momentum.
339  rep = new genf::RKTrackRep(//posM-.5/momM.Mag()*momM,
340  (TVector3)(spacepoints[0]->XYZ()),
341  momM,
342  posErr,
343  momErr,
344  13); // mu- hypothesis
345  // std::cout<<"Track3DKalman: about to do GFTrack. repDim is " << rep->getDim() <<std::endl;
346 
347 
348  genf::GFTrack fitTrack(rep);//initialized with smeared rep
349  // Gonna sort in x cuz I want to essentially transform here to volTPC coords.
350  // volTPC coords, cuz that's what the Geant3/Geane stepper wants, as that's its understanding
351  // from the Geant4 geometry, which it'll use. EC, 7-Jan-2011.
352  int ihit = 0;
353  fptsNo = 0;
354  for (unsigned int point=0;point<spacepoints.size();++point)
355  {
356 
357  TVector3 spt3(spacepoints[point]->XYZ());
358  if (point > 0 )
359  {
360  double eps (0.1);
361  TVector3 magNew(spt3[0],spt3[1],spt3[2]);
362  TVector3 magLast(spacepoints.back()->XYZ()[0],
363  spacepoints.back()->XYZ()[1],
364  spacepoints.back()->XYZ()[2]
365  );
366  if (!(magNew.Mag()>=magLast.Mag()+eps ||
367  magNew.Mag()<=magLast.Mag()-eps)
368  )
369  continue;
370  }
371 
372  if (point%20) // Jump out of loop except on every 20th pt.
373  {
374  //continue;
375  // Icarus paper suggests we may wanna decimate our data in order to give
376  // trackfitter a better idea of multiple-scattering. EC, 7-Jan-2011.
377  //if (std::abs(spt3[0]-spacepoints.at(point-1).XYZ()[0]) < 2) continue;
378  }
379  if (fptsNo<fDimSize)
380  {
381  fshx[fptsNo] = spt3[0];
382  fshy[fptsNo] = spt3[1];
383  fshz[fptsNo] = spt3[2];
384  }
385  fptsNo++;
386 
387 
388  MF_LOG_DEBUG("Track3DKalman: ") << "ihit xyz..." << spt3[0]<<","<< spt3[1]<<","<< spt3[2];
389  fitTrack.addHit(new genf::PointHit(spt3,resolution),
390  1,//dummy detector id
391  ihit++
392  );
393  }
394 
395  // std::cout<<"Track3DKalman about to do GFKalman."<<std::endl;
397  //k.setBlowUpFactor(500); // Instead of 500 out of box. EC, 6-Jan-2011.
398  //k.setInitialDirection(+1); // Instead of 1 out of box. EC, 6-Jan-2011.
399  k.setNumIterations(numIT);
400  // std::cout<<"Track3DKalman back from setNumIterations."<<std::endl;
401  try{
402  // std::cout<<"Track3DKalman about to processTrack."<<std::endl;
403  k.processTrack(&fitTrack);
404  //std::cout<<"Track3DKalman back from processTrack."<<std::endl;
405  }
406  catch(GFException& e){
407  mf::LogError("Track3DKalman: ") << "just caught a GFException."<<std::endl;
408  e.what();
409  mf::LogError("Track3DKalman: ") << "Exceptions won't be further handled ->exit(1) "<<__LINE__;
410 
411  // exit(1);
412  }
413 
414  if(rep->getStatusFlag()==0) // 0 is successful completion
415  {
416  MF_LOG_DEBUG("Track3DKalman: ") << __FILE__ << " " << __LINE__ ;
417  MF_LOG_DEBUG("Track3DKalman: ") << "Track3DKalman.cxx: Original plane:";
418 
419  if(fGenfPRINT) planeG.Print();
420  MF_LOG_DEBUG("Track3DKalman: ") << "Current (fit) reference Plane:";
422 
423  MF_LOG_DEBUG("Track3DKalman: ") << "Track3DKalman.cxx: Last reference Plane:";
424  if(fGenfPRINT) rep->getLastPlane().Print();
425 
426  if(fGenfPRINT)
427  {
428  if(planeG!=rep->getReferencePlane())
429  MF_LOG_DEBUG("Track3DKalman: ") <<"Track3DKalman: Original hit plane (not surprisingly) not current reference Plane!"<<std::endl;
430  }
431 
432  stREC->ResizeTo(rep->getState());
433  *stREC = rep->getState();
434  covREC->ResizeTo(rep->getCov());
435  *covREC = rep->getCov();
436  if(fGenfPRINT)
437  {
438  MF_LOG_DEBUG("Track3DKalman: ") << " Final State and Cov:";
439  stREC->Print();
440  covREC->Print();
441  }
442  chi2 = rep->getChiSqu();
443  ndf = rep->getNDF();
444  nfail = fitTrack.getFailedHits();
445  chi2ndf = chi2/ndf;
446  TVector3 dircoss = (*trackIter)->VertexDirection<TVector3>();
447 
448  for (int ii=0;ii<3;++ii)
449  {
450  fpMCT[ii] = MCMomentum[ii]/MCMomentum.Mag();
451  fpREC[ii] = rep->getReferencePlane().getNormal()[ii];
452  fpRECL[ii] = rep->getLastPlane().getNormal()[ii];
453  fpRECt3D[ii] = dircoss[ii];
454  }
455  fpMCT[3] = MCMomentum.Mag();
456  fpREC[3] = -1.0/(*stREC)[0][0];
457 
458  evtt = (unsigned int) evt.id().event();
459  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] << ".";
460 
461  tree->Fill();
462 
463  // Get the clusters associated to each track in induction and collection view
464  std::vector< art::Ptr<recob::Cluster> > clusters = fmc.at(nTrks);
465  std::vector< art::Ptr<recob::Hit> > hits = fmh.at(nTrks);
466 
467  // Use new Track constructor... EC, 21-Apr-2011.
468  std::vector<TVector3> dircos(spacepoints.size());
469  dircos[0] = TVector3(fpREC);
470  dircos.back() = TVector3(fpRECL);
471 
472  // fill a vector of TVector3 with the space points used to make this track
473  std::vector<TVector3> xyz(spacepoints.size());
474  size_t spStart = spcol->size();
475  for(size_t tv = 0; tv < spacepoints.size(); ++tv){
476  xyz[tv] = TVector3(spacepoints[tv]->XYZ());
477  spcol->push_back(*(spacepoints[tv].get()));
478  }
479  size_t spEnd = spcol->size();
480 
483  recob::Track::Flags_t(xyz.size()), false),
484  0, -1., 0, recob::tracking::SMatrixSym55(), recob::tracking::SMatrixSym55(), tcol->size()-1));
485 
486  // associate the track with its clusters and tracks
487  util::CreateAssn(*this, evt, *tcol, clusters, *tcassn);
488  util::CreateAssn(*this, evt, *tcol, hits, *thassn);
489 
490  // associate the track to the space points
491  util::CreateAssn(*this, evt, *tcol, *spcol, *tspassn, spStart, spEnd);
492 
493  } // getStatusFlag
494  } // spacepoints.size()>0
495 
496  //
497  //std::cout<<"Track3DKalman found "<< tcol->size() <<" 3D track(s)"<<std::endl;
498  if(trackIter!=trackIn.end()) trackIter++;
499 
500  nTrks++;
501 
502  } // end loop over Track3Dreco/SpacePt tracks/groups (whichever) we brought into this event.
503 
504  evt.put(std::move(tcol));
505  evt.put(std::move(spcol));
506  evt.put(std::move(tcassn));
507  evt.put(std::move(thassn));
508  evt.put(std::move(tspassn));
509 
510 } // end method
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:217
std::vector< double > fMomStart
virtual const char * what() const
standard error message handling for exceptions. use like "std::cerr << e.what();" ...
Definition: GFException.cxx:48
TrackTrajectory::Flags_t Flags_t
Definition: Track.h:68
void Print(std::ostream &out=std::cout) const
Definition: GFDetPlane.cxx:247
const TMatrixT< Double_t > & getState() const
MaybeLogger_< ELseverityLevel::ELsev_error, false > LogError
std::vector< double > fMomErr
genf::GFAbsTrackRep * rep
TVector3 getNormal() const
Definition: GFDetPlane.cxx:145
bool isRealData() const
bool getByLabel(std::string const &label, std::string const &instance, Handle< PROD > &result) const
Definition: DataViewImpl.h:633
typename data_t::const_iterator const_iterator
Definition: PtrVector.h:55
const double e
void push_back(Ptr< U > const &p)
Definition: PtrVector.h:435
A trajectory in space reconstructed from hits.
def move(depos, offset)
Definition: depos.py:107
iterator end()
Definition: PtrVector.h:231
void setNumIterations(Int_t i)
Set number of iterations for Kalman Filter.
Definition: GFKalman.h:89
void processTrack(GFTrack *)
Performs fit on a GFTrack.
Definition: GFKalman.cxx:46
static GFFieldManager * getInstance()
ProductID put(std::unique_ptr< PROD > &&edp, std::string const &instance={})
Definition: DataViewImpl.h:686
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.
std::vector< Point_t > convertCollToPoint(std::vector< Point > const &coll)
Definition: TrackingTypes.h:55
const TMatrixT< Double_t > & getCov() const
size_type size() const
Definition: PtrVector.h:302
Exception class for error handling in GENFIT (provides storage for diagnostic information) ...
Definition: GFException.h:48
TMatrixT< Double_t > * stMCT
#define MF_LOG_DEBUG(id)
TMatrixT< Double_t > * stREC
void init(GFAbsBField *b)
set the magntic field here. Magnetic field classes must be derived from GFAbsBField ...
EventNumber_t event() const
Definition: EventID.h:116
CLHEP::HepRandomEngine & fEngine
double getChiSqu() const
std::vector< double > fPosErr
EventID id() const
Definition: Event.cc:34
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:49
QTextStream & endl(QTextStream &s)
TMatrixT< Double_t > * covMCT
std::string GetLArTPCVolumeName(geo::TPCID const &tpcid) const
Return the name of specified LAr TPC volume.
TMatrixT< Double_t > * covREC

Member Data Documentation

Float_t trkf::Track3DKalman::chi2
private

Definition at line 82 of file Track3DKalman_module.cc.

Float_t trkf::Track3DKalman::chi2ndf
private

Definition at line 83 of file Track3DKalman_module.cc.

TMatrixT<Double_t>* trkf::Track3DKalman::covMCT
private

Definition at line 79 of file Track3DKalman_module.cc.

TMatrixT<Double_t>* trkf::Track3DKalman::covREC
private

Definition at line 81 of file Track3DKalman_module.cc.

unsigned int trkf::Track3DKalman::evtt
private

Definition at line 91 of file Track3DKalman_module.cc.

unsigned int trkf::Track3DKalman::fDimSize
private

Definition at line 97 of file Track3DKalman_module.cc.

CLHEP::HepRandomEngine& trkf::Track3DKalman::fEngine
private

Definition at line 104 of file Track3DKalman_module.cc.

std::string trkf::Track3DKalman::fG4ModuleLabel
private

Definition at line 73 of file Track3DKalman_module.cc.

bool trkf::Track3DKalman::fGenfPRINT
private

Definition at line 74 of file Track3DKalman_module.cc.

std::string trkf::Track3DKalman::fGenieGenModuleLabel
private

Definition at line 72 of file Track3DKalman_module.cc.

std::vector<double> trkf::Track3DKalman::fMomErr
private

Definition at line 100 of file Track3DKalman_module.cc.

std::vector<double> trkf::Track3DKalman::fMomStart
private

Definition at line 101 of file Track3DKalman_module.cc.

Float_t* trkf::Track3DKalman::fpMCT
private

Definition at line 88 of file Track3DKalman_module.cc.

std::vector<double> trkf::Track3DKalman::fPosErr
private

Definition at line 99 of file Track3DKalman_module.cc.

Float_t* trkf::Track3DKalman::fpREC
private

Definition at line 87 of file Track3DKalman_module.cc.

Float_t* trkf::Track3DKalman::fpRECL
private

Definition at line 86 of file Track3DKalman_module.cc.

Float_t* trkf::Track3DKalman::fpRECt3D
private

Definition at line 85 of file Track3DKalman_module.cc.

unsigned int trkf::Track3DKalman::fptsNo
private

Definition at line 93 of file Track3DKalman_module.cc.

Float_t* trkf::Track3DKalman::fshx
private

Definition at line 94 of file Track3DKalman_module.cc.

Float_t* trkf::Track3DKalman::fshy
private

Definition at line 95 of file Track3DKalman_module.cc.

Float_t* trkf::Track3DKalman::fshz
private

Definition at line 96 of file Track3DKalman_module.cc.

std::string trkf::Track3DKalman::fSpacePtsModuleLabel
private

Definition at line 71 of file Track3DKalman_module.cc.

int trkf::Track3DKalman::ndf
private

Definition at line 90 of file Track3DKalman_module.cc.

int trkf::Track3DKalman::nfail
private

Definition at line 89 of file Track3DKalman_module.cc.

unsigned int trkf::Track3DKalman::nTrks
private

Definition at line 92 of file Track3DKalman_module.cc.

genf::GFAbsTrackRep* trkf::Track3DKalman::rep
private

Definition at line 103 of file Track3DKalman_module.cc.

genf::GFAbsTrackRep* trkf::Track3DKalman::repMC
private

Definition at line 102 of file Track3DKalman_module.cc.

TMatrixT<Double_t>* trkf::Track3DKalman::stMCT
private

Definition at line 78 of file Track3DKalman_module.cc.

TMatrixT<Double_t>* trkf::Track3DKalman::stREC
private

Definition at line 80 of file Track3DKalman_module.cc.

TTree* trkf::Track3DKalman::tree
private

Definition at line 76 of file Track3DKalman_module.cc.


The documentation for this class was generated from the following file: