All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
BernCRTAna_module.cc
Go to the documentation of this file.
1 #include "art/Framework/Core/EDAnalyzer.h"
2 #include "art/Framework/Core/ModuleMacros.h"
3 #include "art/Framework/Principal/Event.h"
4 #include "art/Framework/Principal/Handle.h"
5 
6 #include "canvas/Utilities/Exception.h"
7 
8 #include "sbndaq-artdaq-core/Overlays/Common/BernCRTFragment.hh"
9 #include "artdaq-core/Data/Fragment.hh"
10 #include "artdaq-core/Data/ContainerFragment.hh"
11 #include "sbndaq-artdaq-core/Overlays/FragmentType.hh"
12 
13 #include "sbndaq-artdaq-core/Overlays/Common/BernCRTTranslator.hh"
14 
15 //#include "art/Framework/Services/Optional/TFileService.h"
16 #include "art_root_io/TFileService.h"
17 #include "TH1F.h"
18 #include "TNtuple.h"
19 
20 #include <algorithm>
21 #include <cassert>
22 #include <cmath>
23 #include <fstream>
24 #include <iomanip>
25 #include <vector>
26 #include <iostream>
27 
28 namespace sbndaq {
29  class BernCRTAna;
30 }
31 
32 class sbndaq::BernCRTAna : public art::EDAnalyzer {
33 
34 public:
35  explicit BernCRTAna(fhicl::ParameterSet const & pset);
36  virtual ~BernCRTAna();
37 
38  virtual void analyze(art::Event const & evt);
39 
40 
41 private:
42  TTree * hits;
43  bool IsSideCRT(icarus::crt::BernCRTTranslator & hit);
44  bool IsTopCRT(icarus::crt::BernCRTTranslator & hit);
45 //data payload
46  uint8_t mac5; //last 8 bits of FEB mac5 address
47  uint16_t flags;
48  uint16_t lostcpu;
49  uint16_t lostfpga;
50  uint32_t ts0;
51  uint32_t ts1;
52  uint16_t adc[32];
53  uint32_t coinc;
54  int subSys; // Top CRT 0, Side CRT 1, Bottom 2
55 
56 
57 //metadata
58  uint64_t run_start_time;
59  uint64_t this_poll_start;
60  uint64_t this_poll_end;
61  uint64_t last_poll_start;
62  uint64_t last_poll_end;
64  uint32_t feb_hits_in_poll;
65  uint32_t feb_hit_number;
66 
67 //information from fragment header
68  uint32_t sequence_id;
71  uint16_t lost_hits;
72  uint16_t hits_in_fragment;
73 };
74 
75 //Define the constructor
76 sbndaq::BernCRTAna::BernCRTAna(fhicl::ParameterSet const & pset)
77  : EDAnalyzer(pset)
78 {
79 
80  art::ServiceHandle<art::TFileService> tfs; //pointer to a file named tfs
81 
82  hits = tfs->make<TTree>("hits", "FEB hits");
83 
84  hits->Branch("mac5", &mac5, "mac5/b");
85  hits->Branch("flags", &flags, "flags/s");
86  hits->Branch("lostcpu", &lostcpu, "lostcpu/s");
87  hits->Branch("lostfpga", &lostfpga, "lostfpga/s");
88  hits->Branch("ts0", &ts0, "ts0/i");
89  hits->Branch("ts1", &ts1, "ts1/i");
90  hits->Branch("adc", &adc, "adc[32]/s");
91  hits->Branch("coinc", &coinc, "coinc/i");
92  hits->Branch("subSys", &subSys, "subSys/i");
93 
94  hits->Branch("run_start_time", &run_start_time, "run_start_time/l");
95  hits->Branch("this_poll_start", &this_poll_start, "this_poll_start/l");
96  hits->Branch("this_poll_end", &this_poll_end, "this_poll_end/l");
97  hits->Branch("last_poll_start", &last_poll_start, "last_poll_start/l");
98  hits->Branch("last_poll_end", &last_poll_end, "last_poll_end/l");
99  hits->Branch("system_clock_deviation", &system_clock_deviation, "system_clock_deviation/I");
100  hits->Branch("feb_hits_in_poll", &feb_hits_in_poll, "feb_hits_in_poll/i");
101  hits->Branch("feb_hit_number", &feb_hit_number, "feb_hit_number/i");
102 
103  hits->Branch("sequence_id", &sequence_id, "sequence_id/i");
104  hits->Branch("fragment_timestamp", &fragment_timestamp, "fragment_timestamp/l");
105 
106  hits->Branch("last_accepted_timestamp", &last_accepted_timestamp, "last_accepted_timestamp/l");
107  hits->Branch("lost_hits", &lost_hits, "lost_hits/s");
108  hits->Branch("hits_in_fragment", &hits_in_fragment, "hits_in_fragment/s");
109 
110 }
111 
112 bool sbndaq::BernCRTAna::IsSideCRT(icarus::crt::BernCRTTranslator & hit) {
113  /**
114  * * Fragment ID described in SBN doc 16111
115  * */
116  return (hit.fragment_ID & 0x3100) == 0x3100;
117 }
118 bool sbndaq::BernCRTAna::IsTopCRT(icarus::crt::BernCRTTranslator & hit) {
119  /**
120  * * * Fragment ID described in SBN doc 16111
121  * * */
122  return (hit.fragment_ID & 0x3200) == 0x3200;
123 }
124 
126 {
127 }
128 
129 void sbndaq::BernCRTAna::analyze(art::Event const & evt) {
130 
131  //WK 09/02/21. Update to BernCRTTranslator in sbndaq_artdaq_core
132  std::vector<icarus::crt::BernCRTTranslator> hit_vector;
133 
134  auto fragmentHandles = evt.getMany<artdaq::Fragments>();
135  for (auto handle : fragmentHandles) {
136  if (!handle.isValid() || handle->size() == 0)
137  continue;
138 
139  auto this_hit_vector = icarus::crt::BernCRTTranslator::getCRTData(*handle);
140 
141  hit_vector.insert(hit_vector.end(),this_hit_vector.begin(),this_hit_vector.end());
142 
143  }
144 
145  for(auto & hit : hit_vector) {
146  // TLOG(TLVL_INFO)<<hit;
147 
148  if(IsSideCRT(hit)){
149  subSys = 1;
150  }
151  else if(IsTopCRT(hit)){
152  subSys = 0;
153  }
154  else{
155  subSys = 2;
156  }
157  fragment_timestamp = hit.timestamp;
158  sequence_id = hit.sequence_id;
159 
160  mac5 = hit.mac5;
161  flags = hit.flags;
162  lostcpu = hit.lostcpu;
163  lostfpga = hit.lostfpga;
164  ts0 = hit.ts0;
165  ts1 = hit.ts1;
166  coinc = hit.coinc;
167 
168  for(int ch=0; ch<32; ch++) adc[ch] = hit.adc[ch];
169 
170  run_start_time = hit.run_start_time;
171  this_poll_start = hit.this_poll_start;
172  this_poll_end = hit.this_poll_end;
173  last_poll_start = hit.last_poll_start;
174  last_poll_end = hit.last_poll_end;
175  system_clock_deviation = hit.system_clock_deviation;
176  feb_hits_in_poll = hit.hits_in_poll;
177  feb_hit_number = hit.feb_hit_number;
178 
179  last_accepted_timestamp = hit.last_accepted_timestamp;
180  lost_hits = hit.lost_hits;
181  hits_in_fragment = hit.hits_in_fragment;
182 
183  hits->Fill();
184  }
185 } //analyze
186 
187 
188 DEFINE_ART_MODULE(sbndaq::BernCRTAna)
189 
BernCRTAna(fhicl::ParameterSet const &pset)
process_name hit
Definition: cheaterreco.fcl:51
bool IsTopCRT(icarus::crt::BernCRTTranslator &hit)
bool IsSideCRT(icarus::crt::BernCRTTranslator &hit)
art::ServiceHandle< art::TFileService > tfs
TCEvent evt
Definition: DataStructs.cxx:8
virtual void analyze(art::Event const &evt)
uint64_t last_accepted_timestamp