LArSoft  v07_13_02
Liquid Argon Software toolkit - http://larsoft.org/
opdet::OpFlashAnaAlg Class Reference

#include "OpFlashAnaAlg.h"

Classes

struct  FlashAnaStruct
 
struct  HitAnaStruct
 

Public Member Functions

 OpFlashAnaAlg ()
 
void SetOpFlashTree (TTree *, bool makeOpDetPEHist=true)
 
void SetOpHitTree (TTree *)
 
void FillOpFlashes (const std::vector< recob::OpFlash > &)
 
void FillOpHits (const std::vector< recob::OpHit > &)
 

Private Member Functions

void FillOpFlashTree (const std::vector< recob::OpFlash > &)
 
void FillOpHitTree (const std::vector< recob::OpHit > &)
 

Private Attributes

FlashAnaStruct fOpFlashAnaStruct
 
bool fMakeOpDetPEHist
 
TTree * fOpFlashTree
 
HitAnaStruct fOpHitAnaStruct
 
TTree * fOpHitTree
 

Detailed Description

Definition at line 24 of file OpFlashAnaAlg.h.

Constructor & Destructor Documentation

opdet::OpFlashAnaAlg::OpFlashAnaAlg ( )
inline

Member Function Documentation

void opdet::OpFlashAnaAlg::FillOpFlashes ( const std::vector< recob::OpFlash > &  flashVector)

Definition at line 28 of file OpFlashAnaAlg.cxx.

References FillOpFlashTree(), and fOpFlashTree.

Referenced by opdet::OpFlashSimpleAna::analyze(), and OpFlashAnaAlg().

29 {
30  if(fOpFlashTree) FillOpFlashTree(flashVector);
31 }
void FillOpFlashTree(const std::vector< recob::OpFlash > &)
void opdet::OpFlashAnaAlg::FillOpFlashTree ( const std::vector< recob::OpFlash > &  flashVector)
private

Definition at line 38 of file OpFlashAnaAlg.cxx.

References opdet::OpFlashAnaAlg::FlashAnaStruct::FlashAbsTime, opdet::OpFlashAnaAlg::FlashAnaStruct::FlashFastToTotal, opdet::OpFlashAnaAlg::FlashAnaStruct::FlashFrame, opdet::OpFlashAnaAlg::FlashAnaStruct::FlashInBeamFrame, opdet::OpFlashAnaAlg::FlashAnaStruct::FlashOnBeamTime, opdet::OpFlashAnaAlg::FlashAnaStruct::FlashOpDetPEHist, opdet::OpFlashAnaAlg::FlashAnaStruct::FlashTime, opdet::OpFlashAnaAlg::FlashAnaStruct::FlashTimeWidth, opdet::OpFlashAnaAlg::FlashAnaStruct::FlashTotalPE, opdet::OpFlashAnaAlg::FlashAnaStruct::FlashY, opdet::OpFlashAnaAlg::FlashAnaStruct::FlashYWidth, opdet::OpFlashAnaAlg::FlashAnaStruct::FlashZ, opdet::OpFlashAnaAlg::FlashAnaStruct::FlashZWidth, fMakeOpDetPEHist, fOpFlashAnaStruct, and fOpFlashTree.

Referenced by FillOpFlashes().

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 }
FlashAnaStruct fOpFlashAnaStruct
Definition: OpFlashAnaAlg.h:60
void opdet::OpFlashAnaAlg::FillOpHits ( const std::vector< recob::OpHit > &  hitVector)

Definition at line 33 of file OpFlashAnaAlg.cxx.

References FillOpHitTree(), and fOpHitTree.

Referenced by opdet::OpFlashSimpleAna::analyze(), and OpFlashAnaAlg().

34 {
35  if(fOpHitTree) FillOpHitTree(hitVector);
36 }
void FillOpHitTree(const std::vector< recob::OpHit > &)
void opdet::OpFlashAnaAlg::FillOpHitTree ( const std::vector< recob::OpHit > &  hitVector)
private

Definition at line 66 of file OpFlashAnaAlg.cxx.

References fOpHitAnaStruct, fOpHitTree, opdet::OpFlashAnaAlg::HitAnaStruct::HitAmplitude, opdet::OpFlashAnaAlg::HitAnaStruct::HitArea, opdet::OpFlashAnaAlg::HitAnaStruct::HitFastToTotal, opdet::OpFlashAnaAlg::HitAnaStruct::HitFrame, opdet::OpFlashAnaAlg::HitAnaStruct::HitOpChannel, opdet::OpFlashAnaAlg::HitAnaStruct::HitPE, opdet::OpFlashAnaAlg::HitAnaStruct::HitPeakTime, opdet::OpFlashAnaAlg::HitAnaStruct::HitPeakTimeAbs, and opdet::OpFlashAnaAlg::HitAnaStruct::HitWidth.

Referenced by FillOpHits().

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 }
HitAnaStruct fOpHitAnaStruct
Definition: OpFlashAnaAlg.h:83
Detector simulation of raw signals on wires.
void opdet::OpFlashAnaAlg::SetOpFlashTree ( TTree *  tree,
bool  makeOpDetPEHist = true 
)

Definition at line 4 of file OpFlashAnaAlg.cxx.

References opdet::OpFlashAnaAlg::FlashAnaStruct::FlashOpDetPEHist, fMakeOpDetPEHist, fOpFlashAnaStruct, fOpFlashTree, and opdet::OpFlashAnaAlg::FlashAnaStruct::LeafList.

Referenced by opdet::OpFlashSimpleAna::beginJob(), and OpFlashAnaAlg().

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 }
FlashAnaStruct fOpFlashAnaStruct
Definition: OpFlashAnaAlg.h:60
void opdet::OpFlashAnaAlg::SetOpHitTree ( TTree *  tree)

Definition at line 20 of file OpFlashAnaAlg.cxx.

References fOpHitAnaStruct, fOpHitTree, and opdet::OpFlashAnaAlg::HitAnaStruct::LeafList.

Referenced by opdet::OpFlashSimpleAna::beginJob(), and OpFlashAnaAlg().

21 {
22  fOpHitTree = tree;
23  fOpHitTree->Branch(fOpHitTree->GetName(),
25  fOpHitAnaStruct.LeafList.c_str());
26 }
HitAnaStruct fOpHitAnaStruct
Definition: OpFlashAnaAlg.h:83

Member Data Documentation

bool opdet::OpFlashAnaAlg::fMakeOpDetPEHist
private

Definition at line 61 of file OpFlashAnaAlg.h.

Referenced by FillOpFlashTree(), OpFlashAnaAlg(), and SetOpFlashTree().

FlashAnaStruct opdet::OpFlashAnaAlg::fOpFlashAnaStruct
private

Definition at line 60 of file OpFlashAnaAlg.h.

Referenced by FillOpFlashTree(), and SetOpFlashTree().

TTree* opdet::OpFlashAnaAlg::fOpFlashTree
private

Definition at line 63 of file OpFlashAnaAlg.h.

Referenced by FillOpFlashes(), FillOpFlashTree(), and SetOpFlashTree().

HitAnaStruct opdet::OpFlashAnaAlg::fOpHitAnaStruct
private

Definition at line 83 of file OpFlashAnaAlg.h.

Referenced by FillOpHitTree(), and SetOpHitTree().

TTree* opdet::OpFlashAnaAlg::fOpHitTree
private

Definition at line 85 of file OpFlashAnaAlg.h.

Referenced by FillOpHits(), FillOpHitTree(), and SetOpHitTree().


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