LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
OpFlashAnaAlg.cxx
Go to the documentation of this file.
2 
3 #include "TTree.h"
4 
5 void opdet::OpFlashAnaAlg::SetOpFlashTree(TTree* tree, bool makeOpDetPEHist)
6 {
7  fOpFlashTree = tree;
8  fMakeOpDetPEHist = makeOpDetPEHist;
9 
10  fOpFlashTree->Branch(
12  if (fMakeOpDetPEHist) {
13  fOpFlashAnaStruct.FlashOpDetPEHist->SetName("hOpDetPEs");
14  fOpFlashAnaStruct.FlashOpDetPEHist->SetTitle("PEs per Optical Detector; OpDet Number; PE");
15  fOpFlashTree->Branch("flash_OpDetPEHist", "TH1D", &(fOpFlashAnaStruct.FlashOpDetPEHist));
16  }
17 }
18 
20 {
21  fOpHitTree = tree;
22  fOpHitTree->Branch(fOpHitTree->GetName(), &fOpHitAnaStruct, fOpHitAnaStruct.LeafList.c_str());
23 }
24 
25 void opdet::OpFlashAnaAlg::FillOpFlashes(const std::vector<recob::OpFlash>& flashVector)
26 {
27  if (fOpFlashTree) FillOpFlashTree(flashVector);
28 }
29 
30 void opdet::OpFlashAnaAlg::FillOpHits(const std::vector<recob::OpHit>& hitVector)
31 {
32  if (fOpHitTree) FillOpHitTree(hitVector);
33 }
34 
35 void opdet::OpFlashAnaAlg::FillOpFlashTree(const std::vector<recob::OpFlash>& flashVector)
36 {
37  for (auto const& flash : flashVector) {
38 
39  fOpFlashAnaStruct.FlashTime = flash.Time();
40  fOpFlashAnaStruct.FlashTimeWidth = flash.TimeWidth();
41  fOpFlashAnaStruct.FlashAbsTime = flash.AbsTime();
42  fOpFlashAnaStruct.FlashFrame = flash.Frame();
43  fOpFlashAnaStruct.FlashY = flash.YCenter();
44  fOpFlashAnaStruct.FlashYWidth = flash.YWidth();
45  fOpFlashAnaStruct.FlashZ = flash.ZCenter();
46  fOpFlashAnaStruct.FlashZWidth = flash.ZWidth();
47  fOpFlashAnaStruct.FlashInBeamFrame = flash.InBeamFrame();
48  fOpFlashAnaStruct.FlashOnBeamTime = flash.OnBeamTime();
49  fOpFlashAnaStruct.FlashFastToTotal = flash.FastToTotal();
50  fOpFlashAnaStruct.FlashTotalPE = flash.TotalPE();
51 
52  if (fMakeOpDetPEHist) {
53  int n_opdets = flash.WireCenters().size();
54  fOpFlashAnaStruct.FlashOpDetPEHist->SetBins(n_opdets, -0.5, (float)n_opdets - 0.5);
55  for (int i = 0; i < n_opdets; i++)
56  fOpFlashAnaStruct.FlashOpDetPEHist->SetBinContent(i + 1, flash.PE(i));
57  }
58 
59  fOpFlashTree->Fill();
60  }
61 }
62 
63 void opdet::OpFlashAnaAlg::FillOpHitTree(const std::vector<recob::OpHit>& hitVector)
64 {
65  for (auto const& hit : hitVector) {
66 
67  fOpHitAnaStruct.HitPeakTime = hit.PeakTime();
68  fOpHitAnaStruct.HitPeakTimeAbs = hit.PeakTimeAbs();
69  fOpHitAnaStruct.HitWidth = hit.Width();
70  fOpHitAnaStruct.HitArea = hit.Area();
71  fOpHitAnaStruct.HitAmplitude = hit.Amplitude();
72  fOpHitAnaStruct.HitFastToTotal = hit.FastToTotal();
73  fOpHitAnaStruct.HitPE = hit.PE();
74  fOpHitAnaStruct.HitFrame = hit.Frame();
75  fOpHitAnaStruct.HitOpChannel = hit.OpChannel();
76 
77  fOpHitTree->Fill();
78  }
79 }
void FillOpFlashes(const std::vector< recob::OpFlash > &)
void FillOpHits(const std::vector< recob::OpHit > &)
void SetOpFlashTree(TTree *, bool makeOpDetPEHist=true)
FlashAnaStruct fOpFlashAnaStruct
Definition: OpFlashAnaAlg.h:61
void SetOpHitTree(TTree *)
HitAnaStruct fOpHitAnaStruct
Definition: OpFlashAnaAlg.h:84
Detector simulation of raw signals on wires.
void FillOpFlashTree(const std::vector< recob::OpFlash > &)
void FillOpHitTree(const std::vector< recob::OpHit > &)