LArSoft  v07_13_02
Liquid Argon Software toolkit - http://larsoft.org/
trkf::Track3Dreco Class Reference
Inheritance diagram for trkf::Track3Dreco:
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

 Track3Dreco (fhicl::ParameterSet const &pset)
 
 ~Track3Dreco ()
 
void reconfigure (fhicl::ParameterSet const &p)
 
void produce (art::Event &evt)
 
void beginJob ()
 
void endJob ()
 
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

int ftmatch
 tolerance for time matching (in time samples) More...
 
double fchi2dof
 tolerance for chi2/dof of cluster fit to function More...
 
std::string fClusterModuleLabel
 label for input cluster collection More...
 

Detailed Description

Definition at line 56 of file Track3Dreco_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::Track3Dreco::Track3Dreco ( fhicl::ParameterSet const &  pset)
explicit

Definition at line 83 of file Track3Dreco_module.cc.

References reconfigure().

84 {
85  this->reconfigure(pset);
86  produces< std::vector<recob::Track> >();
87  produces< std::vector<recob::SpacePoint> >();
88  produces< art::Assns<recob::Track, recob::Cluster> >();
89  produces< art::Assns<recob::Track, recob::SpacePoint> >();
90  produces< art::Assns<recob::SpacePoint, recob::Hit> >();
91  produces< art::Assns<recob::Track, recob::Hit> >();
92 }
void reconfigure(fhicl::ParameterSet const &p)
trkf::Track3Dreco::~Track3Dreco ( )

Definition at line 95 of file Track3Dreco_module.cc.

96 {
97 }

Member Function Documentation

void trkf::Track3Dreco::beginJob ( )
virtual

Reimplemented from art::EDProducer.

Definition at line 107 of file Track3Dreco_module.cc.

108 {
109 
110 
111 }
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::Track3Dreco::endJob ( )
virtual

Reimplemented from art::EDProducer.

Definition at line 113 of file Track3Dreco_module.cc.

114 {
115 }
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::Track3Dreco::produce ( art::Event evt)
virtual
Todo:
: This is very bad practice and should be changed ASAP

Implements art::EDProducer.

Definition at line 118 of file Track3Dreco_module.cc.

References recob::tracking::convertCollToPoint(), recob::tracking::convertCollToVector(), util::CreateAssn(), DEFINE_ART_MODULE, geo::GeometryCore::DetHalfHeight(), detinfo::DetectorProperties::DriftVelocity(), fClusterModuleLabel, ftmatch, art::DataViewImpl::getByLabel(), geo::GeometryCore::GetLArTPCVolumeName(), hits(), geo::kZ, geo::GeometryCore::Plane(), geo::GeometryCore::PlanePitch(), art::PtrVector< T >::push_back(), art::Event::put(), art::PtrVector< T >::size(), t1, t2, geo::WireGeo::ThetaZ(), recob::Cluster::View(), geo::PlaneGeo::Wire(), and geo::GeometryCore::WirePitch().

119 {
120  // get services
122  const detinfo::DetectorProperties* detprop = lar::providerFrom<detinfo::DetectorPropertiesService>();
123 
124  std::unique_ptr<std::vector<recob::Track> > tcol(new std::vector<recob::Track>);
125  std::unique_ptr<std::vector<recob::SpacePoint> > spacepoints(new std::vector<recob::SpacePoint>);
126  std::unique_ptr< art::Assns<recob::Track, recob::Cluster> > cassn (new art::Assns<recob::Track, recob::Cluster>);
127  std::unique_ptr< art::Assns<recob::Track, recob::SpacePoint> > sassn (new art::Assns<recob::Track, recob::SpacePoint>);
128  std::unique_ptr< art::Assns<recob::SpacePoint, recob::Hit> > shassn(new art::Assns<recob::SpacePoint, recob::Hit>);
129  std::unique_ptr< art::Assns<recob::Track, recob::Hit> > hassn (new art::Assns<recob::Track, recob::Hit>);
130 
131  // define TPC parameters
132  TString tpcName = geom->GetLArTPCVolumeName();
133 
134  // double YC = (m_TPCHalfZ-5.)*2.; // TPC height in cm
135  double YC = (geom->DetHalfHeight())*2.; // *ArgoNeuT* TPC active-volume height in cm
136  double Angle = geom->Plane(1).Wire(0).ThetaZ(false)-TMath::Pi()/2.; // wire angle with respect to the vertical direction
137  // Parameters temporary defined here, but possibly to be retrieved somewhere in the code
138  double timetick = 0.198; //time sample in us
139  double presamplings = 60.;
140  const double wireShift=50.; // half the number of wires from the Induction(Collection) plane intersecting with a wire from the Collection(Induction) plane.
141  double plane_pitch = geom->PlanePitch(0,1); //wire plane pitch in cm
142  double wire_pitch = geom->WirePitch(); //wire pitch in cm
143  double Efield_drift = 0.5; // Electric Field in the drift region in kV/cm
144  double Efield_SI = 0.7; // Electric Field between Shield and Induction planes in kV/cm
145  double Efield_IC = 0.9; // Electric Field between Induction and Collection planes in kV/cm
146  double Temperature = 87.6; // LAr Temperature in K
147 
148  double driftvelocity = detprop->DriftVelocity(Efield_drift,Temperature); //drift velocity in the drift
149  //region (cm/us)
150  double driftvelocity_SI = detprop->DriftVelocity(Efield_SI,Temperature); //drift velocity between shield
151  //and induction (cm/us)
152  double driftvelocity_IC = detprop->DriftVelocity(Efield_IC,Temperature); //drift velocity between induction
153  //and collection (cm/us)
154  double timepitch = driftvelocity*timetick; //time sample (cm)
155  double tSI = plane_pitch/driftvelocity_SI/timetick; //drift time between Shield and
156  //Collection planes (time samples)
157  double tIC = plane_pitch/driftvelocity_IC/timetick; //drift time between Induction and
158  //Collection planes (time samples)
159 
160 
161  // get input Cluster object(s).
162  art::Handle< std::vector<recob::Cluster> > clusterListHandle;
163  evt.getByLabel(fClusterModuleLabel,clusterListHandle);
164 
165 
166  // Declare some vectors..
167  // Induction
168  std::vector<double> Iwirefirsts; // in cm
169  std::vector<double> Iwirelasts; // in cm
170  std::vector<double> Itimefirsts; // in cm
171  std::vector<double> Itimelasts; // in cm
172  std::vector<double> Itimefirsts_line; // in cm
173  std::vector<double> Itimelasts_line; // in cm
174  std::vector < std::vector<art::Ptr<recob::Hit> > > IclusHitlists;
175  std::vector<unsigned int> Icluster_count;
176 
177  // Collection
178  std::vector<double> Cwirefirsts; // in cm
179  std::vector<double> Cwirelasts; // in cm
180  std::vector<double> Ctimefirsts; // in cm
181  std::vector<double> Ctimelasts; // in cm
182  std::vector<double> Ctimefirsts_line; // in cm
183  std::vector<double> Ctimelasts_line; // in cm
184  std::vector< std::vector< art::Ptr<recob::Hit> > > CclusHitlists;
185  std::vector<unsigned int> Ccluster_count;
186 
187  // Some variables for the hit
188  float time; //hit time at maximum
189  unsigned int wire = 0; //hit wire number
190  unsigned int plane = 0; //hit plane number
191 
192  size_t startSPIndex = spacepoints->size(); //index for knowing which spacepoints are with which cluster
193  size_t endSPIndex = spacepoints->size(); //index for knowing which spacepoints are with which cluster
194 
195  art::FindManyP<recob::Hit> fmh(clusterListHandle, evt, fClusterModuleLabel);
196 
197  for(size_t ii = 0; ii < clusterListHandle->size(); ++ii){
198 
199  art::Ptr<recob::Cluster> cl(clusterListHandle, ii);
200 
204 
205  // Gaaaaaah! Change me soon!!! But, for now,
206  // let's just chuck one plane's worth of info. EC, 30-Mar-2011.
208  if (cl->View() == geo::kZ) continue;
209 
210  // Can not be const, cuz we're gonna sort 'em.
211  std::vector< art::Ptr<recob::Hit> > hitlist (fmh.at(ii));
212 
213  if(hitlist.size() == 1) continue;//only one Hit in this Cluster...will cause TGraph fit to fail.
214 
215  // sort the hit list to be sure it is in the correct order
216  // using the Hit < operator
217  std::sort(hitlist.begin(), hitlist.end());
218 
219  TGraph *the2Dtrack = new TGraph(hitlist.size());
220 
221  std::vector<double> wires;
222  std::vector<double> times;
223 
224  int np = 0;
225  //loop over cluster hits
226  for(art::PtrVector<recob::Hit>::const_iterator theHit = hitlist.begin(); theHit != hitlist.end(); theHit++){
227  //recover the Hit
228  // recob::Hit* theHit = (recob::Hit*)(*hitIter);
229  time = (*theHit)->PeakTime() ;
230 
231  time -= presamplings;
232 
233  plane = (*theHit)->WireID().Plane;
234  wire = (*theHit)->WireID().Wire;
235 
236  //correct for the distance between wire planes
237  if(plane == 1) time -= tIC; // Collection
238 
239  //transform hit wire and time into cm
240  double wire_cm;
241  if(plane == 0)
242  wire_cm = (double)((wire+3.95) * wire_pitch);
243  else
244  wire_cm = (double)((wire+1.84) * wire_pitch);
245 
246  double time_cm;
247  if(time > tSI) time_cm = (double)( (time-tSI)*timepitch + tSI*driftvelocity_SI*timetick);
248  else time_cm = time*driftvelocity_SI*timetick;
249 
250  wires.push_back(wire_cm);
251  times.push_back(time_cm);
252 
253  the2Dtrack->SetPoint(np,wire_cm,time_cm);
254  np++;
255  }//end of loop over cluster hits
256 
257  // fit the 2Dtrack and get some info to store
258  try{
259  the2Dtrack->Fit("pol1","Q");
260  }
261  catch(...){
262  mf::LogWarning("Track3Dreco") << "The 2D track fit failed";
263  continue;
264  }
265 
266  TF1 *pol1=(TF1*) the2Dtrack->GetFunction("pol1");
267  double par[2];
268  pol1->GetParameters(par);
269  double intercept = par[0];
270  double slope = par[1];
271 
272  double w0 = wires.front(); // first hit wire (cm)
273  double w1 = wires.back(); // last hit wire (cm)
274  double t0 = times.front(); // first hit time (cm)
275  double t1 = times.back(); // last hit time (cm)
276  double t0_line = intercept + (w0)*slope; // time coordinate at wire w0 on the fit line (cm)
277  double t1_line = intercept + (w1)*slope; // time coordinate at wire w1 on the fit line (cm)
278 
279  // actually store the 2Dtrack info
280  switch(plane){
281  case 0:
282  Iwirefirsts.push_back(w0);
283  Iwirelasts.push_back(w1);
284  Itimefirsts.push_back(t0);
285  Itimelasts.push_back(t1);
286  Itimefirsts_line.push_back(t0_line);
287  Itimelasts_line.push_back(t1_line);
288  IclusHitlists.push_back(hitlist);
289  Icluster_count.push_back(ii);
290  break;
291  case 1:
292  Cwirefirsts.push_back(w0);
293  Cwirelasts.push_back(w1);
294  Ctimefirsts.push_back(t0);
295  Ctimelasts.push_back(t1);
296  Ctimefirsts_line.push_back(t0_line);
297  Ctimelasts_line.push_back(t1_line);
298  CclusHitlists.push_back(hitlist);
299  Ccluster_count.push_back(ii);
300  break;
301  }
302  //delete the2Dtrack;
303  delete pol1;
304  }// end of loop over all input clusters
305 
309 
310  //loop over Collection view 2D tracks
311  for(size_t collectionIter = 0; collectionIter < CclusHitlists.size(); ++collectionIter){
312  // Recover previously stored info
313  double Cw0 = Cwirefirsts[collectionIter];
314  double Cw1 = Cwirelasts[collectionIter];
315  //double Ct0 = Ctimefirsts[collectionIter];
316  //double Ct1 = Ctimelasts[collectionIter];
317  double Ct0_line = Ctimefirsts_line[collectionIter];
318  double Ct1_line = Ctimelasts_line[collectionIter];
319  std::vector< art::Ptr<recob::Hit> > hitsCtrk = CclusHitlists[collectionIter];
320 
321  double collLength = TMath::Sqrt( TMath::Power(Ct1_line - Ct0_line,2) + TMath::Power(Cw1 - Cw0,2));
322 
323  //loop over Induction view 2D tracks
324  for(size_t inductionIter = 0; inductionIter < IclusHitlists.size(); ++inductionIter){
325  // Recover previously stored info
326  double Iw0 = Iwirefirsts[inductionIter];
327  double Iw1 = Iwirelasts[inductionIter];
328  //double It0 = Itimefirsts[inductionIter];
329  //double It1 = Itimelasts[inductionIter];
330  double It0_line = Itimefirsts_line[inductionIter];
331  double It1_line = Itimelasts_line[inductionIter];
332  std::vector< art::Ptr<recob::Hit> > hitsItrk = IclusHitlists[inductionIter];
333 
334  double indLength = TMath::Sqrt( TMath::Power(It1_line - It0_line,2) + TMath::Power(Iw1 - Iw0,2));
335 
336  bool forward_match = ((std::abs(Ct0_line-It0_line)<ftmatch*timepitch) &&
337  (std::abs(Ct1_line-It1_line)<ftmatch*timepitch));
338  bool backward_match = ((std::abs(Ct0_line-It1_line)<ftmatch*timepitch) &&
339  (std::abs(Ct1_line-It0_line)<ftmatch*timepitch));
340 
341 
342 
343  // match 2D tracks
344  if(forward_match || backward_match ){
345 
346  // Reconstruct the 3D track
347  TVector3 XYZ0, XYZ1; // track endpoints
348  if(forward_match){
349  XYZ0.SetXYZ(Ct0_line,(Cw0-Iw0)/(2.*TMath::Sin(Angle)),
350  (Cw0+Iw0)/(2.*TMath::Cos(Angle))-YC/2.*TMath::Tan(Angle));
351  XYZ1.SetXYZ(Ct1_line,(Cw1-Iw1)/(2.*TMath::Sin(Angle)),
352  (Cw1+Iw1)/(2.*TMath::Cos(Angle))-YC/2.*TMath::Tan(Angle));
353  }
354  else{
355  XYZ0.SetXYZ(Ct0_line,(Cw0-Iw1)/(2.*TMath::Sin(Angle)),
356  (Cw0+Iw1)/(2.*TMath::Cos(Angle))-YC/2.*TMath::Tan(Angle));
357  XYZ1.SetXYZ(Ct1_line,(Cw1-Iw0)/(2.*TMath::Sin(Angle)),
358  (Cw1+Iw0)/(2.*TMath::Cos(Angle))-YC/2.*TMath::Tan(Angle));
359  }
360 
361  //compute track direction in Local co-ordinate system
362  //WARNING: There is an ambiguity introduced here for the case of backwards-going tracks.
363  //If available, vertex info. could sort this out.
364  TVector3 startpointVec,endpointVec;
365  TVector2 collVtx, indVtx;
366  if(XYZ0.Z() <= XYZ1.Z()){
367  startpointVec.SetXYZ(XYZ0.X(),XYZ0.Y(),XYZ0.Z());
368  endpointVec.SetXYZ(XYZ1.X(),XYZ1.Y(),XYZ1.Z());
369  if(forward_match){
370  collVtx.Set(Ct0_line,Cw0);
371  indVtx.Set(It0_line,Iw0);
372  }
373  else{
374  collVtx.Set(Ct0_line,Cw0);
375  indVtx.Set(It1_line,Iw1);
376  }
377  }
378  else{
379  startpointVec.SetXYZ(XYZ1.X(),XYZ1.Y(),XYZ1.Z());
380  endpointVec.SetXYZ(XYZ0.X(),XYZ0.Y(),XYZ0.Z());
381  if(forward_match){
382  collVtx.Set(Ct1_line,Cw1);
383  indVtx.Set(It1_line,Iw1);
384  }
385  else{
386  collVtx.Set(Ct1_line,Cw1);
387  indVtx.Set(It0_line,Iw0);
388  }
389  }
390 
391  //compute track (normalized) cosine directions in the TPC co-ordinate system
392  TVector3 DirCos = endpointVec - startpointVec;
393 
394  //SetMag casues a crash if the magnitude of the vector is zero
395  try{
396  DirCos.SetMag(1.0);//normalize vector
397  }
398  catch(...){
399  mf::LogWarning("Track3Dreco") <<"The Spacepoint is infinitely small";
400  continue;
401  }
402 
403  art::Ptr <recob::Cluster> cl1(clusterListHandle,Icluster_count[inductionIter]);
404  art::Ptr <recob::Cluster> cl2(clusterListHandle,Ccluster_count[collectionIter]);
405  art::PtrVector<recob::Cluster> clustersPerTrack;
406  clustersPerTrack.push_back(cl1);
407  clustersPerTrack.push_back(cl2);
408 
409 
411  // Match hits
413 
414  std::vector< art::Ptr<recob::Hit> > minhits = hitsCtrk.size() <= hitsItrk.size() ? hitsCtrk : hitsItrk;
415  std::vector< art::Ptr<recob::Hit> > maxhits = hitsItrk.size() <= hitsCtrk.size() ? hitsCtrk : hitsItrk;
416 
417 
418  std::vector<bool> maxhitsMatch(maxhits.size());
419  for(size_t it = 0; it < maxhits.size(); ++it) maxhitsMatch[it] = false;
420 
421  std::vector<recob::Hit*> hits3Dmatched;
422  // For the matching start from the view where the track projection presents less hits
423  unsigned int imaximum = 0;
424  //loop over hits
425 
426  startSPIndex = spacepoints->size();
427 
428  for(size_t imin = 0; imin < minhits.size(); ++imin){
429  //get wire - time coordinate of the hit
430  unsigned int wire,plane1,plane2;
431  wire = minhits[imin]->WireID().Wire;
432  plane1 = minhits[imin]->WireID().Plane;
433 
434  // get the wire-time co-ordinates of the hit to be matched
435  double w1;
436  if(plane1 == 0)
437  w1 = (double)((wire+3.95)*wire_pitch);
438  else
439  w1 = (double)((wire+1.84) * wire_pitch);
440  double temptime1 = minhits[imin]->PeakTime()-presamplings;
441  if(plane1 == 1) temptime1 -= tIC;
442  double t1;
443  if(temptime1>tSI) t1 = (double)( (temptime1-tSI)*timepitch + tSI*driftvelocity_SI*timetick);
444  else t1 = temptime1*driftvelocity_SI*timetick;
445 
446  //get the track origin co-ordinates in the two views
447  TVector2 minVtx2D;
448  plane1==1 ? minVtx2D.Set(collVtx.X(),collVtx.Y()): minVtx2D.Set(indVtx.X(),indVtx.Y());
449  TVector2 maxVtx2D;
450  plane1==1 ? maxVtx2D.Set(indVtx.X(),indVtx.Y()): maxVtx2D.Set(collVtx.X(),collVtx.Y());
451 
452  double ratio = (collLength>indLength) ? collLength/indLength : indLength/collLength;
453 
454  //compute the distance of the hit (imin) from the relative track origin
455  double minDistance = ratio*TMath::Sqrt(TMath::Power(t1-minVtx2D.X(),2)
456  + TMath::Power(w1-minVtx2D.Y(),2));
457 
458  //core matching algorithm
459  double difference = 9999999.;
460 
461  //loop over hits of the other view
462  for(size_t imax = 0; imax < maxhits.size(); ++imax){
463  if(!maxhitsMatch[imax]){
464  //get wire - time coordinate of the hit
465  wire = maxhits[imax]->WireID().Wire;
466  plane2 = maxhits[imax]->WireID().Plane;
467 
468  double w2;
469  if(plane2 == 0)
470  w2 = (double)((wire+3.95)*wire_pitch);
471  else
472  w2 = (double)((wire+1.84)*wire_pitch);
473  double temptime2 = maxhits[imax]->PeakTime()-presamplings;
474  if(plane2 == 1) temptime2 -= tIC;
475  double t2;
476  if(temptime2 > tSI) t2 = (double)( (temptime2-tSI)*timepitch + tSI*driftvelocity_SI*timetick);
477  else t2 = temptime2*driftvelocity_SI*timetick;
478 
479 
480  bool timematch = (std::abs(t1-t2)<ftmatch*timepitch);
481  bool wirematch = (std::abs(w1-w2)<wireShift*wire_pitch);
482 
483  double maxDistance = TMath::Sqrt(TMath::Power(t2-maxVtx2D.X(),2)+TMath::Power(w2-maxVtx2D.Y(),2));
484  if (wirematch && timematch && std::abs(maxDistance-minDistance)<difference) {
485  difference = std::abs(maxDistance-minDistance);
486  imaximum = imax;
487  }
488  }
489  }// end loop over max hits
490  maxhitsMatch[imaximum]=true;
491 
493  if(difference!= 9999999.){
494  sp_hits.push_back(minhits[imin]);
495  sp_hits.push_back(maxhits[imaximum]);
496  }
497 
498  // Get the time-wire co-ordinates of the matched hit
499  wire = maxhits[imaximum]->WireID().Wire;
500  plane2 = maxhits[imaximum]->WireID().Plane;
501 
502  double w1_match;
503  if(plane2 == 0)
504  w1_match = (double)((wire+3.95)*wire_pitch);
505  else
506  w1_match = (double)((wire+1.84)*wire_pitch);
507  double temptime3 = maxhits[imaximum]->PeakTime()-presamplings;
508  if(plane2 == 1) temptime3 -= tIC;
509  double t1_match;
510  if(temptime3 > tSI) t1_match = (double)( (temptime3-tSI)*timepitch + tSI*driftvelocity_SI*timetick);
511  else t1_match = temptime3*driftvelocity_SI*timetick;
512 
513 
514  // create the 3D hit, compute its co-ordinates and add it to the 3D hits list
515  double Ct = plane1==1?t1:t1_match;
516  double Cw = plane1==1?w1:w1_match;
517  double Iw = plane1==1?w1_match:w1;
518 
519  const TVector3 hit3d(Ct,(Cw-Iw)/(2.*TMath::Sin(Angle)),(Cw+Iw)/(2.*TMath::Cos(Angle))-YC/2.*TMath::Tan(Angle));
520  Double_t hitcoord[3];
521  hitcoord[0] = hit3d.X();
522  hitcoord[1] = hit3d.Y();
523  hitcoord[2] = hit3d.Z();
524 
525  Double_t hitcoord_errs[3];
526  for (int i=0; i<3; i++) hitcoord_errs[i]=-1.000;
527 
528  //3d point at end of track
529  recob::SpacePoint mysp(hitcoord, hitcoord_errs, -1., spacepoints->size());
530 
531  spacepoints->push_back(mysp);
532 
533  // associate the hits to the space point
534  util::CreateAssn(*this, evt, *spacepoints, sp_hits, *shassn);
535 
536  }//loop over min-hits
537 
538  endSPIndex = spacepoints->size();
539 
540  // Add the 3D track to the vector of the reconstructed tracks
541  if(spacepoints->size() > startSPIndex || clustersPerTrack.size()>0){
542 
543  std::vector<TVector3> xyz;
544  xyz.push_back(startpointVec);
545  xyz.push_back(endpointVec);
546  std::vector<TVector3> dir_xyz;
547  dir_xyz.push_back(DirCos);
548  dir_xyz.push_back(DirCos);
549 
552  recob::Track::Flags_t(xyz.size()), false),
553  0, -1., 0, recob::tracking::SMatrixSym55(), recob::tracking::SMatrixSym55(), tcol->size());
554  tcol->push_back(the3DTrack);
555 
556  // associate the track with its spacepoints
557  util::CreateAssn(*this, evt, *tcol, *spacepoints, *sassn, startSPIndex, endSPIndex);
558 
559  // associate the track with its clusters
560  util::CreateAssn(*this, evt, *tcol, clustersPerTrack, *cassn);
561 
562  art::FindManyP<recob::Hit> fmhc(clustersPerTrack, evt, fClusterModuleLabel);
563 
564  // get the hits associated with each cluster and associate those with the track
565  for(size_t p = 0; p < clustersPerTrack.size(); ++p){
566  const std::vector< art::Ptr<recob::Hit> >& hits = fmhc.at(p);
567  util::CreateAssn(*this, evt, *tcol, hits, *hassn);
568  }
569 
570  }
571 
572  } //close match 2D tracks
573 
574 
575  }//close loop over Induction view 2D tracks
576 
577  }//close loop over Collection view 2D tracks
578 
579  mf::LogVerbatim("Summary") << std::setfill('-') << std::setw(175) << "-" << std::setfill(' ');
580  mf::LogVerbatim("Summary") << "Track3Dreco Summary:";
581  for(unsigned int i = 0; i<tcol->size(); ++i){
582  mf::LogVerbatim("Summary") << tcol->at(i) ;
583  }
584 
585  evt.put(std::move(tcol));
586  evt.put(std::move(spacepoints));
587  evt.put(std::move(cassn));
588  evt.put(std::move(sassn));
589  evt.put(std::move(hassn));
590  evt.put(std::move(shassn));
591 
592  return;
593 }
code to link reconstructed objects back to the MC truth information
MaybeLogger_< ELseverityLevel::ELsev_info, true > LogVerbatim
PlaneGeo const & Plane(unsigned int const p, unsigned int const tpc=0, unsigned int const cstat=0) const
Returns the specified wire.
TTree * t1
Definition: plottest35.C:26
geo::Length_t PlanePitch(geo::TPCID const &tpcid, geo::PlaneID::PlaneID_t p1=0, geo::PlaneID::PlaneID_t p2=1) const
Returns the distance between two planes.
WireGeo const & Wire(unsigned int iwire) const
Definition: PlaneGeo.cxx:506
ROOT::Math::SMatrix< Double32_t, 5, 5, ROOT::Math::MatRepSym< Double32_t, 5 > > SMatrixSym55
Definition: TrackingTypes.h:85
std::vector< Vector_t > convertCollToVector(std::vector< Vector > const &coll)
Definition: TrackingTypes.h:58
std::string fClusterModuleLabel
label for input cluster collection
TrackTrajectory::Flags_t Flags_t
Definition: Track.h:71
Planes which measure Z direction.
Definition: geo_types.h:79
double ThetaZ() const
Returns angle of wire with respect to z axis in the Y-Z plane in radians.
Definition: WireGeo.h:192
geo::Length_t WirePitch(geo::PlaneID const &planeid) const
Returns the distance between two consecutive wires.
ProductID put(std::unique_ptr< PROD > &&product)
Definition: Event.h:102
void hits()
Definition: readHits.C:15
geo::Length_t DetHalfHeight(geo::TPCID const &tpcid) const
Returns the half height of the active volume of the specified TPC.
void push_back(Ptr< U > const &p)
Definition: PtrVector.h:441
A trajectory in space reconstructed from hits.
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
size_type size() const
Definition: PtrVector.h:308
TTree * t2
Definition: plottest35.C:36
int ftmatch
tolerance for time matching (in time samples)
data_t::const_iterator const_iterator
Definition: PtrVector.h:61
bool getByLabel(std::string const &label, std::string const &productInstanceName, Handle< PROD > &result) const
Definition: DataViewImpl.h:344
virtual double DriftVelocity(double efield=0., double temperature=0.) const =0
MaybeLogger_< ELseverityLevel::ELsev_warning, false > LogWarning
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
std::string GetLArTPCVolumeName(geo::TPCID const &tpcid) const
Return the name of specified LAr TPC volume.
void trkf::Track3Dreco::reconfigure ( fhicl::ParameterSet const &  p)

Definition at line 99 of file Track3Dreco_module.cc.

References fchi2dof, fClusterModuleLabel, ftmatch, and fhicl::ParameterSet::get().

Referenced by Track3Dreco().

100 {
101  fClusterModuleLabel = pset.get< std::string >("ClusterModuleLabel");
102  ftmatch = pset.get< int >("TMatch");
103  fchi2dof = pset.get< double >("Chi2DOFmax");
104 }
std::string fClusterModuleLabel
label for input cluster collection
double fchi2dof
tolerance for chi2/dof of cluster fit to function
int ftmatch
tolerance for time matching (in time samples)
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

double trkf::Track3Dreco::fchi2dof
private

tolerance for chi2/dof of cluster fit to function

Definition at line 72 of file Track3Dreco_module.cc.

Referenced by reconfigure().

std::string trkf::Track3Dreco::fClusterModuleLabel
private

label for input cluster collection

Definition at line 73 of file Track3Dreco_module.cc.

Referenced by produce(), and reconfigure().

int trkf::Track3Dreco::ftmatch
private

tolerance for time matching (in time samples)

Definition at line 71 of file Track3Dreco_module.cc.

Referenced by produce(), and reconfigure().


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