LArSoft  v10_04_05
Liquid Argon Software toolkit - https://larsoft.org/
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
OpHitAna_module.cc
Go to the documentation of this file.
1 // Christie Chiu and Ben Jones, MIT, 2012
2 //
3 // This is an analyzer module which writes the raw optical
4 // detector pulses for each PMT to an output file
5 //
6 
7 // LArSoft includes
11 
12 // Framework includes
18 #include "art_root_io/TFileService.h"
20 #include "fhiclcpp/ParameterSet.h"
22 
23 // ROOT includes
24 #include "TH1D.h"
25 #include "TLorentzVector.h"
26 #include "TTree.h"
27 
28 // C++ Includes
29 #include "math.h"
30 
31 namespace opdet {
32 
33  class OpHitAna : public art::EDAnalyzer {
34  public:
35  // Standard constructor and destructor for an ART module.
37 
38  // The analyzer routine, called once per event.
39  void analyze(const art::Event&);
40 
41  private:
42  // The stuff below is the part you'll most likely have to change to
43  // go from this custom example to your own task.
44 
45  // The parameters we'll read from the .fcl file.
46  std::string fInputModule; // Input tag for OpDet collection
47  float fSampleFreq; // in MHz
48  float fTimeBegin; // in us
49  float fTimeEnd; // in us
50  // short fPEheight; // in ADC counts
51 
52  // Flags to enable or disable output of debugging TH1 / TH2s
56 
57  // Output TTree and its branch variables
58  TTree* fHitTree;
59  Int_t fEventID;
60  Int_t fOpChannel;
61  Float_t fPeakTime;
62  Float_t fNPe;
63  };
64 
65 }
66 
67 // OpHitAna_module.cc
68 
69 // This is a short program required by the ART framework. It enables
70 // a program (OpHitAna, in this case) to be called as a module
71 // from a .fcl file. It is unlikely that you'll have to make any
72 // changes to this file.
73 
74 namespace opdet {
76 }
77 
78 namespace opdet {
79 
80  //-----------------------------------------------------------------------
81  // Constructor
83  {
84 
85  // Indicate that the Input Module comes from .fcl
86  fInputModule = pset.get<std::string>("InputModule");
87 
89  fTimeBegin = odp->TimeBegin();
90  fTimeEnd = odp->TimeEnd();
91  fSampleFreq = odp->SampleFreq();
92 
93  fMakeHistPerChannel = pset.get<bool>("MakeHistPerChannel");
94  fMakeHistAllChannels = pset.get<bool>("MakeHistAllChannels");
95  fMakeHitTree = pset.get<bool>("MakeHitTree");
96 
97  // If required, make TTree for output
98 
99  if (fMakeHitTree) {
101  fHitTree = tfs->make<TTree>("HitTree", "HitTree");
102  fHitTree->Branch("EventID", &fEventID, "EventID/I");
103  fHitTree->Branch("OpChannel", &fOpChannel, "OpChannel/I");
104  fHitTree->Branch("PeakTime", &fPeakTime, "PeakTime/F");
105  fHitTree->Branch("NPe", &fNPe, "NPe/F");
106  }
107  }
108 
109  //-----------------------------------------------------------------------
111  {
112  // Create a handle for our vector of pulses
114 
115  // Create string for histogram name
116  char HistName[50];
117 
118  // Read in HitHandle
119  evt.getByLabel(fInputModule, HitHandle);
120 
121  // Access ART's TFileService, which will handle creating and writing
122  // histograms for us.
124 
125  int NOpChannels = art::ServiceHandle<geo::WireReadout const>()->Get().NOpChannels();
126 
127  std::vector<TH1D*> HitHist;
128 
129  sprintf(HistName, "Event %d AllOpDets", evt.id().event());
130 
131  TH1D* AllHits = nullptr;
132  if (fMakeHistAllChannels) {
133  AllHits = tfs->make<TH1D>(HistName,
134  ";t (ns);",
135  int((fTimeEnd - fTimeBegin) * fSampleFreq),
136  fTimeBegin * 1000.,
137  fTimeEnd * 1000.);
138  }
139 
140  for (int i = 0; i != NOpChannels; ++i) {
141 
142  sprintf(HistName, "Event %d OpDet %i", evt.id().event(), i);
144  HitHist.push_back(tfs->make<TH1D>(HistName,
145  ";t (ns);",
146  int((fTimeEnd - fTimeBegin) * fSampleFreq),
147  fTimeBegin * 1000.,
148  fTimeEnd * 1000.));
149  }
150 
151  fEventID = evt.id().event();
152 
153  // For every OpHit in the vector
154  for (unsigned int i = 0; i < HitHandle->size(); ++i) {
155  // Get OpHit
156  art::Ptr<recob::OpHit> TheHitPtr(HitHandle, i);
157  recob::OpHit TheHit = *TheHitPtr;
158 
159  fOpChannel = TheHit.OpChannel();
160  fNPe = TheHit.PE();
161  fPeakTime = TheHit.PeakTime();
162 
163  if (fMakeHitTree) fHitTree->Fill();
164 
165  if (fMakeHistPerChannel) {
166  if (fOpChannel > int(HitHist.size())) {
167  mf::LogError("OpHitAna") << "Error : Trying to fill channel out of range, " << fOpChannel;
168  }
169  HitHist[fOpChannel]->Fill(fPeakTime, fNPe);
170  }
171 
172  if (fMakeHistAllChannels) AllHits->Fill(fPeakTime, fNPe);
173  }
174  }
175 
176 } // namespace opdet
double PeakTime() const
Definition: OpHit.h:94
cout<< "Opened file "<< fin<< " ixs= "<< ixs<< endl;if(ixs==0) hhh=(TH1F *) fff-> Get("h1")
Definition: AddMC.C:8
MaybeLogger_< ELseverityLevel::ELsev_error, false > LogError
EDAnalyzer(fhicl::ParameterSet const &pset)
Definition: EDAnalyzer.cc:6
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:65
T get(std::string const &key) const
Definition: ParameterSet.h:314
double PE() const
Definition: OpHit.h:122
bool getByLabel(std::string const &label, std::string const &instance, Handle< PROD > &result) const
std::string fInputModule
void analyze(const art::Event &)
EventNumber_t event() const
Definition: EventID.h:116
TCEvent evt
Definition: DataStructs.cxx:8
int OpChannel() const
Definition: OpHit.h:86
EventID id() const
Definition: Event.cc:23
OpHitAna(const fhicl::ParameterSet &)