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