LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
FIFOHistogramAna_module.cc
Go to the documentation of this file.
1 // Ben Jones, MIT, 2013
2 //
3 // This simple ana module writes a root file of pulse histograms,
4 // organized by event and channel number
5 //
6 
7 // LArSoft includes
9 
10 // Framework includes
16 #include "art_root_io/TFileDirectory.h"
17 #include "art_root_io/TFileService.h"
18 #include "fhiclcpp/ParameterSet.h"
19 
20 // ROOT includes
21 #include "TH1.h"
22 
23 // C++ Includes
24 #include <cstring>
25 #include <map>
26 #include <sstream>
27 
28 namespace opdet {
29 
31  public:
32  // Standard constructor and destructor for an ART module.
34 
35  // The analyzer routine, called once per event.
36  void analyze(const art::Event&);
37 
38  private:
39  // The parameters we'll read from the .fcl file.
40  std::string fInputModule; // Input tag for OpDet collection
41  };
42 
43 }
44 
45 namespace opdet {
46 
47  //-----------------------------------------------------------------------
48  // Constructor
50  {
51  // Indicate that the Input Module comes from .fcl
52  fInputModule = pset.get<std::string>("InputModule");
53  }
54 
55  //-----------------------------------------------------------------------
57  {
58 
60 
61  // Create a handle for our vector of pulses
63 
64  // Read in WaveformHandle
65  evt.getByLabel(fInputModule, FIFOChannelHandle);
66 
67  int Run = evt.run();
68  int EID = evt.event();
69 
70  std::stringstream FolderName;
71  FolderName.flush();
72  FolderName << "run" << Run << "_evt" << EID;
73 
74  art::TFileDirectory evtfolder = tfs->mkdir(FolderName.str().c_str());
75 
76  std::map<int, bool> ChanFolderMade;
77  std::map<uint32_t, int> ChanFolderIndex;
78  std::vector<art::TFileDirectory> ChanFolders;
79 
80  for (size_t i = 0; i != FIFOChannelHandle->size(); ++i) {
81  uint32_t Frame = FIFOChannelHandle->at(i).Frame();
82  uint32_t TimeSlice = FIFOChannelHandle->at(i).TimeSlice();
83  uint32_t Channel = FIFOChannelHandle->at(i).ChannelNumber();
84 
85  if (!ChanFolderMade[Channel]) {
86  std::stringstream ChannelLabel;
87  ChannelLabel.flush();
88  ChannelLabel << "chan" << Channel;
89  ChanFolderIndex[Channel] = ChanFolders.size();
90  ChanFolders.push_back(evtfolder.mkdir(ChannelLabel.str().c_str()));
91  ChanFolderMade[Channel] = true;
92  }
93 
94  std::stringstream HistName;
95  HistName.flush();
96  HistName << "frm" << Frame << "_"
97  << "tsl" << TimeSlice;
98 
99  TH1D* ThisHist = ChanFolders[ChanFolderIndex[Channel]].make<TH1D>(
100  HistName.str().c_str(),
101  HistName.str().c_str(),
102  FIFOChannelHandle->at(i).size(),
103  float(TimeSlice) - 0.0001,
104  float(FIFOChannelHandle->at(i).size()) - 0.0001 + TimeSlice);
105 
106  for (size_t j = 0; j != FIFOChannelHandle->at(i).size(); ++j) {
107  ThisHist->Fill(TimeSlice + j, FIFOChannelHandle->at(i).at(j));
108  }
109  }
110  }
111 
112 } // namespace opdet
113 
114 namespace opdet {
116 }
void analyze(const art::Event &)
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
EventNumber_t event() const
Definition: Event.cc:41
bool getByLabel(std::string const &label, std::string const &instance, Handle< PROD > &result) const
FIFOHistogramAna(const fhicl::ParameterSet &)
TCEvent evt
Definition: DataStructs.cxx:8
RunNumber_t run() const
Definition: Event.cc:29