LArSoft  v07_13_02
Liquid Argon Software toolkit - http://larsoft.org/
BeamFlashTrackMatchTaggerAlg.cxx
Go to the documentation of this file.
1 
14 #include <limits>
15 
17  : COSMIC_TYPE_FLASHMATCH(anab::CosmicTagID_t::kFlash_BeamIncompatible),
18  COSMIC_TYPE_OUTSIDEDRIFT(anab::CosmicTagID_t::kOutsideDrift_Partial),
19  DEBUG_FLAG(p.get<bool>("RunDebugMode",false))
20 {
21  this->reconfigure(p);
22 }
23 
25  fMinTrackLength = p.get<float>("MinTrackLength");
26  fMinOpHitPE = p.get<float>("MinOpHitPE",0.1);
27  fMIPdQdx = p.get<float>("MIPdQdx",2.1);
28  fOpDetSaturation = p.get<float>("OpDetSaturation",200.);
29 
30  fSingleChannelCut = p.get<float>("SingleChannelCut");
31  fCumulativeChannelThreshold = p.get<float>("CumulativeChannelThreshold");
32  fCumulativeChannelCut = p.get<unsigned int>("CumulativeChannelCut");
33  fIntegralCut = p.get<float>("IntegralCut");
34 
35  fMakeOutsideDriftTags = p.get<bool>("MakeOutsideDriftTags",false);
36  fNormalizeHypothesisToFlash = p.get<bool>("NormalizeHypothesisToFlash");
37 }
38 
40  TH1F* hist_flash, TH1F* hist_hyp){
41  cTree = tree;
42 
43  cOpDetHist_flash = hist_flash;
44  cOpDetHist_flash->SetNameTitle("opdet_hist_flash","Optical Detector Occupancy, Flash");
45 
46  cOpDetHist_hyp = hist_hyp;
47  cOpDetHist_hyp->SetNameTitle("opdet_hist_hyp","Optical Detector Occupancy, Hypothesis");
48 
50  cTree->Branch("opdet_hyp",&cOpDetVector_hyp);
51  cTree->Branch("opdet_flash",&cOpDetVector_flash);
52  cTree->Branch("opdet_hist_flash","TH1F",cOpDetHist_flash);
53  cTree->Branch("opdet_hist_hyp","TH1F",cOpDetHist_hyp);
54 }
55 
56 void cosmic::BeamFlashTrackMatchTaggerAlg::RunCompatibilityCheck(std::vector<recob::OpFlash> const& flashVector,
57  std::vector<recob::Track> const& trackVector,
58  std::vector<anab::CosmicTag>& cosmicTagVector,
59  std::vector<size_t>& assnTrackTagVector,
60  Providers_t providers,
62  opdet::OpDigiProperties const& opdigip){
63 
64  auto const& geom = *(providers.get<geo::GeometryCore>());
65 
66  std::vector< const recob::OpFlash* > flashesOnBeamTime;
67  for(auto const& flash : flashVector){
68  if(!flash.OnBeamTime()) continue;
69  flashesOnBeamTime.push_back(&flash);
70  }
71 
72  //make sure this association vector is initialized properly
73  assnTrackTagVector.resize(trackVector.size(),std::numeric_limits<size_t>::max());
74  cosmicTagVector.reserve(trackVector.size());
75 
76  for(size_t track_i=0; track_i<trackVector.size(); track_i++){
77 
78  recob::Track const& track(trackVector[track_i]);
79 
80  if(track.Length() < fMinTrackLength) continue;
81 
82  //get the begin and end points of this track
83  TVector3 const& pt_begin = track.LocationAtPoint<TVector3>(0);
84  TVector3 const& pt_end = track.LocationAtPoint<TVector3>(track.NumberTrajectoryPoints()-1);
85  std::vector<float> xyz_begin = { (float)pt_begin.x(), (float)pt_begin.y(), (float)pt_begin.z()};
86  std::vector<float> xyz_end = {(float)pt_end.x(), (float)pt_end.y(), (float)pt_end.z()};
87 
88  //check if this track is outside the drift window, and if it is continue
89  if(!InDriftWindow(pt_begin.x(),pt_end.x(),geom)) {
91  cosmicTagVector.emplace_back(xyz_begin,xyz_end,1.,COSMIC_TYPE_OUTSIDEDRIFT);
92  assnTrackTagVector[track_i] = cosmicTagVector.size()-1;
93  }
94  continue;
95  }
96 
97  //get light hypothesis for track
98  std::vector<float> lightHypothesis = GetMIPHypotheses(track,providers,pvs,opdigip);
99 
100  //check compatibility with beam flash
101  bool compatible=false;
102  for(const recob::OpFlash* flashPointer : flashesOnBeamTime){
103  CompatibilityResultType result = CheckCompatibility(lightHypothesis,flashPointer, geom);
104  if(result==CompatibilityResultType::kCompatible) compatible=true;
105  if(DEBUG_FLAG){
106  PrintTrackProperties(track);
107  PrintFlashProperties(*flashPointer);
108  PrintHypothesisFlashComparison(lightHypothesis,flashPointer,geom,result);
109  }
110  }
111 
112  //make tag
113  float cosmicScore=1.;
114  if(compatible) cosmicScore=0.;
115  cosmicTagVector.emplace_back(xyz_begin,xyz_end,cosmicScore,COSMIC_TYPE_FLASHMATCH);
116  assnTrackTagVector[track_i]=cosmicTagVector.size()-1;
117  }
118 
119 }
120 
121 //this compares the hypothesis to the flash itself.
123  unsigned int const event,
124  std::vector<recob::OpFlash> const& flashVector,
125  std::vector<recob::Track> const& trackVector,
126  Providers_t providers,
128  opdet::OpDigiProperties const& opdigip){
129 
130  auto const& geom = *(providers.get<geo::GeometryCore>());
131 
132  cFlashComparison_p.run = run;
133  cFlashComparison_p.event = event;
134 
135  std::vector< std::pair<unsigned int, const recob::OpFlash*> > flashesOnBeamTime;
136  for(unsigned int i=0; i<flashVector.size(); i++){
137  recob::OpFlash const& flash = flashVector[i];
138  if(!flash.OnBeamTime()) continue;
139  flashesOnBeamTime.push_back(std::make_pair(i,&flash));
140  }
141 
142  for(size_t track_i=0; track_i<trackVector.size(); track_i++){
143 
144  recob::Track const& track(trackVector[track_i]);
145 
146  if(track.Length() < fMinTrackLength) continue;
147 
148  //get the begin and end points of this track
149  TVector3 const& pt_begin = track.LocationAtPoint<TVector3>(0);
150  TVector3 const& pt_end = track.LocationAtPoint<TVector3>(track.NumberTrajectoryPoints()-1);
151  std::vector<float> xyz_begin = { (float)pt_begin.x(), (float)pt_begin.y(), (float)pt_begin.z()};
152  std::vector<float> xyz_end = {(float)pt_end.x(), (float)pt_end.y(), (float)pt_end.z()};
153 
154  //check if this track is outside the drift window, and if it is continue
155  if(!InDriftWindow(pt_begin.x(),pt_end.x(),geom)) continue;
156 
157  cFlashComparison_p.trk_startx = pt_begin.x();
158  cFlashComparison_p.trk_starty = pt_begin.y();
159  cFlashComparison_p.trk_startz = pt_begin.z();
160  cFlashComparison_p.trk_endx = pt_end.x();
161  cFlashComparison_p.trk_endy = pt_end.y();
162  cFlashComparison_p.trk_endz = pt_end.z();
163 
164  //get light hypothesis for track
165  cOpDetVector_hyp = GetMIPHypotheses(track,providers,pvs,opdigip);
166 
167  cFlashComparison_p.hyp_index = track_i;
172  geom);
173 
174  for(auto flash : flashesOnBeamTime){
175  cOpDetVector_flash = std::vector<float>(geom.NOpDets(),0);
177  for(size_t c=0; c<=geom.MaxOpChannel(); c++){
178  if ( geom.IsValidOpChannel( c ) ) {
179  unsigned int OpDet = geom.OpDetFromOpChannel(c);
180  cOpDetVector_flash[OpDet] += flash.second->PE(c);
181  }
182  }
183  for(size_t o=0; o<cOpDetVector_flash.size(); o++)
185 
186  cFlashComparison_p.flash_index = flash.first;
187  cFlashComparison_p.flash_totalPE = flash.second->TotalPE();
188  cFlashComparison_p.flash_y = flash.second->YCenter();
189  cFlashComparison_p.flash_sigmay = flash.second->YWidth();
190  cFlashComparison_p.flash_z = flash.second->ZCenter();
191  cFlashComparison_p.flash_sigmaz = flash.second->ZWidth();
192 
194 
195  cTree->Fill();
196  }//end loop over flashes
197 
198  }//end loop over tracks
199 
200 
201 }
202 
203 //this compares the hypothesis to the flash itself.
205  unsigned int const event,
206  std::vector<recob::OpFlash> const& flashVector,
207  std::vector<simb::MCParticle> const& mcParticleVector,
208  Providers_t providers,
210  opdet::OpDigiProperties const& opdigip){
211 
212  auto const& geom = *(providers.get<geo::GeometryCore>());
213 
214  cFlashComparison_p.run = run;
215  cFlashComparison_p.event = event;
216 
217  std::vector< std::pair<unsigned int, const recob::OpFlash*> > flashesOnBeamTime;
218  for(unsigned int i=0; i<flashVector.size(); i++){
219  recob::OpFlash const& flash = flashVector[i];
220  if(!flash.OnBeamTime()) continue;
221  flashesOnBeamTime.push_back(std::make_pair(i,&flash));
222  }
223 
224  for(size_t particle_i=0; particle_i<mcParticleVector.size(); particle_i++){
225 
226  simb::MCParticle const& particle(mcParticleVector[particle_i]);
227  if(particle.Process().compare("primary")!=0) continue;
228  if(particle.Trajectory().TotalLength() < fMinTrackLength) continue;
229 
230  //get the begin and end points of this track
231  size_t start_i=0, end_i=particle.NumberTrajectoryPoints()-1;
232  bool prev_inside=false;
233  for(size_t pt_i=0; pt_i < particle.NumberTrajectoryPoints(); pt_i++){
234  bool inside = InDetector(particle.Position(pt_i).Vect(),geom);
235  if(inside && !prev_inside) start_i = pt_i;
236  if(!inside && prev_inside) { end_i = pt_i-1; break; }
237  prev_inside = inside;
238  }
239  TVector3 const& pt_begin = particle.Position(start_i).Vect();
240  TVector3 const& pt_end = particle.Position(end_i).Vect();
241  std::vector<float> xyz_begin = { (float)pt_begin.x(), (float)pt_begin.y(), (float)pt_begin.z()};
242  std::vector<float> xyz_end = {(float)pt_end.x(), (float)pt_end.y(), (float)pt_end.z()};
243 
244  //check if this track is outside the drift window, and if it is continue
245  if(!InDriftWindow(pt_begin.x(),pt_end.x(),geom)) continue;
246 
247  cFlashComparison_p.trk_startx = pt_begin.x();
248  cFlashComparison_p.trk_starty = pt_begin.y();
249  cFlashComparison_p.trk_startz = pt_begin.z();
250  cFlashComparison_p.trk_endx = pt_end.x();
251  cFlashComparison_p.trk_endy = pt_end.y();
252  cFlashComparison_p.trk_endz = pt_end.z();
253 
254  //get light hypothesis for track
255  cOpDetVector_hyp = GetMIPHypotheses(particle,start_i,end_i,providers,pvs,opdigip);
256 
257  cFlashComparison_p.hyp_index = particle_i;
262  geom);
263 
264  for(auto flash : flashesOnBeamTime){
265  cOpDetVector_flash = std::vector<float>(geom.NOpDets(),0);
267  //for(size_t c=0; c<geom.NOpChannels(); c++){
268  for(size_t c=0; c<=geom.MaxOpChannel(); c++){
269  if ( geom.IsValidOpChannel(c) ) {
270  unsigned int OpDet = geom.OpDetFromOpChannel(c);
271  cOpDetVector_flash[OpDet] += flash.second->PE(c);
272  }
273  }
274  for(size_t o=0; o<cOpDetVector_flash.size(); o++)
276  cFlashComparison_p.flash_index = flash.first;
277  cFlashComparison_p.flash_totalPE = flash.second->TotalPE();
278  cFlashComparison_p.flash_y = flash.second->YCenter();
279  cFlashComparison_p.flash_sigmay = flash.second->YWidth();
280  cFlashComparison_p.flash_z = flash.second->ZCenter();
281  cFlashComparison_p.flash_sigmaz = flash.second->ZWidth();
282 
284 
285  //cOpDetHist_flash->SetBins(cOpDetVector_flash.size(),0,cOpDetVector_flash.size());
286  for(size_t i=0; i<cOpDetVector_flash.size(); i++)
287  cOpDetHist_flash->SetBinContent(i+1,cOpDetVector_flash[i]);
288 
289  //cOpDetHist_hyp->SetBins(cOpDetVector_hyp.size(),0,cOpDetVector_hyp.size());
290  for(size_t i=0; i<cOpDetVector_hyp.size(); i++)
291  cOpDetHist_hyp->SetBinContent(i+1,cOpDetVector_hyp[i]);
292 
293  for(size_t i=0; i<cOpDetVector_flash.size(); i++){
294  std::cout << "Flash/Hyp " << i << " : "
295  << cOpDetHist_flash->GetBinContent(i+1) << " "
296  << cOpDetHist_hyp->GetBinContent(i+1) << std::endl;
297  }
298 
299  cTree->Fill();
300  }//end loop over flashes
301 
302  }//end loop over tracks
303 
304 
305 }
306 
307 void cosmic::BeamFlashTrackMatchTaggerAlg::FillFlashProperties(std::vector<float> const& opdetVector,
308  float& sum,
309  float& y, float& sigmay,
310  float& z, float& sigmaz,
311  geo::GeometryCore const& geom){
312  y=0; sigmay=0; z=0; sigmaz=0; sum=0;
313  double xyz[3];
314  for(unsigned int opdet=0; opdet<opdetVector.size(); opdet++){
315  sum+=opdetVector[opdet];
316  geom.Cryostat(0).OpDet(opdet).GetCenter(xyz);
317  y += opdetVector[opdet]*xyz[1];
318  z += opdetVector[opdet]*xyz[2];
319  }
320 
321  y /= sum; z /= sum;
322 
323  for(unsigned int opdet=0; opdet<opdetVector.size(); opdet++){
324  geom.Cryostat(0).OpDet(opdet).GetCenter(xyz);
325  sigmay += (opdetVector[opdet]*xyz[1]-y)*(opdetVector[opdet]*xyz[1]-y);
326  sigmaz += (opdetVector[opdet]*xyz[2]-y)*(opdetVector[opdet]*xyz[2]-y);
327  }
328 
329  sigmay = std::sqrt(sigmay)/sum;
330  sigmaz = std::sqrt(sigmaz)/sum;
331 
332 }
333 
335  if(pt.x() < 0 || pt.x() > 2*geom.DetHalfWidth()) return false;
336  if(std::abs(pt.y()) > geom.DetHalfHeight()) return false;
337  if(pt.z() < 0 || pt.z() > geom.DetLength()) return false;
338  return true;
339 }
340 
341 bool cosmic::BeamFlashTrackMatchTaggerAlg::InDriftWindow(double start_x, double end_x, geo::GeometryCore const& geom){
342  if(start_x < 0. || end_x < 0.) return false;
343  if(start_x > 2*geom.DetHalfWidth() || end_x > 2*geom.DetHalfWidth()) return false;
344  return true;
345 }
346 
348  TVector3 const& pt2,
349  std::vector<float> & lightHypothesis,
350  float & totalHypothesisPE,
351  geo::GeometryCore const& geom,
353  float const& PromptMIPScintYield,
354  float XOffset){
355 
356  double xyz_segment[3];
357  xyz_segment[0] = 0.5*(pt2.x()+pt1.x()) + XOffset;
358  xyz_segment[1] = 0.5*(pt2.y()+pt1.y());
359  xyz_segment[2] = 0.5*(pt2.z()+pt1.z());
360 
361  //get the visibility vector
362  const float* PointVisibility = pvs.GetAllVisibilities(xyz_segment);
363 
364  //check pointer, as it may be null if given a y/z outside some range
365  if(!PointVisibility) return;
366 
367  //get the amount of light
368  float LightAmount = PromptMIPScintYield*(pt2-pt1).Mag();
369 
370  for(size_t opdet_i=0; opdet_i<pvs.NOpChannels(); opdet_i++){
371  lightHypothesis[opdet_i] += PointVisibility[opdet_i]*LightAmount;
372  totalHypothesisPE += PointVisibility[opdet_i]*LightAmount;
373 
374  //apply saturation limit
375  if(lightHypothesis[opdet_i]>fOpDetSaturation){
376  totalHypothesisPE -= (lightHypothesis[opdet_i]-fOpDetSaturation);
377  lightHypothesis[opdet_i] = fOpDetSaturation;
378  }
379  }
380 
381 }//end AddLightFromSegment
382 
383 void cosmic::BeamFlashTrackMatchTaggerAlg::NormalizeLightHypothesis(std::vector<float> & lightHypothesis,
384  float const& totalHypothesisPE,
385  geo::GeometryCore const& geom){
386  for(size_t opdet_i=0; opdet_i<geom.NOpDets(); opdet_i++)
387  lightHypothesis[opdet_i] /= totalHypothesisPE;
388 }
389 
390 
391 // Get a hypothesis for the light collected for a track
393  Providers_t providers,
395  opdet::OpDigiProperties const& opdigip,
396  float XOffset)
397 {
398  auto const& geom = *(providers.get<geo::GeometryCore>());
399  auto const& larp = *(providers.get<detinfo::LArProperties>());
400  std::vector<float> lightHypothesis(geom.NOpDets(),0);
401  float totalHypothesisPE=0;
402  const float PromptMIPScintYield = larp.ScintYield()*larp.ScintYieldRatio()*opdigip.QE()*fMIPdQdx;
403 
404  //get QE from ubChannelConfig, which gives per tube, so goes in AddLightFromSegment
405  //VisibleEnergySeparation(step);
406 
407  for(size_t pt=1; pt<track.NumberTrajectoryPoints(); pt++)
408  AddLightFromSegment(track.LocationAtPoint<TVector3>(pt-1),track.LocationAtPoint<TVector3>(pt),
409  lightHypothesis,totalHypothesisPE,
410  geom,pvs,PromptMIPScintYield,
411  XOffset);
412 
413  if(fNormalizeHypothesisToFlash && totalHypothesisPE > std::numeric_limits<float>::epsilon())
414  NormalizeLightHypothesis(lightHypothesis,totalHypothesisPE,geom);
415 
416  return lightHypothesis;
417 
418 }//end GetMIPHypotheses
419 
420 
421 // Get a hypothesis for the light collected for a particle trajectory
423  size_t start_i, size_t end_i,
424  Providers_t providers,
426  opdet::OpDigiProperties const& opdigip,
427  float XOffset)
428 {
429  auto const& geom = *(providers.get<geo::GeometryCore>());
430  auto const& larp = *(providers.get<detinfo::LArProperties>());
431  std::vector<float> lightHypothesis(geom.NOpDets(),0);
432  float totalHypothesisPE=0;
433  const float PromptMIPScintYield = larp.ScintYield()*larp.ScintYieldRatio()*opdigip.QE()*fMIPdQdx;
434 
435  for(size_t pt=start_i+1; pt<=end_i; pt++)
436  AddLightFromSegment(particle.Position(pt-1).Vect(),particle.Position(pt).Vect(),
437  lightHypothesis,totalHypothesisPE,
438  geom,pvs,PromptMIPScintYield,
439  XOffset);
440 
441  if(fNormalizeHypothesisToFlash && totalHypothesisPE > std::numeric_limits<float>::epsilon())
442  NormalizeLightHypothesis(lightHypothesis,totalHypothesisPE,geom);
443 
444  return lightHypothesis;
445 
446 }//end GetMIPHypotheses
447 
448 
449 //---------------------------------------
450 // Check whether a hypothesis can be accomodated in a flash
451 // Flashes fail if 1 bin is far in excess of the observed signal
452 // or if the whole flash intensity is much too large for the hypothesis.
453 // MIP dEdx is assumed for now. Accounting for real dQdx will
454 // improve performance of this algorithm.
455 //---------------------------------------
457 cosmic::BeamFlashTrackMatchTaggerAlg::CheckCompatibility(std::vector<float> const& lightHypothesis,
458  const recob::OpFlash* flashPointer,
459  geo::GeometryCore const& geom)
460 {
461  float hypothesis_integral=0;
462  float flash_integral=0;
463  unsigned int cumulativeChannels=0;
464 
465  std::vector<double> PEbyOpDet(geom.NOpDets(),0);
466  //for (unsigned int c = 0; c < geom.NOpChannels(); c++){
467  for (unsigned int c = 0; c <= geom.MaxOpChannel(); c++){
468  if ( geom.IsValidOpChannel(c) ) {
469  unsigned int o = geom.OpDetFromOpChannel(c);
470  PEbyOpDet[o] += flashPointer->PE(c);
471  }
472  }
473 
474  float hypothesis_scale=1.;
475  if(fNormalizeHypothesisToFlash) hypothesis_scale = flashPointer->TotalPE();
476 
477  for(size_t pmt_i=0; pmt_i<lightHypothesis.size(); pmt_i++){
478 
479  flash_integral += PEbyOpDet[pmt_i];
480 
481  if(lightHypothesis[pmt_i] < std::numeric_limits<float>::epsilon() ) continue;
482  hypothesis_integral += lightHypothesis[pmt_i]*hypothesis_scale;
483 
484  if(PEbyOpDet[pmt_i] < fMinOpHitPE) continue;
485 
486  float diff_scaled = (lightHypothesis[pmt_i]*hypothesis_scale - PEbyOpDet[pmt_i])/std::sqrt(lightHypothesis[pmt_i]*hypothesis_scale);
487 
488  if( diff_scaled > fSingleChannelCut ) return CompatibilityResultType::kSingleChannelCut;
489 
490  if( diff_scaled > fCumulativeChannelThreshold ) cumulativeChannels++;
491  if(cumulativeChannels >= fCumulativeChannelCut) return CompatibilityResultType::kCumulativeChannelCut;
492 
493  }
494 
495  if( (hypothesis_integral - flash_integral)/std::sqrt(hypothesis_integral)
496  > fIntegralCut) return CompatibilityResultType::kIntegralCut;
497 
498  return CompatibilityResultType::kCompatible;
499 }
500 
501 
502 float cosmic::BeamFlashTrackMatchTaggerAlg::CalculateChi2(std::vector<float> const& light_flash,
503  std::vector<float> const& light_track){
504 
505  float chi2=0;
506  for(size_t pmt_i=0; pmt_i<light_flash.size(); pmt_i++){
507 
508  if(light_flash[pmt_i] < fMinOpHitPE) continue;
509 
510  float err2 = 1;
511  if(light_track[pmt_i] > 1) err2 = light_track[pmt_i];
512 
513  chi2 += (light_flash[pmt_i]-light_track[pmt_i])*(light_flash[pmt_i]-light_track[pmt_i]) / err2;
514  }
515 
516  return chi2;
517 }
518 
519 
521 {
522  *output << "----------------------------------------------" << std::endl;
523  *output << "Track properties: ";
524  *output << "\n\tLength=" << track.Length();
525 
526  auto const& pt_begin = track.LocationAtPoint(0);
527  *output << "\n\tBegin Location (x,y,z)=(" << pt_begin.x() << "," << pt_begin.y() << "," << pt_begin.z() << ")";
528 
529  auto const& pt_end = track.LocationAtPoint(track.NumberTrajectoryPoints()-1);
530  *output << "\n\tEnd Location (x,y,z)=(" << pt_end.x() << "," << pt_end.y() << "," << pt_end.z() << ")";
531 
532  *output << "\n\tTrajectoryPoints=" << track.NumberTrajectoryPoints();
533  *output << std::endl;
534  *output << "----------------------------------------------" << std::endl;
535 }
536 
538 {
539  *output << "----------------------------------------------" << std::endl;
540  *output << "Flash properties: ";
541 
542  *output << "\n\tTime=" << flash.Time();
543  *output << "\n\tOnBeamTime=" << flash.OnBeamTime();
544  *output << "\n\ty position (center,width)=(" << flash.YCenter() << "," << flash.YWidth() << ")";
545  *output << "\n\tz position (center,width)=(" << flash.ZCenter() << "," << flash.ZWidth() << ")";
546  *output << "\n\tTotal PE=" << flash.TotalPE();
547 
548  *output << std::endl;
549  *output << "----------------------------------------------" << std::endl;
550 
551 }
552 
553 void cosmic::BeamFlashTrackMatchTaggerAlg::PrintHypothesisFlashComparison(std::vector<float> const& lightHypothesis,
554  const recob::OpFlash* flashPointer,
555  geo::GeometryCore const& geom,
557  std::ostream* output)
558 {
559 
560  *output << "----------------------------------------------" << std::endl;
561  *output << "Hypothesis-flash comparison: ";
562 
563  float hypothesis_integral=0;
564  float flash_integral=0;
565 
566  float hypothesis_scale=1.;
567  if(fNormalizeHypothesisToFlash) hypothesis_scale = flashPointer->TotalPE();
568 
569 
570  std::vector<double> PEbyOpDet(geom.NOpDets(),0);
571  //for (unsigned int c = 0; c < geom.NOpChannels(); c++){
572  for ( unsigned int c = 0; c <= geom.MaxOpChannel(); c++){
573  if ( geom.IsValidOpChannel(c) ) {
574  unsigned int o = geom.OpDetFromOpChannel(c);
575  PEbyOpDet[o] += flashPointer->PE(c);
576  }
577  }
578 
579  for(size_t pmt_i=0; pmt_i<lightHypothesis.size(); pmt_i++){
580 
581  flash_integral += PEbyOpDet[pmt_i];
582 
583  *output << "\n\t pmt_i=" << pmt_i << ", (hypothesis,flash)=("
584  << lightHypothesis[pmt_i]*hypothesis_scale << "," << PEbyOpDet[pmt_i] << ")";
585 
586  if(lightHypothesis[pmt_i] < std::numeric_limits<float>::epsilon() ) continue;
587 
588  *output << " difference="
589  << (lightHypothesis[pmt_i]*hypothesis_scale - PEbyOpDet[pmt_i])/std::sqrt(lightHypothesis[pmt_i]*hypothesis_scale);
590 
591  hypothesis_integral += lightHypothesis[pmt_i]*hypothesis_scale;
592  }
593 
594  *output << "\n\t TOTAL (hypothesis,flash)=("
595  << hypothesis_integral << "," << flash_integral << ")"
596  << " difference=" << (hypothesis_integral - flash_integral)/std::sqrt(hypothesis_integral);
597 
598  *output << std::endl;
599  *output << "End result=" << result << std::endl;
600  *output << "----------------------------------------------" << std::endl;
601 
602 }
unsigned int NumberTrajectoryPoints() const
Definition: MCParticle.h:222
const TLorentzVector & Position(const int i=0) const
Definition: MCParticle.h:223
CompatibilityResultType CheckCompatibility(std::vector< float > const &lightHypothesis, const recob::OpFlash *flashPointer, geo::GeometryCore const &geom)
void reconfigure(fhicl::ParameterSet const &p)
geo::Length_t DetHalfWidth(geo::TPCID const &tpcid) const
Returns the half width of the active volume of the specified TPC.
BeamFlashTrackMatchTaggerAlg(fhicl::ParameterSet const &p)
const simb::MCTrajectory & Trajectory() const
Definition: MCParticle.h:257
enum anab::cosmic_tag_id CosmicTagID_t
Float_t y
Definition: compare.C:6
Point_t const & LocationAtPoint(size_t i) const
Access to track position at different points.
Definition: Track.h:129
Double_t z
Definition: plot.C:279
void RunHypothesisComparison(unsigned int const, unsigned int const, std::vector< recob::OpFlash > const &, std::vector< recob::Track > const &, Providers_t, phot::PhotonVisibilityService const &, opdet::OpDigiProperties const &)
float CalculateChi2(std::vector< float > const &, std::vector< float > const &)
size_t NumberTrajectoryPoints() const
Various functions related to the presence and the number of (valid) points.
Definition: Track.h:105
Provider const * get() const
Returns the provider with the specified type.
Definition: ProviderPack.h:193
void GetCenter(double *xyz, double localz=0.0) const
Definition: OpDetGeo.cxx:47
double PE(unsigned int i) const
Definition: OpFlash.h:87
std::string Process() const
Definition: MCParticle.h:219
void RunCompatibilityCheck(std::vector< recob::OpFlash > const &, std::vector< recob::Track > const &, std::vector< anab::CosmicTag > &, std::vector< size_t > &, Providers_t, phot::PhotonVisibilityService const &, opdet::OpDigiProperties const &)
unsigned int OpDetFromOpChannel(int opChannel) const
Convert unique channel to detector number.
double ZCenter() const
Definition: OpFlash.h:91
std::vector< float > GetMIPHypotheses(recob::Track const &track, Providers_t providers, phot::PhotonVisibilityService const &pvs, opdet::OpDigiProperties const &, float XOffset=0)
double Time() const
Definition: OpFlash.h:83
Int_t max
Definition: plot.C:27
geo::Length_t DetHalfHeight(geo::TPCID const &tpcid) const
Returns the half height of the active volume of the specified TPC.
double Length(size_t p=0) const
Access to various track properties.
Definition: Track.h:170
float const * GetAllVisibilities(double const *xyz, bool wantReflected=false) const
void NormalizeLightHypothesis(std::vector< float > &lightHypothesis, float const &totalHypothesisPE, geo::GeometryCore const &geom)
const OpDetGeo & OpDet(unsigned int iopdet) const
Return the iopdet&#39;th optical detector in the cryostat.
int OnBeamTime() const
Definition: OpFlash.h:97
void AddLightFromSegment(TVector3 const &pt1, TVector3 const &pt2, std::vector< float > &lightHypothesis, float &totalHypothesisPE, geo::GeometryCore const &geom, phot::PhotonVisibilityService const &pvs, float const &PromptMIPScintYield, float XOffset)
T get(std::string const &key) const
Definition: ParameterSet.h:231
void PrintTrackProperties(recob::Track const &, std::ostream *output=&std::cout)
TMarker * pt
Definition: egs.C:25
geo::Length_t DetLength(geo::TPCID const &tpcid) const
Returns the length of the active volume of the specified TPC.
CryostatGeo const & Cryostat(geo::CryostatID const &cryoid) const
Returns the specified cryostat.
unsigned int MaxOpChannel() const
Largest optical channel number.
Description of geometry of one entire detector.
void SetHypothesisComparisonTree(TTree *, TH1F *, TH1F *)
unsigned int NOpDets() const
Number of OpDets in the whole detector.
Encapsulate the geometry of an optical detector.
double YWidth() const
Definition: OpFlash.h:90
double TotalLength() const
void PrintHypothesisFlashComparison(std::vector< float > const &, const recob::OpFlash *, geo::GeometryCore const &geom, CompatibilityResultType, std::ostream *output=&std::cout)
Container for a list of pointers to providers.
Definition: ProviderPack.h:114
bool InDriftWindow(double, double, geo::GeometryCore const &)
void FillFlashProperties(std::vector< float > const &opdetVector, float &, float &, float &, float &, float &, geo::GeometryCore const &geom)
double TotalPE() const
Definition: OpFlash.cxx:56
bool InDetector(TVector3 const &, geo::GeometryCore const &)
double YCenter() const
Definition: OpFlash.h:89
Float_t track
Definition: plot.C:34
double QE() const
Returns quantum efficiency.
bool IsValidOpChannel(int opChannel) const
Is this a valid OpChannel number?
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
double ZWidth() const
Definition: OpFlash.h:92
Event finding and building.
void PrintFlashProperties(recob::OpFlash const &, std::ostream *output=&std::cout)