Jpp  18.3.0-209-g56ce19a
the software that should make you happy
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
JRunAnalyzer.hh
Go to the documentation of this file.
1 #ifndef __JRUNANALYZER__
2 #define __JRUNANALYZER__
3 
4 #include "JSupport/JSupport.hh"
9 #include "JROOT/JManager.hh"
10 #include "JTools/JQuantile.hh"
11 #include "JTools/JRange.hh"
12 #include "JSupport/JTreeScanner.hh"
15 #include "JMath/JConstants.hh"
16 #include "JTrigger/JHitL0.hh"
17 #include "JTrigger/JBuildL0.hh"
18 #include "TFile.h"
19 #include "TH1D.h"
20 #include "JRunHistograms.hh"
21 
22 /**
23  * \author rgruiz, adomi
24  */
25 
27 
28 /**
29  * Class dedicated to the analysis of %KM3NeT runs.
30  *
31  */
32 class JRunAnalyzer {
33 
36  TFile *outputFile;
42 
43 public :
44 
45  /**
46  * Constructor.
47  *
48  * \param file file name
49  * \param detector detector
50  * \param out pointer to output file
51  * \param nTimeslices number of timeslices to be read
52  * \param nSummaryslices number of summary slices to be read
53  * \param nEvents number of events to be read
54  * \param analysislevel option for analysis of PMT or reco data
55  */
57  JLimit_t nTimeslices, JLimit_t nSummaryslices, JLimit_t nEvents, int analysislevel):
58  inputFile (file),
59  router (detector),
60  outputFile (out),
61  numberOfTimeslices (nTimeslices),
62  numberOfSummaryslices (nSummaryslices),
63  numberOfEvents (nEvents),
64  analysis_level (analysislevel)
65  {
66  using namespace JPP;
67  using namespace std;
68 
69  histograms = JRA_Histograms(detector);
70  }
71 
72  /**
73  * Destructor.
74  */
76 
77  /*
78  * Function template used to read events from a tree.
79  *
80  * \param scanner A JTreeScanner
81  * \param
82  */
84 
85  while (scanner.hasNext()) {
86 
87  JDAQEvent event = *(scanner.next());
88 
89  histograms.h_trigger.h_Snapshot_hits -> Fill((Double_t) event.size<JDAQSnapshotHit> ());
90  histograms.h_trigger.h_Triggered_hits -> Fill((Double_t) event.size<JDAQTriggeredHit> ());
91  histograms.h_trigger.h_Triggered_over_Snapshot_hits -> Fill((Double_t) event.size<JDAQTriggeredHit> () / event.size<JDAQSnapshotHit > ());
92  histograms.h_trigger.h_Number_of_overlays -> Fill(event.getOverlays());
93 
94  for (unsigned int i = 0; i != NUMBER_OF_TRIGGER_BITS; ++i) {
95 
96  if (event.hasTriggerBit(i)) {
97 
98  histograms.h_trigger.h_Trigger_bit_event -> Fill((Double_t) i);
99  }
100  }
101 
102  int counter_3dmuon = 0;
103  typedef JRange<double> JRange_t;
104  JRange_t range(JRange_t::DEFAULT_RANGE());
105 
106  for (JDAQEvent::const_iterator<JDAQTriggeredHit> hit = event.begin<JDAQTriggeredHit>(); hit != event.end<JDAQTriggeredHit>(); ++hit) {
107 
108  const JModule& module = router.getModule(hit->getModuleID());
109  const double t1 = getTime(hit->getT(), module.getPMT(hit->getPMT()));
110 
111  range.include(t1);
112 
113  JDAQTriggerMask trigger_mask(event.getTriggerMask(*hit));
114 
115  if(trigger_mask.hasTriggerBit(getTriggerBit<JTrigger3DMuon>())) counter_3dmuon++;
116 
117  if (router.hasModule(hit->getModuleID())) {
118 
119  const JModule& module = router.getModule (hit->getModuleID());
120 
123 
124  int String = module.getString();
125  int Floor = module.getFloor();
126 
127  histograms.h_trigger.h_Triggered_hits_per_module -> Fill(String , Floor);
128 
129  if(trigger_mask.hasTriggerBit(getTriggerBit<JTrigger3DMuon>())){
131  }
132 
133  for (unsigned int i = 0; i != NUMBER_OF_TRIGGER_BITS; ++i) {
134 
135  if (hit -> hasTriggerBit(i)) {
136 
137  histograms.h_trigger.h_Trigger_bit_hit->Fill((Double_t) i);
138  }
139  }
140  }else{
141  FATAL("JModuleRouter trying to access non existing identifier: "<< hit->getModuleID());
142  }
143  }
144 
146 
147  histograms.h_trigger.h_pmt_distribution_triggered_hits->Scale(1. / (Double_t) event.size<JDAQTriggeredHit> ());
148 
149  if(counter_3dmuon != 0) histograms.h_trigger.h_Triggered_hits_3dmuon->Fill(counter_3dmuon);
150 
151  for (JDAQEvent::const_iterator<JDAQSnapshotHit> hit = event.begin<JDAQSnapshotHit>() ; hit != event.end<JDAQSnapshotHit>() ; ++hit){
152 
153  if (router.hasModule(hit->getModuleID())) {
154 
155  const JModule& module = router.getModule (hit->getModuleID());
156 
157  int String = module.getString();
158  int Floor = module.getFloor();
159  int pmt = hit-> getPMT();
160 
161  if(analysis_level == 1){
162 
163  (*histograms.h_trigger.m_Snapshot_hits_per_pmt)[MAKE_STRING("Detector/DU" + to_string(String))] -> Fill(pmt, Floor);
164 
165  }
166 
167  histograms.h_trigger.h_pmt_distribution_snapshot_hits -> Fill(hit->getPMT());
169  histograms.h_trigger.h_Snapshot_hits_per_module -> Fill(String, Floor);
170 
171  }else{
172  FATAL("JModuleRouter trying to access non existing identifier: "<< hit->getModuleID());
173  }
174  }
175  }
176  }
177 
178  /*
179  * Function template used to read reconstructed events from a tree.
180  *
181  * \param scanner A ParallelFileScanner for JDAQEvent and JEvt
182  * \param
183  */
185 
186  TH1D* h_tres = new TH1D("h_tres", ";Time residuals [ns]; Entries", 100, -50, 150);
187  TH1D* h_likelihood = new TH1D ("h_likelihood", " ; Likelihood; Reconstructed Events ", 100, -1000, 1000);
188  TH1D* h_beta0 = new TH1D("h_beta0", "; beta0; Reconstructed Events", 20, 0, 1);
189  TH1D* h_energy = new TH1D ("h_energy", " ; Energy [GeV]; Reconstructed Events ", 65, 0, 9);
190  setLogarithmicX(h_energy);
191  TH1D* h_zenith = new TH1D("h_zenith", "; cosZenith; Reconstructed Events", 20, -1, 1);
192  TH1D* h_azimuth = new TH1D("h_azimuth", "; cosAzimuth; Reconstructed Events", 20, -1, 1);
193  TH1D* h_radial_position = new TH1D ("h_radial_position", "; Radial Position [m]; Reconstructed Events", 60, 0, 4.2);
194  setLogarithmicX(h_radial_position);
195  TH1D* h_z_position = new TH1D ("h_z_position", "; Z Position [m] ; Reconstructed Events", 50, 0, 3.2);
196  setLogarithmicX(h_z_position);
197 
198  while (scanner.hasNext()) {
199 
200  JParallelFileScanner< JTypeList<JDAQEvent, JFIT::JEvt> >::multi_pointer_type ps = scanner.next();
201 
202  JEvt* evt = ps;
203  JDAQEvent* event = ps;
204 
205  if (!evt->empty()) {
206 
207  JEvt::iterator best = evt->begin();
208  h_energy -> Fill(best->getE());
209  h_likelihood -> Fill(best->getQ());
210  h_z_position -> Fill(best->getZ());
211  h_radial_position -> Fill(sqrt(best->getX()*best->getX() + best->getY()*best->getY()));
212 
213  JTrack3D track = getTrack(*best);
214  h_zenith -> Fill(cos(track.getDirection().getTheta()));
215  h_azimuth -> Fill(cos(track.getDirection().getPhi()));
216  h_beta0 -> Fill(best->getW(JGANDALF_BETA0_RAD));
217 
218  std::vector<JHitL0> dataL0;
219  const JBuildL0<JHitL0> buildL0;
220  buildL0(JDAQTimeslice(*event, false), router, back_inserter(dataL0));
221 
222  for (std::vector<JHitL0>::const_iterator hit = dataL0.begin(); hit != dataL0.end(); ++hit) {
223 
224  h_tres -> Fill(hit->getT() - track.getT(hit->getPosition()));
225  }
226  }
227  }
228 
229  f.mkdir("Reco");
230  f.cd("Reco");
231 
232  h_energy->Write();
233  h_likelihood -> Write();
234  h_z_position -> Write();
235  h_radial_position -> Write();
236  h_zenith -> Write();
237  h_azimuth -> Write();
238  h_beta0 -> Write();
239  h_tres -> Write();
240  }
241 
242  struct dom_type :
243  public std::vector< JQuantile>
244  {
247  {}
248  };
249 
250  /*
251  * Function template to read and process time ordered summary slices.
252  *
253  * \param scanner A JTreeScanner from which time ordered summary slices are accessed.
254  */
256 
257  using namespace std;
258  using namespace JPP;
259 
260  map<int, map<int, dom_type> > PMT_rate_quantiles;
261  map<int, map<int, JQuantile> > DOM_rate_quantiles;
262 
263  while (scanner.hasNext()){
264 
265  JDAQSummaryslice slice = *(scanner.next());
266 
267  for (JDAQSummaryslice::const_iterator frame = slice.begin() ; frame != slice.end() ; ++ frame) {
268 
269  if (router.hasModule(frame->getModuleID())) {
270 
271  const JModule& module = router.getModule (frame->getModuleID());
272 
273  int string = module.getString();
274  int floor = module.getFloor ();
275 
276  JDAQFrameStatus status = frame -> getDAQFrameStatus();
277  int nFIFOcount = status.countFIFOStatus();
278  int nHRVcount = status.countHighRateVeto();
279 
280  if (nHRVcount > 0) {
281  histograms.h_summary.h_hrv_per_dom->Fill(string , floor, nHRVcount);
282  }
283 
284  if (status.testDAQStatus() == false) {
285  histograms.h_summary.h_daq_status_per_dom->Fill(string , floor, (1./distance(scanner.begin(), scanner.end()))*100);
286  }
287 
288  histograms.h_summary.h_fifo_per_dom->Fill(string , floor , nFIFOcount);
289 
290  if(analysis_level == 1){
291 
292  TH2D* h2 = (*histograms.h_summary.m_summary_rate_distribution)[MAKE_STRING("Detector/DU" + to_string(string) + "/F" + to_string(floor))];
293 
294  double rate = 0;
295  const double factor = 1.0e-3; // [kHz]
296 
297  for (int i = 0 ; i < NUMBER_OF_PMTS ; i++){
298 
299  rate += frame->getRate(i, factor);
300 
301  h2->Fill(i , frame->getRate(i, factor));
302 
303  PMT_rate_quantiles[string][floor][i].put(frame->getRate(i, factor));
304 
305  histograms.h_summary.h_pmt_rate_distribution->Fill(frame->getRate(i, factor), frame->getWeight(i, factor));
306 
307  }
308 
310 
311  DOM_rate_quantiles[string][floor].put(rate);
312 
313  } else {
314 
315  double rate = 0;
316  const double factor = 1.0e-3; // [kHz]
317 
318  for (int i = 0 ; i < NUMBER_OF_PMTS ; i++){
319 
320  rate += frame->getRate(i, factor);
321 
322  histograms.h_summary.h_pmt_rate_distribution->Fill(frame->getRate(i, factor), frame->getWeight(i, factor));
323  }
324 
326 
327  DOM_rate_quantiles[string][floor].put(rate);
328 
329  }
330 
331  }else{
332  FATAL("JModuleRouter trying to access non existing identifier: "<< frame->getModuleID());
333  }
334  }
335  }
336 
337  for (const auto& i1 : DOM_rate_quantiles) {
338  for (const auto& i2 : i1.second) {
339  if (i2.second.getCount() > 0) {
340  histograms.h_summary.h_rate_summary->Fill(i1.first, i2.first, i2.second.getMean());
341  }
342  }
343  }
344 
345  if (analysis_level == 1){
346  for (const auto& i1 : PMT_rate_quantiles){
347 
348  TH2D* h2 = (*histograms.h_summary.m_mean_summary_rate) [MAKE_STRING("Detector/DU" + to_string(i1.first))];
349  TH1D* h1 = (*histograms.h_summary.m_mean_summary_rate_distribution)[MAKE_STRING("Detector/DU" + to_string(i1.first))];
350 
351  for (const auto& i2 : i1.second) {
352  for (int i3 = 0 ; i3 != NUMBER_OF_PMTS ; ++i3) {
353  if (i2.second[i3].getCount() > 0){
354  h2->Fill(i3, i2.first, i2.second[i3].getMean());
355  h1->Fill(i2.second[i3].getMean());
356  }
357  }
358  }
359  }
360  }
361  }
362 
363  /*
364  * Function template to read and process time ordered timeslices.
365  *
366  * \param scanner A JTreeScanner from which time ordered timeslices are accessed.
367  */
368  template <class T>
370 
372 
373  std::map < int , std::map <int , JQuantile> > DOM_rate_quantiles;
374 
375  double inverseFrameTimeSec = 1. / (1.0e-9 * getFrameTime());
376 
377  while (scanner.hasNext()){
378 
379  T slice = *(scanner.next());
380 
381  for(auto & frame : slice) {
382  if (router.hasModule(frame.getModuleID())) {
383 
384  const JModule& module = router.getModule (frame.getModuleID());
385  double rate = frame.numberOfHits * inverseFrameTimeSec;
386  int string = module.getString();
387  int floor = module.getFloor ();
388 
389  DOM_rate_quantiles[string][floor].put(rate);
390 
391  vector <int> pmt_hit_count (NUMBER_OF_PMTS , 0) ;
392 
393  if(analysis_level == 1){
394 
395  TH2D* h2 = (*histograms.h_timeslice.m_pmt_tot_distributions[ts_type])[MAKE_STRING("Detector/DU" + to_string(string) + "/F" + to_string(floor) + "/" + to_string(frame.getModuleID()))];
396 
397  for (JDAQSuperFrame::const_iterator hit = frame.begin() ; hit != frame.end() ; ++hit){
398 
399  h2 -> Fill(hit->getPMT() , hit->getToT()) ;
400 
401  pmt_hit_count[hit->getPMT()]++;
402 
403  }
404 
405  for (int pmt = 0 ; pmt != NUMBER_OF_PMTS ; ++pmt) {
406 
407  (*histograms.h_timeslice.m_pmt_rate_distributions[ts_type])[MAKE_STRING("Detector/DU" + to_string(string) + "/F" + to_string(floor))] -> Fill(pmt , 1e-3 * pmt_hit_count[pmt] * inverseFrameTimeSec);
408 
409  }
410 
411  }
412 
413  }else{
414  FATAL("JModuleRouter trying to access non existing identifier: "<< frame.getModuleID());
415  }
416  }
417  }
418 
419  for (std::map< int , std::map<int,JQuantile>>::const_iterator i = DOM_rate_quantiles.begin() ; i != DOM_rate_quantiles.end() ; i++){
420 
421  for (std::map<int , JQuantile>::const_iterator j = i->second.begin() ; j != i->second.end() ; j++){
422 
423  if (j->second.getCount() > 0) histograms.h_timeslice.h_dom_mean_rates[ts_type] -> Fill (i->first , j->first , j->second.getMean() ) ;
424  }
425  }
426  }
427 
428  /*
429  * Checks through the use of a JTreeScanner, if summary slices are present in the run file.
430  * In case they are, the corresponding tree is iterated.
431  *
432  */
434 
436  if (scanner.hasNext()) {
438  scanner.rewind();
439  iterateSummarysliceTree(scanner);
440  }
441  }
442 
443  /*
444  * Function template that checks through the use of a JTreeScanner, if objects of different timeslice classes are present in the run file.
445  * In case they are, the corresponding tree is iterated.
446  *
447  */
448  template <class T>
450 
452  if(scanner.hasNext()) {
453 
455  scanner.rewind();
456  iterateTimesliceTree(scanner);
457  }
458  }
459 
460  /*
461  * Function template that checks through the use of a JTreeScanner, if JDAQEvent objects are present in the run file.
462  * In case they are, the corresponding tree is read.
463  *
464  */
465  void readEvents() {
466 
468 
469  if(scanner.hasNext()) {
470 
472  scanner.rewind();
473  iterateEventTree(scanner);
474  }
475  }
476 
477  /*
478  * Function template that checks through the use of a JTreeScanner, if JEvt objects are present in the run file.
479  * In case they are, the corresponding tree is read.
480  */
481  void readRecoEvents() {
482 
484 
485  if(scanner.hasNext()) {
486 
487  scanner.rewind();
488  iterateRecoEventTree(scanner, *outputFile);
489  }
490  }
491 
492  /*
493  * Writes the histograms to a root file.
494  * \param f The root file.
495  */
496  void writeToFile(TFile *f,int analysis_level){
497 
498  f->mkdir("Detector");
499  f->cd("Detector");
500 
507 
509 
510  for (typename vector < JManager < string , TH2D >* >::const_iterator i = histograms.h_timeslice.m_pmt_tot_distributions.begin();
512 
513  if ((*i)){
514 
515  for (typename JManager < string , TH2D >::const_iterator j = (*i) -> begin() ; j != (*i) -> end() ; ++j){
516 
517  j->second->Scale(1., "width");
518  }
519  }
520  }
521 
524 
526 
529 
531 
534 
535  f->mkdir ( MAKE_STRING ("JDAQEvent").c_str());
536  f->cd ("JDAQEvent");
537 
554 
555  }
556 };
557 
558 #endif
TH1D * h_pmt_distribution_snapshot_hits
TH2D * h_Snapshot_hits_per_module
DAQ key hit.
Definition: JDAQKeyHit.hh:19
void iterateTimesliceTree(JTreeScanner< T > &scanner)
ROOT TTree parameter settings of various packages.
static const unsigned int NUMBER_OF_TRIGGER_BITS
Number of trigger bits.
int countFIFOStatus() const
Count FIFO status.
TH1D * h_pmt_distribution_triggered_hits
vector< JManager< string, TH2D > * > m_pmt_tot_distributions
JManager< string, TH1D > * m_mean_summary_rate_distribution
int getFloor() const
Get floor number.
Definition: JLocation.hh:145
const JModule & getModule(const JObjectID &id) const
Get module parameters.
Data structure for a composite optical module.
Definition: JModule.hh:67
Auxiliary class for trigger mask.
then usage $script[< detector identifier >< run range >]< QA/QCfile > nExample script to produce data quality plots nWhen a detector identifier and run range are data are downloaded from the database nand subsequently stored in the given QA QC file
Definition: JDataQuality.sh:19
virtual const multi_pointer_type & next() override
Get next element.
range_type & include(argument_type x)
Include given value to range.
Definition: JRange.hh:397
JTrack3E getTrack(const Trk &track)
Get track.
std::vector< T >::difference_type distance(typename std::vector< T >::const_iterator first, typename PhysicsEvent::const_iterator< T > second)
Specialisation of STL distance.
vector< TH2D * > h_dom_mean_rates
Template specialisation of L0 builder for JHitL0 data type.
Definition: JBuildL0.hh:102
TFile * outputFile
Definition: JRunAnalyzer.hh:36
JRA_Histograms histograms
Definition: JRunAnalyzer.hh:37
TH2D * h_Triggered_hits_per_module
Detector data structure.
Definition: JDetector.hh:89
Template const_iterator.
Definition: JDAQEvent.hh:66
Router for direct addressing of module data in detector data structure.
then set_variable singlesRate set_variable doublesRate set_variable numberOfSlices echo Generating random background echo Singles rate
void readRecoEvents()
o $QUALITY_ROOT d $DEBUG!CHECK_EXIT_CODE JPlot1D f
Definition: JDataQuality.sh:76
void writeToFile(TFile *f, int analysis_level)
General purpose class for parallel reading of objects from a single file or multiple files...
JLimit_t numberOfTimeslices
Definition: JRunAnalyzer.hh:38
Dynamic ROOT object management.
void initialize_summary_histograms()
double getTime(const Hit &hit)
Get true time of hit.
double getPhi() const
Get phi angle.
Definition: JVersor3D.hh:144
JSingleFileScanner inputFile
Definition: JRunAnalyzer.hh:34
void iterateSummarysliceTree(JTreeScanner< JDAQSummaryslice > &scanner)
JManager< string, TH2D > * m_mean_summary_rate
Class dedicated to the analysis of KM3NeT runs.
Definition: JRunAnalyzer.hh:32
Template definition for direct access of elements in ROOT TChain.
JVersor3D getDirection(const JVector3D &pos) const
Get photon direction of Cherenkov light on PMT.
Definition: JTrack3D.hh:147
TH1D * h_dom_rate_distribution
#define MAKE_STRING(A)
Make string.
Definition: JPrint.hh:127
void readTimesliceData()
Basic data structure for L0 hit.
T getLength() const
Get length (difference between upper and lower limit).
Definition: JRange.hh:289
Type list.
Definition: JTypeList.hh:22
void iterateRecoEventTree(JParallelFileScanner< JTypeList< JDAQEvent, JFIT::JEvt > > &scanner, TFile &f)
Auxiliary class for defining the range of iterations of objects.
Definition: JLimit.hh:41
Type definition of range.
Definition: JHead.hh:41
Auxiliary class to manage set of compatible ROOT objects (e.g. histograms) using unique keys...
Definition: JManager.hh:43
TimesliceHistograms h_timeslice
~JRunAnalyzer()
Destructor.
Definition: JRunAnalyzer.hh:75
JModuleRouter router
Definition: JRunAnalyzer.hh:35
Detector file.
Definition: JHead.hh:226
Acoustic event fit.
void setLogarithmicX(TList *list)
Make x-axis of objects in list logarithmic (e.g. after using log10()).
TH1D * h_tot_distribution_snapshot_hits
Hit data structure.
Definition: JDAQHit.hh:34
Mathematical constants.
TriggerHistograms h_trigger
double getFrameTime()
Get frame time duration.
Definition: JDAQClock.hh:162
do set_variable OUTPUT_DIRECTORY $WORKDIR T
double getTheta() const
Get theta angle.
Definition: JVersor3D.hh:128
Parallel scanning of objects from a single file or multiple files according a format that follows fro...
void Write_manager_table_in_key_dir(TFile &f, vector< JManager< T, V > * > table)
then awk string
Support methods.
Data time slice.
static const int JGANDALF_BETA0_RAD
KM3NeT Data Definitions v3.3.0-2-g5cc95cf https://git.km3net.de/common/km3net-dataformat.
TH1D * h_pmt_rate_distribution
JRunAnalyzer(const JSingleFileScanner<> &file, const JDetector &detector, TFile *out, JLimit_t nTimeslices, JLimit_t nSummaryslices, JLimit_t nEvents, int analysislevel)
Constructor.
Definition: JRunAnalyzer.hh:56
void Write_histogram_table_to_file(TFile &f, string dirname, vector< vector< T * > > table)
void readSummaryData()
const JPMT & getPMT(const int index) const
Get PMT.
Definition: JModule.hh:172
double getT(const JVector3D &pos) const
Get arrival time of Cherenkov light at given position.
Definition: JTrack3D.hh:126
void iterateEventTree(JTreeScanner< JDAQEvent > &scanner)
Definition: JRunAnalyzer.hh:83
vector< JManager< string, TH1D > * > m_mean_ToT_distribution
#define FATAL(A)
Definition: JMessage.hh:67
Scanning of objects from multiple files according a format that follows from the extension of each fi...
SummaryHistograms h_summary
Direct access to module in detector data structure.
z range($ZMAX-$ZMIN)< $MINIMAL_DZ." fi fi typeset -Z 4 STRING typeset -Z 2 FLOOR JPlot1D -f $
int getString() const
Get string number.
Definition: JLocation.hh:134
Auxiliary data structure for average.
Definition: JKatoomba_t.hh:76
Auxiliary class to define a range between two values.
int countHighRateVeto() const
Count high-rate veto status.
JManager< string, TH2D > * m_summary_rate_distribution
JLimit_t numberOfEvents
Definition: JRunAnalyzer.hh:40
void Write_manager_in_key_dir(TFile &f, JManager< T, V > *manager)
General purpose string class.
Definition: JHead.hh:152
std::string to_string(const T &value)
Convert value to string.
vector< JManager< string, TH2D > * > m_mean_ToT
vector< JManager< string, TH2D > * > m_pmt_rate_distributions
bool hasModule(const JObjectID &id) const
Has module.
TH1D * h_tot_distribution_triggered_hits
TH1D * h_Triggered_hits_3dmuon
JLimit_t numberOfSummaryslices
Definition: JRunAnalyzer.hh:39
Object reading from a list of files.
int j
Definition: JPolint.hh:792
TH2D * h_Triggered_hits_3dmuon_per_module
Indexing of data type in type list.
Definition: JTypeList.hh:310
static const int NUMBER_OF_PMTS
Total number of PMTs in module.
Definition: JDAQ.hh:26
bool testDAQStatus() const
Test DAQ status of packets.
void readEvents()
void initialize_trigger_histograms()
JManager< string, TH2D > * m_Snapshot_hits_per_pmt
void initialize_timeslice_histograms()
TH1D * h_Triggered_over_Snapshot_hits