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
opdet::SBNDFlashFinder Class Reference
Inheritance diagram for opdet::SBNDFlashFinder:

Public Member Functions

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

Private Member Functions

::lightana::LiteOpHitArray_t GetAssociatedLiteHits (::lightana::LiteOpFlash_t lite_flash,::lightana::LiteOpHitArray_t lite_hits_v)
 

Private Attributes

::lightana::FlashFinderManager _mgr
 
::lightana::PECalib _pecalib
 
std::vector< std::string > _hit_producers
 
std::string _ophit_input_time
 
bool _use_t0tool
 
bool _correct_light_propagation
 
std::unique_ptr
< lightana::FlashGeoBase
_flashgeo
 
std::unique_ptr
< lightana::FlashT0Base
_flasht0calculator
 
std::unique_ptr
< lightana::DriftEstimatorBase
_driftestimator
 

Detailed Description

Definition at line 37 of file SBNDFlashFinder_module.cc.

Constructor & Destructor Documentation

opdet::SBNDFlashFinder::SBNDFlashFinder ( fhicl::ParameterSet const &  p)
explicit

Definition at line 75 of file SBNDFlashFinder_module.cc.

76  : EDProducer{p}
77  // Initialize member data here.
78  {
79  _hit_producers = p.get<std::vector<std::string>>("OpHitProducers");
80 
81  auto const flash_algo = p.get<std::string>("FlashFinderAlgo");
82  auto const flash_pset = p.get<lightana::Config_t>("AlgoConfig");
83  auto algo_ptr = ::lightana::FlashAlgoFactory::get().create(flash_algo,flash_algo);
84  algo_ptr->Configure(flash_pset);
85  _mgr.SetFlashAlgo(algo_ptr);
86  _pecalib.Configure(p.get<lightana::Config_t>("PECalib"));
87  _ophit_input_time = p.get<std::string>("OpHitInputTime", "PeakTime");
88  _use_t0tool = p.get<bool>("UseT0Tool", false);
89  _correct_light_propagation = p.get<bool>("CorrectLightPropagation", false);
90 
91  auto const flashgeo_pset = p.get<lightana::Config_t>("FlashGeoConfig");
92  _flashgeo = art::make_tool<lightana::FlashGeoBase>(flashgeo_pset);
93 
94  if(_use_t0tool){
95  auto const flasht0_pset = p.get<lightana::Config_t>("FlashT0Config");
96  _flasht0calculator = art::make_tool<lightana::FlashT0Base>(flasht0_pset);
97  }
98 
99  if(_correct_light_propagation){
100  auto const driftestimator_pset = p.get<lightana::Config_t>("DriftEstimatorConfig");
101  _driftestimator = art::make_tool<lightana::DriftEstimatorBase>(driftestimator_pset);
102  }
103 
104  produces< std::vector<recob::OpFlash> >();
105  produces< art::Assns <recob::OpHit, recob::OpFlash> >();
106  }
std::unique_ptr< lightana::DriftEstimatorBase > _driftestimator
::lightana::PECalib _pecalib
pdgs p
Definition: selectors.fcl:22
std::unique_ptr< lightana::FlashGeoBase > _flashgeo
::lightana::FlashFinderManager _mgr
std::vector< std::string > _hit_producers
std::unique_ptr< lightana::FlashT0Base > _flasht0calculator
static FlashAlgoFactory & get()
Static sharable instance getter.
FlashAlgoBase * create(const std::string name, const std::string instance_name)
Factory creation method (should be called by clients, possibly you!)
virtual void Configure(const Config_t &p)=0
opdet::SBNDFlashFinder::SBNDFlashFinder ( SBNDFlashFinder const &  )
delete
opdet::SBNDFlashFinder::SBNDFlashFinder ( SBNDFlashFinder &&  )
delete

Member Function Documentation

lightana::LiteOpHitArray_t opdet::SBNDFlashFinder::GetAssociatedLiteHits ( ::lightana::LiteOpFlash_t  lite_flash,
::lightana::LiteOpHitArray_t  lite_hits_v 
)
private

Definition at line 202 of file SBNDFlashFinder_module.cc.

204  {
205 
206  ::lightana::LiteOpHitArray_t flash_hits_v;
207 
208  for(auto const& hitidx : lite_flash.asshit_idx) {
209  flash_hits_v.emplace_back(std::move(lite_hits_v.at(hitidx)));
210  }
211 
212  return flash_hits_v;
213  }
std::vector< lightana::LiteOpHit_t > LiteOpHitArray_t
SBNDFlashFinder& opdet::SBNDFlashFinder::operator= ( SBNDFlashFinder const &  )
delete
SBNDFlashFinder& opdet::SBNDFlashFinder::operator= ( SBNDFlashFinder &&  )
delete
void opdet::SBNDFlashFinder::produce ( art::Event &  e)
override

Definition at line 108 of file SBNDFlashFinder_module.cc.

109  {
110 
111  // produce OpFlash data-product to be filled within module
112  std::unique_ptr< std::vector<recob::OpFlash> > opflashes(new std::vector<recob::OpFlash>);
113  std::unique_ptr< art::Assns <recob::OpHit, recob::OpFlash> > flash2hit_assn_v
114  (new art::Assns<recob::OpHit, recob::OpFlash>);
115 
116  std::vector<art::Ptr<recob::OpHit>> ophit_v;
117 
119  double trigger_time=1.1e20;
120 
121  for (auto producer : _hit_producers) {
122 
123  // load OpHits previously created
124  art::Handle<std::vector<recob::OpHit> > ophit_h;
125  e.getByLabel(producer, ophit_h);
126 
127  // make sure hits look good
128  if(!ophit_h.isValid()) {
129  std::cerr << "\033[93m[ERROR]\033[00m ... could not locate OpHit!" << std::endl;
130  throw std::exception();
131  }
132 
133  std::vector<art::Ptr<recob::OpHit>> temp_v;
134  art::fill_ptr_vector(temp_v, ophit_h);
135  ophit_v.insert(ophit_v.end(), temp_v.begin(), temp_v.end());
136  }
137 
138  for(auto const oph : ophit_v) {
140  if(trigger_time > 1.e20) trigger_time = oph->PeakTimeAbs() - oph->PeakTime();
141 
142  if(_ophit_input_time=="RiseTime") loph.peak_time = oph->StartTime()+oph->RiseTime();
143  else if(_ophit_input_time=="StartTime") loph.peak_time = oph->StartTime();
144  else loph.peak_time = oph->PeakTime();
145 
146  size_t opdet = ::lightana::OpDetFromOpChannel(oph->OpChannel());
147  loph.pe = _pecalib.Calibrate(opdet,oph->Area());
148  loph.channel = oph->OpChannel();
149  ophits.emplace_back(std::move(loph));
150  }
151 
152  auto const flash_v = _mgr.RecoFlash(ophits);
153 
154  for(const auto& lflash : flash_v) {
155 
156  // Get Flash Barycenter
157  double Ycenter, Zcenter, Ywidth, Zwidth;
158  _flashgeo->GetFlashLocation(lflash.channel_pe, Ycenter, Zcenter, Ywidth, Zwidth);
159 
160  // Get flasht0
161  double flasht0 = lflash.time;
162 
163  // Refine t0 calculation
164  if(_use_t0tool)
165  flasht0 = _flasht0calculator->GetFlashT0(lflash.time, GetAssociatedLiteHits(lflash, ophits));
166 
167  // Estimate drift location of the interaction and
168  // make t0 unbias (lght propagation time correction)
170 
171  double drift_distance = _driftestimator->GetDriftPosition( lflash.channel_pe );
172  double propagation_time = _driftestimator->GetPropagationTime( drift_distance );
173  flasht0 = flasht0-propagation_time * 1e-3;
174 
175  drift_distance = (lflash.tpc==0 ? -drift_distance : drift_distance);
176 
177  recob::OpFlash flash(flasht0, lflash.time_err, trigger_time + flasht0,
178  (trigger_time + flasht0) / 1600., lflash.channel_pe,
179  0, 0, 1, // this are just default values
180  drift_distance, -1, Ycenter, Ywidth, Zcenter, Zwidth);
181  opflashes->emplace_back(std::move(flash));
182  }
183  else{
184  recob::OpFlash flash(flasht0, lflash.time_err, trigger_time + flasht0,
185  (trigger_time + flasht0) / 1600., lflash.channel_pe,
186  0, 0, 1, // this are just default values
187  100., -1., Ycenter, Ywidth, Zcenter, Zwidth);
188  opflashes->emplace_back(std::move(flash));
189  }
190 
191  // Create OpHit association
192  for(auto const& hitidx : lflash.asshit_idx) {
193  const art::Ptr<recob::OpHit> hit_ptr(ophit_v.at(hitidx));
194  util::CreateAssn(*this, e, *opflashes, hit_ptr, *flash2hit_assn_v);
195  }
196  }
197 
198  e.put(std::move(opflashes));
199  e.put(std::move(flash2hit_assn_v));
200  }
std::unique_ptr< lightana::DriftEstimatorBase > _driftestimator
BEGIN_PROLOG could also be cerr
::lightana::PECalib _pecalib
timescale_traits< TriggerTimeCategory >::time_point_t trigger_time
A point in time on the trigger time scale.
std::vector< lightana::LiteOpHit_t > LiteOpHitArray_t
double Calibrate(const size_t opdet, const double area) const
std::unique_ptr< lightana::FlashGeoBase > _flashgeo
::lightana::FlashFinderManager _mgr
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.
std::vector< std::string > _hit_producers
std::unique_ptr< lightana::FlashT0Base > _flasht0calculator
LiteOpFlashArray_t RecoFlash(const LiteOpHitArray_t &ophits) const
do i e
::lightana::LiteOpHitArray_t GetAssociatedLiteHits(::lightana::LiteOpFlash_t lite_flash,::lightana::LiteOpHitArray_t lite_hits_v)

Member Data Documentation

bool opdet::SBNDFlashFinder::_correct_light_propagation
private

Definition at line 59 of file SBNDFlashFinder_module.cc.

std::unique_ptr<lightana::DriftEstimatorBase> opdet::SBNDFlashFinder::_driftestimator
private

Definition at line 68 of file SBNDFlashFinder_module.cc.

std::unique_ptr<lightana::FlashGeoBase> opdet::SBNDFlashFinder::_flashgeo
private

Definition at line 62 of file SBNDFlashFinder_module.cc.

std::unique_ptr<lightana::FlashT0Base> opdet::SBNDFlashFinder::_flasht0calculator
private

Definition at line 65 of file SBNDFlashFinder_module.cc.

std::vector<std::string> opdet::SBNDFlashFinder::_hit_producers
private

Definition at line 56 of file SBNDFlashFinder_module.cc.

::lightana::FlashFinderManager opdet::SBNDFlashFinder::_mgr
private

Definition at line 54 of file SBNDFlashFinder_module.cc.

std::string opdet::SBNDFlashFinder::_ophit_input_time
private

Definition at line 57 of file SBNDFlashFinder_module.cc.

::lightana::PECalib opdet::SBNDFlashFinder::_pecalib
private

Definition at line 55 of file SBNDFlashFinder_module.cc.

bool opdet::SBNDFlashFinder::_use_t0tool
private

Definition at line 58 of file SBNDFlashFinder_module.cc.


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