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
lbne::PhotonCounterT0Matching Class Reference
Inheritance diagram for lbne::PhotonCounterT0Matching:

Public Member Functions

 PhotonCounterT0Matching (fhicl::ParameterSet const &p)
 
 PhotonCounterT0Matching (PhotonCounterT0Matching const &)=delete
 
 PhotonCounterT0Matching (PhotonCounterT0Matching &&)=delete
 
PhotonCounterT0Matchingoperator= (PhotonCounterT0Matching const &)=delete
 
PhotonCounterT0Matchingoperator= (PhotonCounterT0Matching &&)=delete
 
void produce (art::Event &e) override
 
void beginJob () override
 

Private Member Functions

void TrackProp (double TrackStart_X, double TrackEnd_X, double &TrackLength_X, double &TrackCentre_X, double TrackStart_Y, double TrackEnd_Y, double &TrackLength_Y, double &TrackCentre_Y, double TrackStart_Z, double TrackEnd_Z, double &TrackLength_Z, double &TrackCentre_Z, double trkTimeStart, double trkTimeEnd, double &trkTimeLengh, double &trkTimeCentre, double &TrackLength)
 
double DistFromPoint (double StartY, double EndY, double StartZ, double EndZ, double PointY, double PointZ)
 

Private Attributes

std::string fTrackModuleLabel
 
std::string fShowerModuleLabel
 
std::string fHitsModuleLabel
 
std::string fFlashModuleLabel
 
std::string fTruthT0ModuleLabel
 
double fPredictedXConstant
 
double fPredictedXPower = 1
 
double fPredictedExpConstant
 
double fPredictedExpGradient
 
double fDriftWindowSize
 
double fWeightOfDeltaYZ
 
double fMatchCriteria
 
double fPEThreshold
 
bool fVerbosity
 
double TrackLength_X
 
double TrackCentre_X
 
double BestTrackCentre_X
 
double TrackLength_Y
 
double TrackCentre_Y
 
double TrackLength_Z
 
double TrackCentre_Z
 
double trkTimeStart
 
double trkTimeEnd
 
double trkTimeLengh
 
double trkTimeCentre
 
double BesttrkTimeCentre
 
double TrackLength
 
double BestTrackLength
 
double PredictedX
 
double BestPredictedX
 
double TimeSepPredX
 
double BestTimeSepPredX
 
double DeltaPredX
 
double BestDeltaPredX
 
double minYZSep
 
double BestminYZSep
 
double FitParam
 
double BestFitParam
 
double FlashTime
 
double BestFlashTime
 
double TimeSep
 
double BestTimeSep
 
int BestFlash
 
int FlashTriggerType = 1
 
double YZSep
 
double MCTruthT0
 
TTree * fTree
 
TH2D * hPredX_T
 
TH2D * hPredX_PE
 
TH2D * hPredX_T_PE
 
TH2D * hdeltaX_deltaYZ
 
TH2D * hdeltaYZ_Length
 
TH2D * hFitParam_Length
 
TH2D * hPhotonT0_MCT0
 
TH1D * hT0_diff_full
 
TH1D * hT0_diff_zoom
 

Detailed Description

Definition at line 65 of file PhotonCounterT0Matching_module.cc.

Constructor & Destructor Documentation

lbne::PhotonCounterT0Matching::PhotonCounterT0Matching ( fhicl::ParameterSet const &  p)
explicit

Definition at line 157 of file PhotonCounterT0Matching_module.cc.

157  : EDProducer{p}
158 {
159  // Call appropriate produces<>() functions here.
160  produces<std::vector<anab::T0>>();
161  produces<art::Assns<recob::Track, anab::T0>>();
162  produces<art::Assns<recob::Shower, anab::T0>>();
163 
164  fTrackModuleLabel = (p.get<std::string>("TrackModuleLabel"));
165  fShowerModuleLabel = (p.get<std::string>("ShowerModuleLabel"));
166  fHitsModuleLabel = (p.get<std::string>("HitsModuleLabel"));
167  fFlashModuleLabel = (p.get<std::string>("FlashModuleLabel"));
168  fTruthT0ModuleLabel = (p.get<std::string>("TruthT0ModuleLabel"));
169 
170  fPredictedXConstant = (p.get<double>("PredictedXConstant"));
171  fPredictedExpConstant = (p.get<double>("PredictedExpConstant"));
172  fPredictedExpGradient = (p.get<double>("PredictedExpGradient"));
173 
174  fDriftWindowSize = (p.get<double>("DriftWindowSize"));
175  fWeightOfDeltaYZ = (p.get<double>("WeightOfDeltaYZ"));
176  fMatchCriteria = (p.get<double>("MatchCriteria"));
177  fPEThreshold = (p.get<double>("PEThreshold"));
178 
179  fVerbosity = (p.get<bool>("Verbose", false));
180 }
pdgs p
Definition: selectors.fcl:22
lbne::PhotonCounterT0Matching::PhotonCounterT0Matching ( PhotonCounterT0Matching const &  )
delete
lbne::PhotonCounterT0Matching::PhotonCounterT0Matching ( PhotonCounterT0Matching &&  )
delete

Member Function Documentation

void lbne::PhotonCounterT0Matching::beginJob ( )
override

Definition at line 183 of file PhotonCounterT0Matching_module.cc.

184 {
185  // Implementation of optional member function here.
186  art::ServiceHandle<art::TFileService const> tfs;
187  fTree = tfs->make<TTree>("PhotonCounterT0Matching", "PhotonCounterT0");
188  fTree->Branch("TrackCentre_X", &BestTrackCentre_X, "TrackCentre_X/D");
189  fTree->Branch("PredictedX", &BestPredictedX, "PredictedX/D");
190  fTree->Branch("TrackTimeCent", &BesttrkTimeCentre, "TimeSepPredX/D");
191  fTree->Branch("FlashTime", &BestFlashTime, "FlashTime/D");
192  fTree->Branch("TimeSep", &BestTimeSep, "TimeSep/D");
193  fTree->Branch("TimeSepPredX", &BestTimeSepPredX, "TimeSepPredX/D");
194  fTree->Branch("minYZSep", &BestminYZSep, "minYZSep/D");
195  fTree->Branch("FitParam", &BestFitParam, "FitParam/D");
196  fTree->Branch("MCTruthT0", &MCTruthT0, "MCTruthT0/D");
197 
198  hPredX_T = tfs->make<TH2D>("hPredX_T",
199  "Predicted X from timing information against reconstructed X; "
200  "Reconstructed X (cm); Predicted X (cm)",
201  30,
202  0,
203  300,
204  30,
205  0,
206  300);
207  hPredX_PE = tfs->make<TH2D>("hPredX_PE",
208  "Predicted X from PE information against reconstructed X; "
209  "Reconstructed X (cm); Predicted X (cm)",
210  30,
211  0,
212  300,
213  30,
214  0,
215  300);
216  hPredX_T_PE = tfs->make<TH2D>("hPredX_T_PE",
217  "Predicted X position from time and PE information; Predicted X "
218  "from timing information (cm); Predicted X from PE information",
219  60,
220  0,
221  300,
222  60,
223  0,
224  300);
225  hdeltaX_deltaYZ = tfs->make<TH2D>(
226  "hdeltaX_deltaYZ",
227  "Difference between X predicted from PE's and T agaisnt distance of flash from track in YZ; "
228  "Difference in X predicted from PE's and T (cm); Distance of flash from track in YZ (cm)",
229  40,
230  0,
231  200,
232  40,
233  0,
234  100);
235  hdeltaYZ_Length = tfs->make<TH2D>("hdeltaYZ_Length",
236  "Distance of flash from track against track length; Distance "
237  "from flash to track (cm); Track length (cm)",
238  20,
239  0,
240  100,
241  60,
242  0,
243  300);
245  tfs->make<TH2D>("hFitParam_Length",
246  "How fit correlates with track length; Fit correlation; Track Length (cm)",
247  50,
248  0,
249  250,
250  30,
251  0,
252  300);
253  hPhotonT0_MCT0 = tfs->make<TH2D>("hPhotonT0_MCT0",
254  "Comparing Photon Counter reconstructed T0 against MCTruth T0; "
255  "Photon Counter T0 (us); MCTruthT0 T0 (us)",
256  1760,
257  -1600,
258  16000,
259  1760,
260  -1600,
261  16000);
262  hT0_diff_full = tfs->make<TH1D>(
263  "hT0_diff_full",
264  "Difference between MCTruth T0 and photon detector T0; Time difference (us); Number",
265  500,
266  -2500,
267  2500);
268  hT0_diff_zoom = tfs->make<TH1D>(
269  "hT0_diff_zoom",
270  "Difference between MCTruth T0 and photon detector T0; Time difference (us); Number",
271  320,
272  -1.6,
273  1.6);
274 }
art::ServiceHandle< art::TFileService > tfs
double lbne::PhotonCounterT0Matching::DistFromPoint ( double  StartY,
double  EndY,
double  StartZ,
double  EndZ,
double  PointY,
double  PointZ 
)
private

Calculate the distance between the centre of the flash and the centre of a line connecting two adjacent space points.

Definition at line 527 of file PhotonCounterT0Matching_module.cc.

533 {
534  ///Calculate the distance between the centre of the flash and the centre of a line connecting two adjacent space points.
535  double Length = hypot(fabs(EndY - StartY), fabs(EndZ - StartZ));
536  double distance =
537  ((PointZ - StartZ) * (EndY - StartY) - (PointY - StartY) * (EndZ - StartZ)) / Length;
538  return fabs(distance);
539 }
float Length(const PFPStruct &pfp)
Definition: PFPUtils.cxx:3304
double distance(geo::Point_t const &point, CathodeDesc_t const &cathode)
Returns the distance of a point from the cathode.
PhotonCounterT0Matching& lbne::PhotonCounterT0Matching::operator= ( PhotonCounterT0Matching const &  )
delete
PhotonCounterT0Matching& lbne::PhotonCounterT0Matching::operator= ( PhotonCounterT0Matching &&  )
delete
void lbne::PhotonCounterT0Matching::produce ( art::Event &  e)
override

Definition at line 277 of file PhotonCounterT0Matching_module.cc.

278 {
279  // Access art services...
280  art::ServiceHandle<geo::Geometry const> geom;
281  auto const clock_data =
282  art::ServiceHandle<detinfo::DetectorClocksService const>()->DataFor(evt);
283  auto const detprop =
284  art::ServiceHandle<detinfo::DetectorPropertiesService const>()->DataFor(evt, clock_data);
285 
286  //TrackList handle
287  art::Handle<std::vector<recob::Track>> trackListHandle;
288  std::vector<art::Ptr<recob::Track>> tracklist;
289  if (evt.getByLabel(fTrackModuleLabel, trackListHandle))
290  art::fill_ptr_vector(tracklist, trackListHandle);
291 
292  //ShowerList handle
293  art::Handle<std::vector<recob::Shower>> showerListHandle;
294  std::vector<art::Ptr<recob::Shower>> showerlist;
295  if (evt.getByLabel(fShowerModuleLabel, showerListHandle))
296  art::fill_ptr_vector(showerlist, showerListHandle);
297 
298  //HitList Handle
299  art::Handle<std::vector<recob::Hit>> hitListHandle;
300  std::vector<art::Ptr<recob::Hit>> hitlist;
301  if (evt.getByLabel(fHitsModuleLabel, hitListHandle)) art::fill_ptr_vector(hitlist, hitListHandle);
302 
303  //FlashList Handle
304  art::Handle<std::vector<recob::OpFlash>> flashListHandle;
305  std::vector<art::Ptr<recob::OpFlash>> flashlist;
306  if (evt.getByLabel(fFlashModuleLabel, flashListHandle))
307  art::fill_ptr_vector(flashlist, flashListHandle);
308 
309  // Create anab::T0 objects and make association with recob::Track
310 
311  std::unique_ptr<std::vector<anab::T0>> T0col(new std::vector<anab::T0>);
312  std::unique_ptr<art::Assns<recob::Track, anab::T0>> Trackassn(
313  new art::Assns<recob::Track, anab::T0>);
314  std::unique_ptr<art::Assns<recob::Shower, anab::T0>> Showerassn(
315  new art::Assns<recob::Shower, anab::T0>);
316 
317  if (trackListHandle.isValid() && flashListHandle.isValid()) {
318  //Access tracks and hits
319  art::FindManyP<recob::Hit> fmtht(trackListHandle, evt, fTrackModuleLabel);
320  art::FindMany<anab::T0> fmtruth(trackListHandle, evt, fTruthT0ModuleLabel);
321 
322  size_t NTracks = tracklist.size();
323  size_t NFlashes = flashlist.size();
324 
325  if (fVerbosity)
326  std::cout << "There were " << NTracks << " tracks and " << NFlashes
327  << " flashes in this event." << std::endl;
328 
329  // Now to access PhotonCounter for each track...
330  for (size_t iTrk = 0; iTrk < NTracks; ++iTrk) {
331  if (fVerbosity) std::cout << "\n New Track " << (int)iTrk << std::endl;
332  // Reset Variables.
335  bool ValidTrack = false;
336 
337  // Work out Properties of the track.
338  recob::Track::Point_t trackStart, trackEnd;
339  std::tie(trackStart, trackEnd) = tracklist[iTrk]->Extent();
340  std::vector<art::Ptr<recob::Hit>> allHits = fmtht.at(iTrk);
341  size_t nHits = allHits.size();
342  trkTimeStart = allHits[nHits - 1]->PeakTime() /
343  clock_data.TPCClock().Frequency(); // Got in ticks, now in us!
344  trkTimeEnd =
345  allHits[0]->PeakTime() / clock_data.TPCClock().Frequency(); // Got in ticks, now in us!
346  TrackProp(trackStart.X(),
347  trackEnd.X(),
350  trackStart.Y(),
351  trackEnd.Y(),
354  trackStart.Z(),
355  trackEnd.Z(),
358  trkTimeStart,
359  trkTimeEnd,
360  trkTimeLengh,
361  trkTimeCentre, // times in us!
362  TrackLength);
363 
364  // Some cout statement about track properties.
365  if (fVerbosity) {
366  std::cout << trackStart.X() << " " << trackEnd.X() << " " << TrackLength_X << " "
367  << TrackCentre_X << "\n"
368  << trackStart.Y() << " " << trackEnd.Y() << " " << TrackLength_Y << " "
369  << TrackCentre_Y << "\n"
370  << trackStart.Z() << " " << trackEnd.Z() << " " << TrackLength_Z << " "
371  << TrackCentre_Z << "\n"
372  << trkTimeStart << " " << trkTimeEnd << " " << trkTimeLengh << " "
373  << trkTimeCentre << std::endl;
374  }
375  // ----- Loop over flashes ------
376  for (size_t iFlash = 0; iFlash < NFlashes; ++iFlash) {
377  //Reset some flash specific quantities
378  YZSep = minYZSep = 9999;
379  FlashTime = TimeSep = 9999;
381  // Check flash could be caused by track...
382  FlashTime = flashlist[iFlash]->Time(); // Got in us!
383  TimeSep = trkTimeCentre - FlashTime; // Time in us!
384  if (TimeSep < 0 || TimeSep > (fDriftWindowSize / clock_data.TPCClock().Frequency()))
385  continue; // Times compared in us!
386 
387  // Check flash has enough PE's to satisfy our threshold
388  if (flashlist[iFlash]->TotalPE() < fPEThreshold) continue;
389 
390  // Work out some quantities for this flash...
391  // PredictedX = ( A / x^n ) + exp ( B + Cx )
392  PredictedX =
393  (fPredictedXConstant / pow(flashlist[iFlash]->TotalPE(), fPredictedXPower)) +
394  (exp(fPredictedExpConstant + (fPredictedExpGradient * flashlist[iFlash]->TotalPE())));
395  TimeSepPredX = TimeSep * detprop.DriftVelocity(); // us * cm/us = cm!
397  // Dependant on each point...
398  for (size_t Point = 1; Point < tracklist[iTrk]->NumberTrajectoryPoints(); ++Point) {
399  auto NewPoint = tracklist[iTrk]->LocationAtPoint(Point);
400  auto PrevPoint = tracklist[iTrk]->LocationAtPoint(Point - 1);
401  YZSep = DistFromPoint(NewPoint.Y(),
402  PrevPoint.Y(),
403  NewPoint.Z(),
404  PrevPoint.Z(),
405  flashlist[iFlash]->YCenter(),
406  flashlist[iFlash]->ZCenter());
407  if (Point == 1) minYZSep = YZSep;
408  if (YZSep < minYZSep) minYZSep = YZSep;
409  }
410 
411  // Determine how well matched this track is......
412  if (fMatchCriteria == 0)
413  FitParam =
414  pow(((DeltaPredX * DeltaPredX) + (minYZSep * minYZSep * fWeightOfDeltaYZ)), 0.5);
415  else if (fMatchCriteria == 1)
416  FitParam = minYZSep;
417  else if (fMatchCriteria == 2)
419 
420  //----FLASH INFO-----
421  if (fVerbosity) {
422  std::cout << "\nFlash " << (int)iFlash << " " << TrackCentre_X << ", " << TimeSepPredX
423  << " - " << PredictedX << " = " << DeltaPredX << ", " << minYZSep << " -> "
424  << FitParam << std::endl;
425  }
426  //----Select best flash------
427  if (FitParam < BestFitParam) {
428  ValidTrack = true;
429  BestFlash = (int)iFlash;
438  BestFlashTime = FlashTime;
440  } // Find best Flash
441  } // Loop over Flashes
442 
443  // ---- Now Make association and fill TTree/Histos with the best matched flash.....
444  if (ValidTrack) {
445 
446  // -- Fill Histos --
453  // ------ Compare Photon Matched to MCTruth Matched -------
454  if (fmtruth.isValid()) {
455  std::vector<const anab::T0*> T0s = fmtruth.at((int)iTrk);
456  for (size_t i = 0; i < T0s.size(); ++i) {
457  MCTruthT0 = T0s[i]->Time() / 1e3; // Got in ns, now in us!!
458  hPhotonT0_MCT0->Fill(BestFlashTime, MCTruthT0);
459  hT0_diff_full->Fill(MCTruthT0 - BestFlashTime);
460  hT0_diff_zoom->Fill(MCTruthT0 - BestFlashTime);
461  }
462  }
463  // -- Fill TTree --
464  fTree->Fill();
465  //Make Association
466  T0col->push_back(anab::T0(
467  BestFlashTime * 1e3, FlashTriggerType, (int)BestFlash, (*T0col).size(), BestFitParam));
468  util::CreateAssn(*this, evt, *T0col, tracklist[iTrk], *Trackassn);
469  } // Valid Track
470  } // Loop over tracks
471  }
472  evt.put(std::move(T0col));
473  evt.put(std::move(Trackassn));
474  evt.put(std::move(Showerassn));
475 
476 } // Produce
Definition: T0.h:16
double DistFromPoint(double StartY, double EndY, double StartZ, double EndZ, double PointY, double PointZ)
std::tuple< double, double, const reco::ClusterHit3D * > Point
Definitions used by the VoronoiDiagram algorithm.
Definition: DCEL.h:44
void TrackProp(double TrackStart_X, double TrackEnd_X, double &TrackLength_X, double &TrackCentre_X, double TrackStart_Y, double TrackEnd_Y, double &TrackLength_Y, double &TrackCentre_Y, double TrackStart_Z, double TrackEnd_Z, double &TrackLength_Z, double &TrackCentre_Z, double trkTimeStart, double trkTimeEnd, double &trkTimeLengh, double &trkTimeCentre, double &TrackLength)
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.
tracking::Point_t Point_t
TCEvent evt
Definition: DataStructs.cxx:8
BEGIN_PROLOG could also be cout
void lbne::PhotonCounterT0Matching::TrackProp ( double  TrackStart_X,
double  TrackEnd_X,
double &  TrackLength_X,
double &  TrackCentre_X,
double  TrackStart_Y,
double  TrackEnd_Y,
double &  TrackLength_Y,
double &  TrackCentre_Y,
double  TrackStart_Z,
double  TrackEnd_Z,
double &  TrackLength_Z,
double &  TrackCentre_Z,
double  trkTimeStart,
double  trkTimeEnd,
double &  trkTimeLengh,
double &  trkTimeCentre,
double &  TrackLength 
)
private

Calculate central values for track X, Y, Z and time, as well as lengths and overall track length.

Definition at line 479 of file PhotonCounterT0Matching_module.cc.

496 {
497  ///Calculate central values for track X, Y, Z and time, as well as lengths and overall track length.
498  TrackLength_X = fabs(TrackEnd_X - TrackStart_X);
499  if (TrackStart_X < TrackEnd_X)
500  TrackCentre_X = TrackStart_X + 0.5 * TrackLength_X;
501  else
502  TrackCentre_X = TrackStart_X - 0.5 * TrackLength_X;
503 
504  TrackLength_Y = fabs(TrackEnd_Y - TrackStart_Y);
505  if (TrackStart_Y < TrackEnd_Y)
506  TrackCentre_Y = TrackStart_Y + 0.5 * TrackLength_Y;
507  else
508  TrackCentre_Y = TrackStart_Y - 0.5 * TrackLength_Y;
509 
510  TrackLength_Z = fabs(TrackEnd_Z - TrackStart_Z);
511  if (TrackStart_Z < TrackEnd_Z)
512  TrackCentre_Z = TrackStart_Z + 0.5 * TrackLength_Z;
513  else
514  TrackCentre_Z = TrackStart_Z - 0.5 * TrackLength_Z;
515 
517  trkTimeCentre = trkTimeStart + 0.5 * trkTimeLengh;
518 
519  TrackLength = pow(pow((TrackEnd_X - TrackStart_X), 2) + pow((TrackEnd_Y - TrackStart_Y), 2) +
520  pow((TrackEnd_Z - TrackStart_Z), 2),
521  0.5);
522 
523  return;
524 }

Member Data Documentation

double lbne::PhotonCounterT0Matching::BestDeltaPredX
private

Definition at line 135 of file PhotonCounterT0Matching_module.cc.

double lbne::PhotonCounterT0Matching::BestFitParam
private

Definition at line 137 of file PhotonCounterT0Matching_module.cc.

int lbne::PhotonCounterT0Matching::BestFlash
private

Definition at line 140 of file PhotonCounterT0Matching_module.cc.

double lbne::PhotonCounterT0Matching::BestFlashTime
private

Definition at line 138 of file PhotonCounterT0Matching_module.cc.

double lbne::PhotonCounterT0Matching::BestminYZSep
private

Definition at line 136 of file PhotonCounterT0Matching_module.cc.

double lbne::PhotonCounterT0Matching::BestPredictedX
private

Definition at line 133 of file PhotonCounterT0Matching_module.cc.

double lbne::PhotonCounterT0Matching::BestTimeSep
private

Definition at line 139 of file PhotonCounterT0Matching_module.cc.

double lbne::PhotonCounterT0Matching::BestTimeSepPredX
private

Definition at line 134 of file PhotonCounterT0Matching_module.cc.

double lbne::PhotonCounterT0Matching::BestTrackCentre_X
private

Definition at line 126 of file PhotonCounterT0Matching_module.cc.

double lbne::PhotonCounterT0Matching::BestTrackLength
private

Definition at line 132 of file PhotonCounterT0Matching_module.cc.

double lbne::PhotonCounterT0Matching::BesttrkTimeCentre
private

Definition at line 131 of file PhotonCounterT0Matching_module.cc.

double lbne::PhotonCounterT0Matching::DeltaPredX
private

Definition at line 135 of file PhotonCounterT0Matching_module.cc.

double lbne::PhotonCounterT0Matching::fDriftWindowSize
private

Definition at line 119 of file PhotonCounterT0Matching_module.cc.

std::string lbne::PhotonCounterT0Matching::fFlashModuleLabel
private

Definition at line 113 of file PhotonCounterT0Matching_module.cc.

std::string lbne::PhotonCounterT0Matching::fHitsModuleLabel
private

Definition at line 112 of file PhotonCounterT0Matching_module.cc.

double lbne::PhotonCounterT0Matching::FitParam
private

Definition at line 137 of file PhotonCounterT0Matching_module.cc.

double lbne::PhotonCounterT0Matching::FlashTime
private

Definition at line 138 of file PhotonCounterT0Matching_module.cc.

int lbne::PhotonCounterT0Matching::FlashTriggerType = 1
private

Definition at line 141 of file PhotonCounterT0Matching_module.cc.

double lbne::PhotonCounterT0Matching::fMatchCriteria
private

Definition at line 121 of file PhotonCounterT0Matching_module.cc.

double lbne::PhotonCounterT0Matching::fPEThreshold
private

Definition at line 122 of file PhotonCounterT0Matching_module.cc.

double lbne::PhotonCounterT0Matching::fPredictedExpConstant
private

Definition at line 117 of file PhotonCounterT0Matching_module.cc.

double lbne::PhotonCounterT0Matching::fPredictedExpGradient
private

Definition at line 118 of file PhotonCounterT0Matching_module.cc.

double lbne::PhotonCounterT0Matching::fPredictedXConstant
private

Definition at line 115 of file PhotonCounterT0Matching_module.cc.

double lbne::PhotonCounterT0Matching::fPredictedXPower = 1
private

Definition at line 116 of file PhotonCounterT0Matching_module.cc.

std::string lbne::PhotonCounterT0Matching::fShowerModuleLabel
private

Definition at line 111 of file PhotonCounterT0Matching_module.cc.

std::string lbne::PhotonCounterT0Matching::fTrackModuleLabel
private

Definition at line 110 of file PhotonCounterT0Matching_module.cc.

TTree* lbne::PhotonCounterT0Matching::fTree
private

Definition at line 145 of file PhotonCounterT0Matching_module.cc.

std::string lbne::PhotonCounterT0Matching::fTruthT0ModuleLabel
private

Definition at line 114 of file PhotonCounterT0Matching_module.cc.

bool lbne::PhotonCounterT0Matching::fVerbosity
private

Definition at line 123 of file PhotonCounterT0Matching_module.cc.

double lbne::PhotonCounterT0Matching::fWeightOfDeltaYZ
private

Definition at line 120 of file PhotonCounterT0Matching_module.cc.

TH2D* lbne::PhotonCounterT0Matching::hdeltaX_deltaYZ
private

Definition at line 149 of file PhotonCounterT0Matching_module.cc.

TH2D* lbne::PhotonCounterT0Matching::hdeltaYZ_Length
private

Definition at line 150 of file PhotonCounterT0Matching_module.cc.

TH2D* lbne::PhotonCounterT0Matching::hFitParam_Length
private

Definition at line 151 of file PhotonCounterT0Matching_module.cc.

TH2D* lbne::PhotonCounterT0Matching::hPhotonT0_MCT0
private

Definition at line 152 of file PhotonCounterT0Matching_module.cc.

TH2D* lbne::PhotonCounterT0Matching::hPredX_PE
private

Definition at line 147 of file PhotonCounterT0Matching_module.cc.

TH2D* lbne::PhotonCounterT0Matching::hPredX_T
private

Definition at line 146 of file PhotonCounterT0Matching_module.cc.

TH2D* lbne::PhotonCounterT0Matching::hPredX_T_PE
private

Definition at line 148 of file PhotonCounterT0Matching_module.cc.

TH1D* lbne::PhotonCounterT0Matching::hT0_diff_full
private

Definition at line 153 of file PhotonCounterT0Matching_module.cc.

TH1D* lbne::PhotonCounterT0Matching::hT0_diff_zoom
private

Definition at line 154 of file PhotonCounterT0Matching_module.cc.

double lbne::PhotonCounterT0Matching::MCTruthT0
private

Definition at line 143 of file PhotonCounterT0Matching_module.cc.

double lbne::PhotonCounterT0Matching::minYZSep
private

Definition at line 136 of file PhotonCounterT0Matching_module.cc.

double lbne::PhotonCounterT0Matching::PredictedX
private

Definition at line 133 of file PhotonCounterT0Matching_module.cc.

double lbne::PhotonCounterT0Matching::TimeSep
private

Definition at line 139 of file PhotonCounterT0Matching_module.cc.

double lbne::PhotonCounterT0Matching::TimeSepPredX
private

Definition at line 134 of file PhotonCounterT0Matching_module.cc.

double lbne::PhotonCounterT0Matching::TrackCentre_X
private

Definition at line 126 of file PhotonCounterT0Matching_module.cc.

double lbne::PhotonCounterT0Matching::TrackCentre_Y
private

Definition at line 127 of file PhotonCounterT0Matching_module.cc.

double lbne::PhotonCounterT0Matching::TrackCentre_Z
private

Definition at line 128 of file PhotonCounterT0Matching_module.cc.

double lbne::PhotonCounterT0Matching::TrackLength
private

Definition at line 132 of file PhotonCounterT0Matching_module.cc.

double lbne::PhotonCounterT0Matching::TrackLength_X
private

Definition at line 126 of file PhotonCounterT0Matching_module.cc.

double lbne::PhotonCounterT0Matching::TrackLength_Y
private

Definition at line 127 of file PhotonCounterT0Matching_module.cc.

double lbne::PhotonCounterT0Matching::TrackLength_Z
private

Definition at line 128 of file PhotonCounterT0Matching_module.cc.

double lbne::PhotonCounterT0Matching::trkTimeCentre
private

Definition at line 131 of file PhotonCounterT0Matching_module.cc.

double lbne::PhotonCounterT0Matching::trkTimeEnd
private

Definition at line 129 of file PhotonCounterT0Matching_module.cc.

double lbne::PhotonCounterT0Matching::trkTimeLengh
private

Definition at line 129 of file PhotonCounterT0Matching_module.cc.

double lbne::PhotonCounterT0Matching::trkTimeStart
private

Definition at line 129 of file PhotonCounterT0Matching_module.cc.

double lbne::PhotonCounterT0Matching::YZSep
private

Definition at line 143 of file PhotonCounterT0Matching_module.cc.


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