LArSoft  v07_13_02
Liquid Argon Software toolkit - http://larsoft.org/
trkf::Track3DKalman Class Reference
Inheritance diagram for trkf::Track3DKalman:
art::EDProducer art::ProducerBase art::Consumer art::EngineCreator art::ProductRegistryHelper

Public Types

using ModuleType = EDProducer
 
using WorkerType = WorkerT< EDProducer >
 
template<typename UserConfig , typename KeysToIgnore = void>
using Table = ProducerBase::Table< UserConfig, KeysToIgnore >
 

Public Member Functions

 Track3DKalman (fhicl::ParameterSet const &pset)
 
 ~Track3DKalman ()
 
void produce (art::Event &evt)
 
void beginJob ()
 
void endJob ()
 
void reconfigure (fhicl::ParameterSet const &p)
 
template<typename PROD , BranchType B = InEvent>
ProductID getProductID (std::string const &instanceName={}) const
 
template<typename PROD , BranchType B>
ProductID getProductID (ModuleDescription const &moduleDescription, std::string const &instanceName) const
 
bool modifiesEvent () const
 
template<typename T , BranchType = InEvent>
ProductToken< T > consumes (InputTag const &)
 
template<typename T , art::BranchType BT>
art::ProductToken< T > consumes (InputTag const &it)
 
template<typename T , BranchType = InEvent>
void consumesMany ()
 
template<typename Element , BranchType = InEvent>
ViewToken< Element > consumesView (InputTag const &)
 
template<typename T , art::BranchType BT>
art::ViewToken< T > consumesView (InputTag const &it)
 
template<typename T , BranchType = InEvent>
ProductToken< T > mayConsume (InputTag const &)
 
template<typename T , art::BranchType BT>
art::ProductToken< T > mayConsume (InputTag const &it)
 
template<typename T , BranchType = InEvent>
void mayConsumeMany ()
 
template<typename Element , BranchType = InEvent>
ViewToken< Element > mayConsumeView (InputTag const &)
 
template<typename T , art::BranchType BT>
art::ViewToken< T > mayConsumeView (InputTag const &it)
 
base_engine_tcreateEngine (seed_t seed)
 
base_engine_tcreateEngine (seed_t seed, std::string const &kind_of_engine_to_make)
 
base_engine_tcreateEngine (seed_t seed, std::string const &kind_of_engine_to_make, label_t const &engine_label)
 
seed_t get_seed_value (fhicl::ParameterSet const &pset, char const key[]="seed", seed_t const implicit_seed=-1)
 

Static Public Member Functions

static cet::exempt_ptr< Consumernon_module_context ()
 

Protected Member Functions

CurrentProcessingContext const * currentContext () const
 
void validateConsumedProduct (BranchType const bt, ProductInfo const &pi)
 
void prepareForJob (fhicl::ParameterSet const &pset)
 
void showMissingConsumes () const
 

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
 

Detailed Description

Definition at line 82 of file Track3DKalman_module.cc.

Member Typedef Documentation

using art::EDProducer::ModuleType = EDProducer
inherited

Definition at line 34 of file EDProducer.h.

template<typename UserConfig , typename KeysToIgnore = void>
using art::EDProducer::Table = ProducerBase::Table<UserConfig, KeysToIgnore>
inherited

Definition at line 43 of file EDProducer.h.

using art::EDProducer::WorkerType = WorkerT<EDProducer>
inherited

Definition at line 35 of file EDProducer.h.

Constructor & Destructor Documentation

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

Definition at line 147 of file Track3DKalman_module.cc.

References art::EngineCreator::createEngine(), and reconfigure().

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 }
base_engine_t & createEngine(seed_t seed)
void reconfigure(fhicl::ParameterSet const &p)
trkf::Track3DKalman::~Track3DKalman ( )

Definition at line 179 of file Track3DKalman_module.cc.

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 }

Member Function Documentation

void trkf::Track3DKalman::beginJob ( )
virtual

Reimplemented from art::EDProducer.

Definition at line 196 of file Track3DKalman_module.cc.

References chi2, chi2ndf, covMCT, covREC, evtt, fDimSize, fpMCT, fpREC, fpRECL, fpRECt3D, fptsNo, fshx, fshy, fshz, art::TFileDirectory::make(), ndf, nfail, nTrks, stMCT, stREC, and tree.

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 }
T * make(ARGS...args) const
TMatrixT< Double_t > * stMCT
TMatrixT< Double_t > * stREC
TMatrixT< Double_t > * covMCT
TMatrixT< Double_t > * covREC
template<typename T , BranchType = InEvent>
ProductToken<T> art::Consumer::consumes ( InputTag const &  )
inherited
template<typename T , art::BranchType BT>
art::ProductToken<T> art::Consumer::consumes ( InputTag const &  it)
inherited

Definition at line 147 of file Consumer.h.

References art::InputTag::instance(), art::InputTag::label(), and art::InputTag::process().

148 {
149  if (!moduleContext_)
150  return ProductToken<T>::invalid();
151 
152  consumables_[BT].emplace_back(ConsumableType::Product,
153  TypeID{typeid(T)},
154  it.label(),
155  it.instance(),
156  it.process());
157  return ProductToken<T>{it};
158 }
static ProductToken< T > invalid()
Definition: ProductToken.h:47
ConsumableProducts consumables_
Definition: Consumer.h:138
bool moduleContext_
Definition: Consumer.h:136
template<typename T , art::BranchType BT>
void art::Consumer::consumesMany ( )
inherited

Definition at line 162 of file Consumer.h.

163 {
164  if (!moduleContext_)
165  return;
166 
167  consumables_[BT].emplace_back(ConsumableType::Many, TypeID{typeid(T)});
168 }
ConsumableProducts consumables_
Definition: Consumer.h:138
bool moduleContext_
Definition: Consumer.h:136
template<typename Element , BranchType = InEvent>
ViewToken<Element> art::Consumer::consumesView ( InputTag const &  )
inherited
template<typename T , art::BranchType BT>
art::ViewToken<T> art::Consumer::consumesView ( InputTag const &  it)
inherited

Definition at line 172 of file Consumer.h.

References art::InputTag::instance(), art::InputTag::label(), and art::InputTag::process().

173 {
174  if (!moduleContext_)
175  return ViewToken<T>::invalid();
176 
177  consumables_[BT].emplace_back(ConsumableType::ViewElement,
178  TypeID{typeid(T)},
179  it.label(),
180  it.instance(),
181  it.process());
182  return ViewToken<T>{it};
183 }
static ViewToken< Element > invalid()
Definition: ProductToken.h:75
ConsumableProducts consumables_
Definition: Consumer.h:138
bool moduleContext_
Definition: Consumer.h:136
EngineCreator::base_engine_t & EngineCreator::createEngine ( seed_t  seed,
std::string const &  kind_of_engine_to_make 
)
inherited

Definition at line 32 of file EngineCreator.cc.

References art::EngineCreator::rng().

34 {
35  return rng()->createEngine(
36  placeholder_schedule_id(), seed, kind_of_engine_to_make);
37 }
long seed
Definition: chem4.cc:68
static art::ServiceHandle< art::RandomNumberGenerator > & rng()
EngineCreator::base_engine_t & EngineCreator::createEngine ( seed_t  seed,
std::string const &  kind_of_engine_to_make,
label_t const &  engine_label 
)
inherited

Definition at line 40 of file EngineCreator.cc.

References art::EngineCreator::rng().

43 {
44  return rng()->createEngine(
45  placeholder_schedule_id(), seed, kind_of_engine_to_make, engine_label);
46 }
long seed
Definition: chem4.cc:68
static art::ServiceHandle< art::RandomNumberGenerator > & rng()
CurrentProcessingContext const * art::EDProducer::currentContext ( ) const
protectedinherited

Definition at line 120 of file EDProducer.cc.

References art::EDProducer::current_context_.

121  {
122  return current_context_.get();
123  }
CPC_exempt_ptr current_context_
Definition: EDProducer.h:116
void trkf::Track3DKalman::endJob ( )
virtual

Reimplemented from art::EDProducer.

Definition at line 257 of file Track3DKalman_module.cc.

References rep, and repMC.

258 {
259  if (!rep) delete rep;
260  if (!repMC) delete repMC;
261 }
genf::GFAbsTrackRep * rep
genf::GFAbsTrackRep * repMC
EngineCreator::seed_t EngineCreator::get_seed_value ( fhicl::ParameterSet const &  pset,
char const  key[] = "seed",
seed_t const  implicit_seed = -1 
)
inherited

Definition at line 49 of file EngineCreator.cc.

References fhicl::ParameterSet::get().

Referenced by art::MixFilter< T >::initEngine_().

52 {
53  auto const& explicit_seeds = pset.get<std::vector<int>>(key, {});
54  return explicit_seeds.empty() ? implicit_seed : explicit_seeds.front();
55 }
template<typename PROD , BranchType B>
ProductID art::EDProducer::getProductID ( std::string const &  instanceName = {}) const
inlineinherited

Definition at line 123 of file EDProducer.h.

References art::EDProducer::moduleDescription_.

124  {
125  return ProducerBase::getProductID<PROD, B>(moduleDescription_,
126  instanceName);
127  }
ModuleDescription moduleDescription_
Definition: EDProducer.h:115
template<typename PROD , BranchType B>
ProductID art::ProducerBase::getProductID ( ModuleDescription const &  moduleDescription,
std::string const &  instanceName 
) const
inherited

Definition at line 56 of file ProducerBase.h.

References B, and art::ModuleDescription::moduleLabel().

Referenced by art::ProducerBase::modifiesEvent().

58  {
59  auto const& pd =
60  get_ProductDescription<PROD>(B, md.moduleLabel(), instanceName);
61  return pd.productID();
62  }
Int_t B
Definition: plot.C:25
template<typename T , BranchType = InEvent>
ProductToken<T> art::Consumer::mayConsume ( InputTag const &  )
inherited
template<typename T , art::BranchType BT>
art::ProductToken<T> art::Consumer::mayConsume ( InputTag const &  it)
inherited

Definition at line 190 of file Consumer.h.

References art::InputTag::instance(), art::InputTag::label(), and art::InputTag::process().

191 {
192  if (!moduleContext_)
193  return ProductToken<T>::invalid();
194 
195  consumables_[BT].emplace_back(ConsumableType::Product,
196  TypeID{typeid(T)},
197  it.label(),
198  it.instance(),
199  it.process());
200  return ProductToken<T>{it};
201 }
static ProductToken< T > invalid()
Definition: ProductToken.h:47
ConsumableProducts consumables_
Definition: Consumer.h:138
bool moduleContext_
Definition: Consumer.h:136
template<typename T , art::BranchType BT>
void art::Consumer::mayConsumeMany ( )
inherited

Definition at line 205 of file Consumer.h.

206 {
207  if (!moduleContext_)
208  return;
209 
210  consumables_[BT].emplace_back(ConsumableType::Many, TypeID{typeid(T)});
211 }
ConsumableProducts consumables_
Definition: Consumer.h:138
bool moduleContext_
Definition: Consumer.h:136
template<typename Element , BranchType = InEvent>
ViewToken<Element> art::Consumer::mayConsumeView ( InputTag const &  )
inherited
template<typename T , art::BranchType BT>
art::ViewToken<T> art::Consumer::mayConsumeView ( InputTag const &  it)
inherited

Definition at line 215 of file Consumer.h.

References art::InputTag::instance(), art::InputTag::label(), and art::InputTag::process().

216 {
217  if (!moduleContext_)
218  return ViewToken<T>::invalid();
219 
220  consumables_[BT].emplace_back(ConsumableType::ViewElement,
221  TypeID{typeid(T)},
222  it.label(),
223  it.instance(),
224  it.process());
225  return ViewToken<T>{it};
226 }
static ViewToken< Element > invalid()
Definition: ProductToken.h:75
ConsumableProducts consumables_
Definition: Consumer.h:138
bool moduleContext_
Definition: Consumer.h:136
bool art::ProducerBase::modifiesEvent ( ) const
inlineinherited

Definition at line 40 of file ProducerBase.h.

References art::ProducerBase::getProductID().

41  {
42  return true;
43  }
void art::Consumer::prepareForJob ( fhicl::ParameterSet const &  pset)
protectedinherited

Definition at line 89 of file Consumer.cc.

References fhicl::ParameterSet::get_if_present().

Referenced by art::EDProducer::doBeginJob(), art::EDFilter::doBeginJob(), and art::EDAnalyzer::doBeginJob().

90 {
91  if (!moduleContext_)
92  return;
93 
94  pset.get_if_present("errorOnMissingConsumes", requireConsumes_);
95  for (auto& consumablesPerBranch : consumables_) {
96  cet::sort_all(consumablesPerBranch);
97  }
98 }
bool requireConsumes_
Definition: Consumer.h:137
ConsumableProducts consumables_
Definition: Consumer.h:138
bool moduleContext_
Definition: Consumer.h:136
void trkf::Track3DKalman::produce ( art::Event evt)
virtual
Todo:
: remove this statement, there is no place for checks on isRealData in reconstruction code

Implements art::EDProducer.

Definition at line 265 of file Track3DKalman_module.cc.

References genf::GFTrack::addHit(), art::PtrVector< T >::begin(), chi2, chi2ndf, recob::tracking::convertCollToPoint(), recob::tracking::convertCollToVector(), covMCT, covREC, util::CreateAssn(), DEFINE_ART_MODULE, e, art::PtrVector< T >::end(), art::EventID::event(), evtt, fDimSize, fGenfPRINT, fGenieGenModuleLabel, fMomErr, fMomStart, fpMCT, fPosErr, fpREC, fpRECL, fpRECt3D, fptsNo, fshx, fshy, fshz, fSpacePtsModuleLabel, art::DataViewImpl::getByLabel(), genf::GFAbsTrackRep::getChiSqu(), genf::GFAbsTrackRep::getCov(), art::RandomNumberGenerator::getEngine(), genf::GFTrack::getFailedHits(), genf::GFFieldManager::getInstance(), geo::GeometryCore::GetLArTPCVolumeName(), genf::GFAbsTrackRep::getLastPlane(), genf::GFAbsTrackRep::getNDF(), genf::GFDetPlane::getNormal(), simb::MCTruth::GetParticle(), genf::GFAbsTrackRep::getReferencePlane(), genf::GFAbsTrackRep::getState(), genf::GFAbsTrackRep::getStatusFlag(), hits(), art::Event::id(), genf::GFFieldManager::init(), art::Event::isRealData(), LOG_DEBUG, ndf, nfail, simb::MCTruth::NParticles(), nTrks, part, genf::GFDetPlane::Print(), genf::GFKalman::processTrack(), art::PtrVector< T >::push_back(), art::Event::put(), rep, repMC, art::EngineCreator::rng(), genf::GFKalman::setNumIterations(), art::PtrVector< T >::size(), sp_sort_3dz(), stMCT, stREC, track, tree, and GFException::what().

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  {
352  simb::MCParticle part(mc->GetParticle(jj));
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
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
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
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
std::vector< double > fMomErr
genf::GFAbsTrackRep * rep
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
void push_back(Ptr< U > const &p)
Definition: PtrVector.h:441
A trajectory in space reconstructed from hits.
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
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.
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:308
data_t::const_iterator const_iterator
Definition: PtrVector.h:61
static art::ServiceHandle< art::RandomNumberGenerator > & rng()
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)
EventNumber_t event() const
Definition: EventID.h:117
double getChiSqu() const
Float_t e
Definition: plot.C:34
Float_t track
Definition: plot.C:34
std::vector< double > fPosErr
EventID id() const
Definition: Event.h:56
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
void trkf::Track3DKalman::reconfigure ( fhicl::ParameterSet const &  p)

Definition at line 164 of file Track3DKalman_module.cc.

References fG4ModuleLabel, fGenfPRINT, fGenieGenModuleLabel, fMomErr, fMomStart, fPosErr, fSpacePtsModuleLabel, and fhicl::ParameterSet::get().

Referenced by Track3DKalman().

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 }
std::vector< double > fMomStart
std::vector< double > fMomErr
std::vector< double > fPosErr
void art::Consumer::showMissingConsumes ( ) const
protectedinherited

Definition at line 125 of file Consumer.cc.

Referenced by art::EDProducer::doEndJob(), art::EDFilter::doEndJob(), art::EDAnalyzer::doEndJob(), and art::RootOutput::endJob().

126 {
127  if (!moduleContext_)
128  return;
129 
130  // If none of the branches have missing consumes statements, exit early.
131  if (std::all_of(cbegin(missingConsumes_),
132  cend(missingConsumes_),
133  [](auto const& perBranch) { return perBranch.empty(); }))
134  return;
135 
136  constexpr cet::HorizontalRule rule{60};
137  mf::LogPrint log{"MTdiagnostics"};
138  log << '\n'
139  << rule('=') << '\n'
140  << "The following consumes (or mayConsume) statements are missing from\n"
141  << module_context(moduleDescription_) << '\n'
142  << rule('-') << '\n';
143 
144  cet::for_all_with_index(
145  missingConsumes_, [&log](std::size_t const i, auto const& perBranch) {
146  for (auto const& pi : perBranch) {
147  log << " "
148  << assemble_consumes_statement(static_cast<BranchType>(i), pi)
149  << '\n';
150  }
151  });
152  log << rule('=');
153 }
cet::exempt_ptr< ModuleDescription const > moduleDescription_
Definition: Consumer.h:140
constexpr T pi()
Returns the constant pi (up to 35 decimal digits of precision)
bool moduleContext_
Definition: Consumer.h:136
ConsumableProductSets missingConsumes_
Definition: Consumer.h:139
void art::Consumer::validateConsumedProduct ( BranchType const  bt,
ProductInfo const &  pi 
)
protectedinherited

Definition at line 101 of file Consumer.cc.

References art::errors::ProductRegistrationFailure.

103 {
104  // Early exits if consumes tracking has been disabled or if the
105  // consumed product is an allowed consumable.
106  if (!moduleContext_)
107  return;
108 
109  if (cet::binary_search_all(consumables_[bt], pi))
110  return;
111 
112  if (requireConsumes_) {
114  "Consumer: an error occurred during validation of a "
115  "retrieved product\n\n")
116  << "The following consumes (or mayConsume) statement is missing from\n"
117  << module_context(moduleDescription_) << ":\n\n"
118  << " " << assemble_consumes_statement(bt, pi) << "\n\n";
119  }
120 
121  missingConsumes_[bt].insert(pi);
122 }
cet::exempt_ptr< ModuleDescription const > moduleDescription_
Definition: Consumer.h:140
bool requireConsumes_
Definition: Consumer.h:137
cet::coded_exception< errors::ErrorCodes, ExceptionDetail::translate > Exception
Definition: Exception.h:66
constexpr T pi()
Returns the constant pi (up to 35 decimal digits of precision)
ConsumableProducts consumables_
Definition: Consumer.h:138
bool moduleContext_
Definition: Consumer.h:136
ConsumableProductSets missingConsumes_
Definition: Consumer.h:139

Member Data Documentation

Float_t trkf::Track3DKalman::chi2
private

Definition at line 109 of file Track3DKalman_module.cc.

Referenced by beginJob(), and produce().

Float_t trkf::Track3DKalman::chi2ndf
private

Definition at line 110 of file Track3DKalman_module.cc.

Referenced by beginJob(), and produce().

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

Definition at line 106 of file Track3DKalman_module.cc.

Referenced by beginJob(), and produce().

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

Definition at line 108 of file Track3DKalman_module.cc.

Referenced by beginJob(), and produce().

unsigned int trkf::Track3DKalman::evtt
private

Definition at line 118 of file Track3DKalman_module.cc.

Referenced by beginJob(), and produce().

unsigned int trkf::Track3DKalman::fDimSize
private

Definition at line 124 of file Track3DKalman_module.cc.

Referenced by beginJob(), and produce().

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

Definition at line 99 of file Track3DKalman_module.cc.

Referenced by reconfigure().

bool trkf::Track3DKalman::fGenfPRINT
private

Definition at line 100 of file Track3DKalman_module.cc.

Referenced by produce(), and reconfigure().

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

Definition at line 98 of file Track3DKalman_module.cc.

Referenced by produce(), and reconfigure().

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

Definition at line 127 of file Track3DKalman_module.cc.

Referenced by produce(), and reconfigure().

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

Definition at line 128 of file Track3DKalman_module.cc.

Referenced by produce(), and reconfigure().

Float_t* trkf::Track3DKalman::fpMCT
private

Definition at line 115 of file Track3DKalman_module.cc.

Referenced by beginJob(), and produce().

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

Definition at line 126 of file Track3DKalman_module.cc.

Referenced by produce(), and reconfigure().

Float_t* trkf::Track3DKalman::fpREC
private

Definition at line 114 of file Track3DKalman_module.cc.

Referenced by beginJob(), and produce().

Float_t* trkf::Track3DKalman::fpRECL
private

Definition at line 113 of file Track3DKalman_module.cc.

Referenced by beginJob(), and produce().

Float_t* trkf::Track3DKalman::fpRECt3D
private

Definition at line 112 of file Track3DKalman_module.cc.

Referenced by beginJob(), and produce().

unsigned int trkf::Track3DKalman::fptsNo
private

Definition at line 120 of file Track3DKalman_module.cc.

Referenced by beginJob(), and produce().

Float_t* trkf::Track3DKalman::fshx
private

Definition at line 121 of file Track3DKalman_module.cc.

Referenced by beginJob(), and produce().

Float_t* trkf::Track3DKalman::fshy
private

Definition at line 122 of file Track3DKalman_module.cc.

Referenced by beginJob(), and produce().

Float_t* trkf::Track3DKalman::fshz
private

Definition at line 123 of file Track3DKalman_module.cc.

Referenced by beginJob(), and produce().

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

Definition at line 97 of file Track3DKalman_module.cc.

Referenced by produce(), and reconfigure().

int trkf::Track3DKalman::ndf
private

Definition at line 117 of file Track3DKalman_module.cc.

Referenced by beginJob(), and produce().

int trkf::Track3DKalman::nfail
private

Definition at line 116 of file Track3DKalman_module.cc.

Referenced by beginJob(), and produce().

unsigned int trkf::Track3DKalman::nTrks
private

Definition at line 119 of file Track3DKalman_module.cc.

Referenced by beginJob(), and produce().

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

Definition at line 130 of file Track3DKalman_module.cc.

Referenced by endJob(), and produce().

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

Definition at line 129 of file Track3DKalman_module.cc.

Referenced by endJob(), and produce().

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

Definition at line 105 of file Track3DKalman_module.cc.

Referenced by beginJob(), and produce().

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

Definition at line 107 of file Track3DKalman_module.cc.

Referenced by beginJob(), and produce().

TTree* trkf::Track3DKalman::tree
private

Definition at line 103 of file Track3DKalman_module.cc.

Referenced by beginJob(), and produce().


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