All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Private Member Functions | Private Attributes | List of all members
SBNDOpT0Finder Class Reference
Inheritance diagram for SBNDOpT0Finder:

Public Member Functions

 SBNDOpT0Finder (fhicl::ParameterSet const &p)
 
 SBNDOpT0Finder (SBNDOpT0Finder const &)=delete
 
 SBNDOpT0Finder (SBNDOpT0Finder &&)=delete
 
SBNDOpT0Finderoperator= (SBNDOpT0Finder const &)=delete
 
SBNDOpT0Finderoperator= (SBNDOpT0Finder &&)=delete
 
void produce (art::Event &e) override
 

Private Member Functions

void DoMatch (art::Event &e, int tpc, std::unique_ptr< std::vector< anab::T0 >> &t0_v, std::unique_ptr< art::Assns< recob::Slice, anab::T0 >> &slice_t0_assn_v, std::unique_ptr< art::Assns< recob::OpFlash, anab::T0 >> &flash_t0_assn_v)
 Performs the matching in a specified tpc. More...
 
bool ConstructLightClusters (art::Event &e, unsigned int tpc)
 Constructs all the LightClusters (TPC Objects) in a specified TPC. More...
 
float GetNPhotons (const float charge, const art::Ptr< recob::PFParticle > &pfp)
 Returns the number of photons given charge and PFParticle. More...
 
std::vector< int > PDNamesToList (std::vector< std::string >)
 Convert from a list of PDS names to a list of op channels. More...
 
std::vector< int > GetUncoatedPTMList (std::vector< int > ch_to_use)
 Returns a list of uncoated PMTs that are a subset of those in ch_to_use. More...
 

Private Attributes

std::unique_ptr
< SemiAnalyticalModel
_semi_model
 
fhicl::ParameterSet _vuv_params
 
fhicl::ParameterSet _vis_params
 
::flashmatch::FlashMatchManager _mgr
 The flash matching manager. More...
 
std::vector
< flashmatch::FlashMatch_t
_result_v
 Matching result will be stored here. More...
 
std::vector< std::string > _opflash_producer_v
 The OpFlash producers (to be set) More...
 
std::vector< unsigned int > _tpc_v
 TPC number per OpFlash producer (to be set) More...
 
std::string _slice_producer
 The Slice producer (to be set) More...
 
double _flash_trange_start
 The time start from where to include flashes (to be set) More...
 
double _flash_trange_end
 The time stop from where to stop including flashes (to be set) More...
 
float _charge_to_n_photons_track
 The conversion factor betweeen hit integral and photons (to be set) More...
 
float _charge_to_n_photons_shower
 The conversion factor betweeen hit integral and photons (to be set) More...
 
std::vector< std::string > _photo_detectors
 The photodetector to use (to be set) More...
 
std::vector< int > _opch_to_use
 List of opch to use (will be infered from _photo_detectors) More...
 
std::vector< int > _uncoated_pmts
 List of uncoated opch to use (will be infered from _opch_to_use) More...
 
opdet::sbndPDMapAlg _pds_map
 map for photon detector types More...
 
std::vector
< flashmatch::QCluster_t
_light_cluster_v
 Vector that contains all the TPC objects. More...
 
std::map< int, art::Ptr
< recob::Slice > > 
_clusterid_to_slice
 
std::map< int, art::Ptr
< recob::OpFlash > > 
_flashid_to_opflash
 Will contain map tpc object id -> Slice. More...
 
TTree * _tree1
 Will contain map flash id -> OpFlash. More...
 
int _run
 
int _subrun
 
int _event
 
int _tpc
 
int _matchid
 
int _flashid
 
int _tpcid
 
double _t0
 
double _score
 
double _tpc_xmin
 
double _qll_xmin
 
double _hypo_pe
 
double _flash_pe
 
std::vector< double > _flash_spec
 
std::vector< double > _hypo_spec
 
TTree * _tree2
 
std::vector< float > _dep_x
 
std::vector< float > _dep_y
 
std::vector< float > _dep_z
 
std::vector< float > _dep_charge
 
std::vector< float > _dep_n_photons
 
std::vector< int > _dep_slice
 

Detailed Description

Definition at line 56 of file SBNDOpT0Finder_module.cc.

Constructor & Destructor Documentation

SBNDOpT0Finder::SBNDOpT0Finder ( fhicl::ParameterSet const &  p)
explicit

Definition at line 137 of file SBNDOpT0Finder_module.cc.

138  : EDProducer{p}
139 {
140  produces<std::vector<anab::T0>>();
141  produces<art::Assns<recob::Slice, anab::T0>>();
142  produces<art::Assns<recob::OpFlash, anab::T0>>();
143 
144  ::art::ServiceHandle<geo::Geometry> geo;
145 
146  _vuv_params = p.get<fhicl::ParameterSet>("VUVHits");
147  _vis_params = p.get<fhicl::ParameterSet>("VIVHits");
148  _semi_model = std::make_unique<SemiAnalyticalModel>(_vuv_params, _vis_params, true, false);
149 
150  _opflash_producer_v = p.get<std::vector<std::string>>("OpFlashProducers");
151  _tpc_v = p.get<std::vector<unsigned int>>("TPCs");
152  _slice_producer = p.get<std::string>("SliceProducer");
153 
154  _flash_trange_start = p.get<double>("FlashVetoTimeStart", 0);
155  _flash_trange_end = p.get<double>("FlashVetoTimeEnd", 2);
156 
157  _photo_detectors = p.get<std::vector<std::string>>("PhotoDetectors");
158  _opch_to_use = this->PDNamesToList(_photo_detectors);
159  _uncoated_pmts = this->GetUncoatedPTMList(_opch_to_use);
160 
161  _charge_to_n_photons_track = p.get<float>("ChargeToNPhotonsTrack");
162  _charge_to_n_photons_shower = p.get<float>("ChargeToNPhotonsShower");
163 
164  if (_tpc_v.size() != _opflash_producer_v.size()) {
165  throw cet::exception("SBNDOpT0Finder")
166  << "TPC vector and OpFlash producer vector don't have the same size, check your fcl params.";
167  }
168 
169  _mgr.Configure(p.get<flashmatch::Config_t>("FlashMatchConfig"));
170 
171  _mgr.SetChannelMask(_opch_to_use);
172 
174 
175  _mgr.SetSemiAnalyticalModel(std::move(_semi_model));
176 
177  _flash_spec.resize(geo->NOpDets(), 0.);
178  _hypo_spec.resize(geo->NOpDets(), 0.);
179 
180  art::ServiceHandle<art::TFileService> fs;
181 
182  _tree1 = fs->make<TTree>("slice_deposition_tree","");
183  _tree1->Branch("run", &_run, "run/I");
184  _tree1->Branch("subrun", &_subrun, "subrun/I");
185  _tree1->Branch("event", &_event, "event/I");
186  _tree1->Branch("dep_slice", "std::vector<int>", &_dep_slice);
187  _tree1->Branch("dep_x", "std::vector<float>", &_dep_x);
188  _tree1->Branch("dep_y", "std::vector<float>", &_dep_y);
189  _tree1->Branch("dep_z", "std::vector<float>", &_dep_z);
190  _tree1->Branch("dep_charge", "std::vector<float>", &_dep_charge);
191  _tree1->Branch("dep_n_photons", "std::vector<float>", &_dep_n_photons);
192 
193  _tree2 = fs->make<TTree>("flash_match_tree","");
194  _tree2->Branch("run", &_run, "run/I");
195  _tree2->Branch("subrun", &_subrun, "subrun/I");
196  _tree2->Branch("event", &_event, "event/I");
197  _tree2->Branch("tpc", &_tpc, "tpc/I");
198  _tree2->Branch("matchid", &_matchid, "matchid/I");
199  _tree2->Branch("tpcid", &_tpcid, "tpcid/I");
200  _tree2->Branch("flashid", &_flashid, "flashid/I");
201  _tree2->Branch("tpc_xmin", &_tpc_xmin, "tpc_xmin/D");
202  _tree2->Branch("qll_xmin", &_qll_xmin, "qll_xmin/D");
203  _tree2->Branch("t0", &_t0, "t0/D");
204  _tree2->Branch("score", &_score, "score/D");
205  _tree2->Branch("hypo_pe", &_hypo_pe, "hypo_pe/D");
206  _tree2->Branch("flash_pe", &_flash_pe, "flash_pe/D");
207  _tree2->Branch("hypo_spec", "std::vector<double>", &_hypo_spec);
208  _tree2->Branch("flash_spec", "std::vector<double>", &_flash_spec);
209 }
double _flash_trange_end
The time stop from where to stop including flashes (to be set)
std::vector< float > _dep_z
std::vector< float > _dep_n_photons
std::unique_ptr< SemiAnalyticalModel > _semi_model
std::vector< int > _uncoated_pmts
List of uncoated opch to use (will be infered from _opch_to_use)
void SetSemiAnalyticalModel(std::unique_ptr< SemiAnalyticalModel > model)
std::vector< int > _dep_slice
std::vector< std::string > _opflash_producer_v
The OpFlash producers (to be set)
pdgs p
Definition: selectors.fcl:22
std::vector< float > _dep_charge
std::vector< std::string > _photo_detectors
The photodetector to use (to be set)
fhicl::ParameterSet Config_t
Configuration object.
Definition: FMWKInterface.h:31
std::vector< double > _flash_spec
::flashmatch::FlashMatchManager _mgr
The flash matching manager.
double _flash_trange_start
The time start from where to include flashes (to be set)
void SetChannelMask(std::vector< int >)
Sets the op channels to be used for matching.
void SetUncoatedPMTs(std::vector< int > ch_uncoated)
Sets the channels sensitive to visible light.
std::vector< unsigned int > _tpc_v
TPC number per OpFlash producer (to be set)
std::vector< int > PDNamesToList(std::vector< std::string >)
Convert from a list of PDS names to a list of op channels.
std::string _slice_producer
The Slice producer (to be set)
std::vector< int > GetUncoatedPTMList(std::vector< int > ch_to_use)
Returns a list of uncoated PMTs that are a subset of those in ch_to_use.
std::vector< float > _dep_x
float _charge_to_n_photons_track
The conversion factor betweeen hit integral and photons (to be set)
std::vector< float > _dep_y
std::vector< double > _hypo_spec
fhicl::ParameterSet _vis_params
TTree * _tree1
Will contain map flash id -&gt; OpFlash.
std::vector< int > _opch_to_use
List of opch to use (will be infered from _photo_detectors)
void Configure(const Config_t &cfg)
Configuration.
float _charge_to_n_photons_shower
The conversion factor betweeen hit integral and photons (to be set)
fhicl::ParameterSet _vuv_params
SBNDOpT0Finder::SBNDOpT0Finder ( SBNDOpT0Finder const &  )
delete
SBNDOpT0Finder::SBNDOpT0Finder ( SBNDOpT0Finder &&  )
delete

Member Function Documentation

bool SBNDOpT0Finder::ConstructLightClusters ( art::Event &  e,
unsigned int  tpc 
)
private

Constructs all the LightClusters (TPC Objects) in a specified TPC.

Definition at line 402 of file SBNDOpT0Finder_module.cc.

402  {
403  // One slice is one QCluster_t.
404  // Start from a slice, get all the PFParticles, from there get all the spacepoints, from
405  // there get all the hits on the collection plane.
406  // Use the charge on the collection plane to estimate the light, and the 3D spacepoint
407  // position for the 3D location.
408 
409  _light_cluster_v.clear();
410 
411  ::art::Handle<std::vector<recob::Slice>> slice_h;
412  e.getByLabel(_slice_producer, slice_h);
413  if(!slice_h.isValid() || slice_h->empty()) {
414  mf::LogWarning("SBNDOpT0Finder") << "Don't have good Slices." << std::endl;
415  return false;
416  }
417 
418  ::art::Handle<std::vector<recob::PFParticle>> pfp_h;
419  e.getByLabel(_slice_producer, pfp_h);
420  if(!pfp_h.isValid() || pfp_h->empty()) {
421  mf::LogWarning("SBNDOpT0Finder") << "Don't have good PFParticle." << std::endl;
422  return false;
423  }
424 
425  ::art::Handle<std::vector<recob::SpacePoint>> spacepoint_h;
426  e.getByLabel(_slice_producer, spacepoint_h);
427  if(!spacepoint_h.isValid() || spacepoint_h->empty()) {
428  mf::LogWarning("SBNDOpT0Finder") << "Don't have good SpacePoint." << std::endl;
429  return false;
430  }
431 
432  // Construct the vector of Slices
433  std::vector<art::Ptr<recob::Slice>> slice_v;
434  art::fill_ptr_vector(slice_v, slice_h);
435 
436  // Get the associations between slice->pfp->spacepoint->hit
437  art::FindManyP<recob::PFParticle> slice_to_pfps (slice_h, e, _slice_producer);
438  art::FindManyP<recob::SpacePoint> pfp_to_spacepoints (pfp_h, e, _slice_producer);
439  art::FindManyP<recob::Hit> spacepoint_to_hits (spacepoint_h, e, _slice_producer);
440 
441  // Loop over the Slices
442  for (size_t n_slice = 0; n_slice < slice_h->size(); n_slice++) {
443 
444  flashmatch::QCluster_t light_cluster;
445 
446  _dep_slice.clear();
447  _dep_x.clear();
448  _dep_y.clear();
449  _dep_z.clear();
450  _dep_charge.clear();
451  _dep_n_photons.clear();
452 
453  // Get the associated PFParticles
454  std::vector<art::Ptr<recob::PFParticle>> pfp_v = slice_to_pfps.at(n_slice);
455 
456  for (size_t n_pfp = 0; n_pfp < pfp_v.size(); n_pfp++) {
457 
458  auto pfp = pfp_v[n_pfp];
459 
460  // Get the associated SpacePoints
461  std::vector<art::Ptr<recob::SpacePoint>> spacepoint_v = pfp_to_spacepoints.at(pfp.key());
462 
463  for (size_t n_spacepoint = 0; n_spacepoint < spacepoint_v.size(); n_spacepoint++) {
464 
465  auto spacepoint = spacepoint_v[n_spacepoint];
466 
467  // Get the associated hits
468  std::vector<art::Ptr<recob::Hit>> hit_v = spacepoint_to_hits.at(spacepoint.key());
469 
470  for (size_t n_hit = 0; n_hit < hit_v.size(); n_hit++) {
471 
472  auto hit = hit_v[n_hit];
473 
474  // Only select hits from the collection plane
475  if (hit->View() != geo::kZ) {
476  continue;
477  }
478 
479  // Only use hits (and so spacepoints) that are in the specified TPC
480  if (hit->WireID().TPC != tpc) {
481  continue;
482  }
483 
484  const auto &position(spacepoint->XYZ());
485  const auto charge(hit->Integral());
486 
487  // Emplace this point with charge to the light cluster
488  light_cluster.emplace_back(position[0],
489  position[1],
490  position[2],
491  GetNPhotons(charge, pfp));
492 
493  // Also save the quantites for the output tree
494  _dep_slice.push_back(_light_cluster_v.size());
495  _dep_x.push_back(position[0]);
496  _dep_y.push_back(position[1]);
497  _dep_z.push_back(position[2]);
498  _dep_charge.push_back(charge);
499  _dep_n_photons.push_back(GetNPhotons(charge, pfp));
500  }
501  } // End loop over Spacepoints
502  } // End loop over PFParticle
503 
504  _tree1->Fill();
505 
506  // Don't include clusters with zero points
507  if (!light_cluster.size()) {
508  continue;
509  }
510 
511  // Save the light cluster, and remember the correspondance from index to slice
512  _clusterid_to_slice[_light_cluster_v.size()] = slice_v.at(n_slice);
513  _light_cluster_v.emplace_back(light_cluster);
514 
515  } // End loop over Slices
516 
517  return true;
518 }
std::vector< float > _dep_z
std::vector< float > _dep_n_photons
std::vector< int > _dep_slice
float GetNPhotons(const float charge, const art::Ptr< recob::PFParticle > &pfp)
Returns the number of photons given charge and PFParticle.
std::vector< float > _dep_charge
std::map< int, art::Ptr< recob::Slice > > _clusterid_to_slice
Planes which measure Z direction.
Definition: geo_types.h:132
process_name hit
Definition: cheaterreco.fcl:51
Collection of charge deposition 3D point (cluster)
std::string _slice_producer
The Slice producer (to be set)
std::vector< float > _dep_x
std::vector< float > _dep_y
TTree * _tree1
Will contain map flash id -&gt; OpFlash.
do i e
std::vector< flashmatch::QCluster_t > _light_cluster_v
Vector that contains all the TPC objects.
void SBNDOpT0Finder::DoMatch ( art::Event &  e,
int  tpc,
std::unique_ptr< std::vector< anab::T0 >> &  t0_v,
std::unique_ptr< art::Assns< recob::Slice, anab::T0 >> &  slice_t0_assn_v,
std::unique_ptr< art::Assns< recob::OpFlash, anab::T0 >> &  flash_t0_assn_v 
)
private

Performs the matching in a specified tpc.

Definition at line 250 of file SBNDOpT0Finder_module.cc.

254  {
255 
256  _flashid_to_opflash.clear();
257  _clusterid_to_slice.clear();
258 
259  auto const & flash_h = e.getValidHandle<std::vector<recob::OpFlash>>(_opflash_producer_v[tpc]);
260  if(!flash_h.isValid() || flash_h->empty()) {
261  mf::LogInfo("SBNDOpT0Finder") << "Don't have good flashes from producer "
262  << _opflash_producer_v[tpc] << std::endl;
263  return;
264  }
265 
266  // Construct the vector of OpFlashes
267  std::vector<art::Ptr<recob::OpFlash>> flash_v;
268  art::fill_ptr_vector(flash_v, flash_h);
269 
270  ::art::ServiceHandle<geo::Geometry> geo;
271 
272  int n_flashes = 0;
273  std::vector<::flashmatch::Flash_t> all_flashes;
274 
275  for (size_t n = 0; n < flash_v.size(); n++) {
276 
277  auto const& flash = *flash_v[n];
278 
279  mf::LogDebug("SBNDOpT0Finder") << "Flash time from " << _opflash_producer_v[tpc] << ": " << flash.Time() << std::endl;
280 
281  if(flash.Time() < _flash_trange_start || _flash_trange_end < flash.Time()) {
282  continue;
283  }
284 
285  _flashid_to_opflash[n_flashes] = flash_v[n];
286 
287  n_flashes++;
288 
289  // Construct a Flash_t
291  f.x = f.x_err = 0;
292  f.pe_v.resize(geo->NOpDets());
293  f.pe_err_v.resize(geo->NOpDets());
294  for (unsigned int op_ch = 0; op_ch < f.pe_v.size(); op_ch++) {
295  unsigned int opdet = geo->OpDetFromOpChannel(op_ch);
296  if (std::find(_opch_to_use.begin(), _opch_to_use.end(), op_ch) == _opch_to_use.end()) {
297  f.pe_v[opdet] = 0;
298  f.pe_err_v[opdet] = 0;
299  } else {
300  f.pe_v[opdet] = flash.PE(op_ch);
301  f.pe_err_v[opdet] = sqrt(flash.PE(op_ch));
302  }
303  }
304  f.y = flash.YCenter();
305  f.z = flash.ZCenter();
306  f.y_err = flash.YWidth();
307  f.z_err = flash.ZWidth();
308  f.time = flash.Time();
309  f.idx = n_flashes-1;
310  all_flashes.resize(n_flashes);
311  all_flashes[n_flashes-1] = f;
312 
313  } // flash loop
314 
315  // Don't waste time if there are no flashes
316  if (n_flashes == 0) {
317  mf::LogInfo("SBNDOpT0Finder") << "Zero good flashes in this event." << std::endl;
318  return;
319  }
320 
321  // Get all the ligh clusters
322  // auto light_cluster_v = GetLighClusters(e);
323  if (!ConstructLightClusters(e, tpc)) {
324  mf::LogInfo("SBNDOpT0Finder") << "Cannot construct Light Clusters." << std::endl;
325  return;
326  }
327 
328  // Don't waste time if there are no clusters
329  if (!_light_cluster_v.size()) {
330  mf::LogInfo("SBNDOpT0Finder") << "No slices to work with in TPC " << tpc << "." << std::endl;
331  return;
332  }
333 
334  // Emplace flashes to Flash Matching Manager
335  for (auto f : all_flashes) {
336  _mgr.Emplace(std::move(f));
337  }
338 
339  // Emplace clusters to Flash Matching Manager
340  for (auto lc : _light_cluster_v) {
341  _mgr.Emplace(std::move(lc));
342  }
343 
344  // Run the matching
345  _result_v = _mgr.Match();
346 
347  // Loop over the matching results
348  for(_matchid = 0; _matchid < (int)(_result_v.size()); ++_matchid) {
349 
350  auto const& match = _result_v[_matchid];
351 
352  _tpcid = match.tpc_id;
353  _flashid = match.flash_id;
354  _score = match.score;
355  _qll_xmin = match.tpc_point.x;
356 
357  mf::LogInfo("SBNDOpT0Finder") << "Matched TPC object " << _tpcid
358  << " with flash number " << _flashid
359  << " in TPC " << tpc
360  << " -> score: " << _score
361  << ", qll xmin: " << _qll_xmin << std::endl;
362 
363  // Get the minimum x position of the TPC Object
364  _tpc_xmin = 1.e4;
365  for(auto const& pt : _mgr.QClusterArray()[_tpcid]) {
366  if(pt.x < _tpc_xmin) _tpc_xmin = pt.x;
367  }
368 
369  // Get the matched flash time, the t0
370  auto const& flash = _mgr.FlashArray()[_flashid];
371  _t0 = flash.time;
372 
373  // Save the reconstructed flash and hypothesis flash PE spectrum
374  if(_hypo_spec.size() != match.hypothesis.size()) {
375  throw cet::exception("SBNDOpT0Finder") << "Hypothesis size mismatch!";
376  }
377  for(size_t pmt=0; pmt<_hypo_spec.size(); ++pmt) _hypo_spec[pmt] = match.hypothesis[pmt];
378  for(size_t pmt=0; pmt<_hypo_spec.size(); ++pmt) _flash_spec[pmt] = flash.pe_v[pmt];
379 
380  // Also save the total number of photoelectrons
381  _flash_pe = 0.;
382  _hypo_pe = 0.;
383  for(auto const& v : _hypo_spec) _hypo_pe += v;
384  for(auto const& v : _flash_spec) _flash_pe += v;
385 
386  _tree2->Fill();
387 
388  // Construct the anab::T0 dataproduc to put in the Event
389  auto t0 = anab::T0(_t0, // "Time": The recontructed flash time, or t0
390  _flash_pe, // "TriggerType": placing the reconstructed total PE instead
391  _tpcid, // "TriggerBits": placing the tpc id instead
392  _flashid, // "ID": placing the flash id instead
393  _score); // "TriggerConfidence": Matching score
394 
395  t0_v->push_back(t0);
396  util::CreateAssn(*this, e, *t0_v, _clusterid_to_slice[_tpcid], *slice_t0_assn_v);
397  util::CreateAssn(*this, e, *t0_v, _flashid_to_opflash[_flashid], *flash_t0_assn_v);
398  }
399 
400 }
double _flash_trange_end
The time stop from where to stop including flashes (to be set)
std::vector< flashmatch::FlashMatch_t > Match()
std::map< int, art::Ptr< recob::OpFlash > > _flashid_to_opflash
Will contain map tpc object id -&gt; Slice.
std::vector< std::string > _opflash_producer_v
The OpFlash producers (to be set)
std::map< int, art::Ptr< recob::Slice > > _clusterid_to_slice
std::vector< double > _flash_spec
::flashmatch::FlashMatchManager _mgr
The flash matching manager.
double _flash_trange_start
The time start from where to include flashes (to be set)
Struct to represent an optical flash.
const QClusterArray_t & QClusterArray() const
Access to an input: TPC objects in the form of QClusterArray_t.
std::vector< double > pe_v
PE distribution over photo-detectors.
std::vector< double > _hypo_spec
bool CreateAssn(art::Event &evt, std::vector< T > const &a, art::Ptr< U > const &b, art::Assns< U, T > &assn, std::string a_instance, size_t index=UINT_MAX)
Creates a single one-to-one association.
const FlashArray_t & FlashArray() const
Access to an input: PMT objects in the form of FlashArray_t.
bool ConstructLightClusters(art::Event &e, unsigned int tpc)
Constructs all the LightClusters (TPC Objects) in a specified TPC.
std::vector< flashmatch::FlashMatch_t > _result_v
Matching result will be stored here.
std::vector< int > _opch_to_use
List of opch to use (will be infered from _photo_detectors)
do i e
std::vector< double > pe_err_v
PE value error.
double z
Flash position.
void Emplace(flashmatch::QCluster_t &&obj)
Emplacer of a TPC object (hidden from ROOT5 CINT)
double time
Flash timing, a candidate T0.
double z_err
Flash position error.
process_name opdaq physics producers generator physics producers generator physics producers generator physics producers generator physics producers generator physics producers generator physics producers generator physics producers generator T0
Definition: gen_protons.fcl:45
std::vector< flashmatch::QCluster_t > _light_cluster_v
Vector that contains all the TPC objects.
float SBNDOpT0Finder::GetNPhotons ( const float  charge,
const art::Ptr< recob::PFParticle > &  pfp 
)
private

Returns the number of photons given charge and PFParticle.

Definition at line 520 of file SBNDOpT0Finder_module.cc.

521  {
524 }
float _charge_to_n_photons_track
The conversion factor betweeen hit integral and photons (to be set)
static bool IsTrack(const art::Ptr< recob::PFParticle > particle)
Determine whether a particle has been reconstructed as track-like.
float _charge_to_n_photons_shower
The conversion factor betweeen hit integral and photons (to be set)
std::vector< int > SBNDOpT0Finder::GetUncoatedPTMList ( std::vector< int >  ch_to_use)
private

Returns a list of uncoated PMTs that are a subset of those in ch_to_use.

Definition at line 539 of file SBNDOpT0Finder_module.cc.

539  {
540  std::vector<int> out_v;
541 
542  for (auto ch : ch_to_use) {
543  if (_pds_map.isPDType(ch, "pmt_uncoated")) {
544  out_v.push_back(ch);
545  }
546  }
547 
548  return out_v;
549 }
bool isPDType(size_t ch, std::string pdname) const override
opdet::sbndPDMapAlg _pds_map
map for photon detector types
SBNDOpT0Finder& SBNDOpT0Finder::operator= ( SBNDOpT0Finder const &  )
delete
SBNDOpT0Finder& SBNDOpT0Finder::operator= ( SBNDOpT0Finder &&  )
delete
std::vector< int > SBNDOpT0Finder::PDNamesToList ( std::vector< std::string >  pd_names)
private

Convert from a list of PDS names to a list of op channels.

Definition at line 526 of file SBNDOpT0Finder_module.cc.

526  {
527 
528  std::vector<int> out_ch_v;
529 
530  for (auto name : pd_names) {
531  auto ch_v = _pds_map.getChannelsOfType(name);
532  out_ch_v.insert(out_ch_v.end(), ch_v.begin(), ch_v.end());
533  }
534 
535  return out_ch_v;
536 
537 }
opdet::sbndPDMapAlg _pds_map
map for photon detector types
std::vector< int > getChannelsOfType(std::string pdname) const
then echo fcl name
void SBNDOpT0Finder::produce ( art::Event &  e)
override

Definition at line 211 of file SBNDOpT0Finder_module.cc.

212 {
213  std::unique_ptr<std::vector<anab::T0>> t0_v (new std::vector<anab::T0>);
214  std::unique_ptr< art::Assns<recob::Slice, anab::T0>> slice_t0_assn_v (new art::Assns<recob::Slice, anab::T0>);
215  std::unique_ptr< art::Assns<recob::OpFlash, anab::T0>> flash_t0_assn_v (new art::Assns<recob::OpFlash, anab::T0>);
216 
217  // _mgr.PrintConfig();
218 
219  _run = e.id().run();
220  _subrun = e.id().subRun();
221  _event = e.id().event();
222 
223  // Loop over the specified TPCs
224  for (auto tpc : _tpc_v) {
225 
226  mf::LogInfo("SBNDOpT0Finder") << "Performing matching in TPC " << tpc << std::endl;
227 
228  // Reset the manager and the result vector
229  _mgr.Reset();
230  _result_v.clear();
231  _tpc = tpc;
232 
233  // Tell the manager what TPC and cryostat we are going to be doing
234  // the matching in. For SBND, the cryostat is always zero.
235  _mgr.SetTPCCryo(tpc, 0);
236 
237  // Perform the matching in the specified TPC
238  DoMatch(e, tpc, t0_v, slice_t0_assn_v, flash_t0_assn_v);
239 
240  }
241 
242  // Finally, place the anab::T0 vector and the associations in the Event
243  e.put(std::move(t0_v));
244  e.put(std::move(slice_t0_assn_v));
245  e.put(std::move(flash_t0_assn_v));
246 
247  return;
248 }
void SetTPCCryo(int tpc, int _cryo)
Sets the TPC and Cryo numbers.
::flashmatch::FlashMatchManager _mgr
The flash matching manager.
std::vector< unsigned int > _tpc_v
TPC number per OpFlash producer (to be set)
void Reset()
Clears locally kept TPC object (QClusterArray_t) and flash (FlashArray_t), both provided by a user...
std::vector< flashmatch::FlashMatch_t > _result_v
Matching result will be stored here.
void DoMatch(art::Event &e, int tpc, std::unique_ptr< std::vector< anab::T0 >> &t0_v, std::unique_ptr< art::Assns< recob::Slice, anab::T0 >> &slice_t0_assn_v, std::unique_ptr< art::Assns< recob::OpFlash, anab::T0 >> &flash_t0_assn_v)
Performs the matching in a specified tpc.
do i e

Member Data Documentation

float SBNDOpT0Finder::_charge_to_n_photons_shower
private

The conversion factor betweeen hit integral and photons (to be set)

Definition at line 107 of file SBNDOpT0Finder_module.cc.

float SBNDOpT0Finder::_charge_to_n_photons_track
private

The conversion factor betweeen hit integral and photons (to be set)

Definition at line 106 of file SBNDOpT0Finder_module.cc.

std::map<int, art::Ptr<recob::Slice> > SBNDOpT0Finder::_clusterid_to_slice
private

Definition at line 118 of file SBNDOpT0Finder_module.cc.

std::vector<float> SBNDOpT0Finder::_dep_charge
private

Definition at line 132 of file SBNDOpT0Finder_module.cc.

std::vector<float> SBNDOpT0Finder::_dep_n_photons
private

Definition at line 132 of file SBNDOpT0Finder_module.cc.

std::vector<int> SBNDOpT0Finder::_dep_slice
private

Definition at line 133 of file SBNDOpT0Finder_module.cc.

std::vector<float> SBNDOpT0Finder::_dep_x
private

Definition at line 132 of file SBNDOpT0Finder_module.cc.

std::vector<float> SBNDOpT0Finder::_dep_y
private

Definition at line 132 of file SBNDOpT0Finder_module.cc.

std::vector<float> SBNDOpT0Finder::_dep_z
private

Definition at line 132 of file SBNDOpT0Finder_module.cc.

int SBNDOpT0Finder::_event
private

Definition at line 122 of file SBNDOpT0Finder_module.cc.

double SBNDOpT0Finder::_flash_pe
private

Definition at line 127 of file SBNDOpT0Finder_module.cc.

std::vector<double> SBNDOpT0Finder::_flash_spec
private

Definition at line 128 of file SBNDOpT0Finder_module.cc.

double SBNDOpT0Finder::_flash_trange_end
private

The time stop from where to stop including flashes (to be set)

Definition at line 104 of file SBNDOpT0Finder_module.cc.

double SBNDOpT0Finder::_flash_trange_start
private

The time start from where to include flashes (to be set)

Definition at line 103 of file SBNDOpT0Finder_module.cc.

int SBNDOpT0Finder::_flashid
private

Definition at line 124 of file SBNDOpT0Finder_module.cc.

std::map<int, art::Ptr<recob::OpFlash> > SBNDOpT0Finder::_flashid_to_opflash
private

Will contain map tpc object id -> Slice.

Definition at line 119 of file SBNDOpT0Finder_module.cc.

double SBNDOpT0Finder::_hypo_pe
private

Definition at line 127 of file SBNDOpT0Finder_module.cc.

std::vector<double> SBNDOpT0Finder::_hypo_spec
private

Definition at line 129 of file SBNDOpT0Finder_module.cc.

std::vector<flashmatch::QCluster_t> SBNDOpT0Finder::_light_cluster_v
private

Vector that contains all the TPC objects.

Definition at line 116 of file SBNDOpT0Finder_module.cc.

int SBNDOpT0Finder::_matchid
private

Definition at line 124 of file SBNDOpT0Finder_module.cc.

::flashmatch::FlashMatchManager SBNDOpT0Finder::_mgr
private

The flash matching manager.

Definition at line 96 of file SBNDOpT0Finder_module.cc.

std::vector<int> SBNDOpT0Finder::_opch_to_use
private

List of opch to use (will be infered from _photo_detectors)

Definition at line 110 of file SBNDOpT0Finder_module.cc.

std::vector<std::string> SBNDOpT0Finder::_opflash_producer_v
private

The OpFlash producers (to be set)

Definition at line 99 of file SBNDOpT0Finder_module.cc.

opdet::sbndPDMapAlg SBNDOpT0Finder::_pds_map
private

map for photon detector types

Definition at line 113 of file SBNDOpT0Finder_module.cc.

std::vector<std::string> SBNDOpT0Finder::_photo_detectors
private

The photodetector to use (to be set)

Definition at line 109 of file SBNDOpT0Finder_module.cc.

double SBNDOpT0Finder::_qll_xmin
private

Definition at line 126 of file SBNDOpT0Finder_module.cc.

std::vector<flashmatch::FlashMatch_t> SBNDOpT0Finder::_result_v
private

Matching result will be stored here.

Definition at line 97 of file SBNDOpT0Finder_module.cc.

int SBNDOpT0Finder::_run
private

Definition at line 122 of file SBNDOpT0Finder_module.cc.

double SBNDOpT0Finder::_score
private

Definition at line 125 of file SBNDOpT0Finder_module.cc.

std::unique_ptr<SemiAnalyticalModel> SBNDOpT0Finder::_semi_model
private

Definition at line 92 of file SBNDOpT0Finder_module.cc.

std::string SBNDOpT0Finder::_slice_producer
private

The Slice producer (to be set)

Definition at line 101 of file SBNDOpT0Finder_module.cc.

int SBNDOpT0Finder::_subrun
private

Definition at line 122 of file SBNDOpT0Finder_module.cc.

double SBNDOpT0Finder::_t0
private

Definition at line 125 of file SBNDOpT0Finder_module.cc.

int SBNDOpT0Finder::_tpc
private

Definition at line 123 of file SBNDOpT0Finder_module.cc.

std::vector<unsigned int> SBNDOpT0Finder::_tpc_v
private

TPC number per OpFlash producer (to be set)

Definition at line 100 of file SBNDOpT0Finder_module.cc.

double SBNDOpT0Finder::_tpc_xmin
private

Definition at line 126 of file SBNDOpT0Finder_module.cc.

int SBNDOpT0Finder::_tpcid
private

Definition at line 124 of file SBNDOpT0Finder_module.cc.

TTree* SBNDOpT0Finder::_tree1
private

Will contain map flash id -> OpFlash.

Definition at line 121 of file SBNDOpT0Finder_module.cc.

TTree* SBNDOpT0Finder::_tree2
private

Definition at line 131 of file SBNDOpT0Finder_module.cc.

std::vector<int> SBNDOpT0Finder::_uncoated_pmts
private

List of uncoated opch to use (will be infered from _opch_to_use)

Definition at line 111 of file SBNDOpT0Finder_module.cc.

fhicl::ParameterSet SBNDOpT0Finder::_vis_params
private

Definition at line 94 of file SBNDOpT0Finder_module.cc.

fhicl::ParameterSet SBNDOpT0Finder::_vuv_params
private

Definition at line 93 of file SBNDOpT0Finder_module.cc.


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