Jpp 20.0.0-rc.1
the software that should make you happy
Loading...
Searching...
No Matches
JRunHistograms.hh
Go to the documentation of this file.
1#ifndef __JRUNHISTOGRAMS__
2#define __JRUNHISTOGRAMS__
3
4/**
5 * \author rgruiz, adomi
6 */
9
12
15
19
21
22#include "JROOT/JManager.hh"
23
24#include "TH1D.h"
25#include "TH2D.h"
26#include "TProfile2D.h"
27#include "TAxis.h"
28#include "TCanvas.h"
29#include "TPaveText.h"
30#include "TStyle.h"
31#include "TString.h"
32#include "TMath.h"
33#include "TColor.h"
34#include "TDirectory.h"
35#include "TPRegexp.h"
36#include "TObjArray.h"
37#include "TObjString.h"
38
39using namespace std ;
40using namespace KM3NETDAQ ;
41using namespace JLANG ;
42using namespace JPP ;
43using namespace JSUPPORT ;
44
45
46/*
47 * factory function to generate TH2D with x axis being DUs, with x labels set to DU numbers
48 */
49
50
51TH2D* h2_detector_map_factory(const char *name, const char *title, std::set<int> & du_ids , int modules_per_string){
52 TH2D * h2 = new TH2D (name, title,
53 du_ids.size(), -0.5 , du_ids.size() - 0.5 ,
54 modules_per_string , 0.5 , modules_per_string + 0.5);
55
56 for (Int_t ix = 1; ix <= h2->GetXaxis()->GetNbins(); ++ix) {
57 h2->GetXaxis()->SetBinLabel(ix, MAKE_CSTRING(*std::next(du_ids.begin(), ix-1)));
58 }
59
60 return h2;
61};
62
63
64/*
65 * Structure to store histograms obtained from the JDAQSummary TTree.
66 */
67
76
77 /* One histogram for each DU */
78 JManager < string , TH2D >* m_mean_summary_rate;
79 JManager < string , TH1D >* m_mean_summary_rate_distribution;
80
81 /* One histogram for each module */
82 JManager < string , TH2D >* m_summary_rate_distribution;
83
95
96 /*
97 * Initializes the histograms for summary slices
98 */
99 void initialize(std::set<int> & du_ids , int modules_per_string){
100
101 h_fifo_per_dom = h2_detector_map_factory("h_fifo_per_dom", " FIFO ; String ; Floor ; Number of slices with FIFO almost full x number of PMTs ", du_ids, modules_per_string);
102
103 h_daq_status_per_dom = h2_detector_map_factory("h_daq_status_per_dom",
104 " Fraction of wrong DAQ Status [%]; String ; Floor ; Fraction of slices with wrong DAQ status of packets [%]", du_ids, modules_per_string);
105
106 h_hrv_per_dom = h2_detector_map_factory("h_hrv_per_dom", "HRV ; String ; Floor ; Number of slices x number of PMTs in HRV", du_ids, modules_per_string);
107
108 h_rate_summary = h2_detector_map_factory("h_rate_summary", "Summary slices ; String ; Floor ; Mean rate over all summary slices [kHz]", du_ids, modules_per_string);
109
110 h_pmt_rate_distribution = new TH1D ("h_pmt_rate_distribution", "PMT rate distribution from summary slices ; rate [kHz] ; Counts",
112
113 h_dom_rate_distribution = new TH1D ("h_dom_rate_distribution", "DOM rate distribution from summary slices ; rate [kHz] ; Counts", 50 , log10(50) , log10(2000));
115
116 TH1D* h = new TH1D("%/h_mean_summary_rate_distribution", " ; rate [kHz] ; # PMTs", 40 , 0 , log10(1000));
117 setLogarithmicX (h);
118 m_mean_summary_rate_distribution = new JManager < string , TH1D > (h);
119
120 m_mean_summary_rate = new JManager < string , TH2D > (new TH2D("%/h_mean_summary_rate", " ; TDC Channel ; Floor ; rate [kHz]",
121 NUMBER_OF_PMTS , -0.5 , NUMBER_OF_PMTS - 0.5,
122 modules_per_string , 0.5 , 0.5 + modules_per_string));
123
124 TH2D* h_summary_rate_distribution = new TH2D ("%/h_pmt_rate_distributions_Summaryslice", "Summaryslice ; TDC channel ; rate [kHz] ; counts",
125 NUMBER_OF_PMTS , -0.5 , NUMBER_OF_PMTS - 0.5, 100 , -1 , log10(10000));
126 setLogarithmicY (h_summary_rate_distribution);
127 m_summary_rate_distribution = new JManager < string , TH2D > (h_summary_rate_distribution);
128
129 }
130
131};
132
134
139
140 /* One histogram per timeslice type. I decided not to use a JManager here because the range of each histogram could be different for each timeslice type. */
141
142 vector < TH2D* > h_dom_mean_rates;
143
144 /* One JManager per timeslice type. Each JManager hosts a histogram for each DU */
145
146 vector < JManager < string , TH2D >* > m_mean_ToT;
147 vector < JManager < string , TH1D >* > m_mean_ToT_distribution;
148
149 /* One JManager per timeslice type. Each manager hosts a histogram for an optical module. The key is expected to follow the pattern SXXFXX */
150
151 vector < JManager < string , TH2D >* > m_pmt_tot_distributions;
152 vector < JManager < string , TH2D >* > m_pmt_rate_distributions;
153
155 min_logdt (0),
156 max_logdt (9),
157 nbins_logdt (150),
158 nbins_time (200)
159 {
160 int number_of_timeslice_types = JLength<JDAQTimesliceTypes_t>::value ;
161
162 h_dom_mean_rates .resize (number_of_timeslice_types , NULL);
163 m_mean_ToT .resize (number_of_timeslice_types , NULL);
164 m_mean_ToT_distribution .resize (number_of_timeslice_types , NULL);
165 m_pmt_tot_distributions .resize (number_of_timeslice_types , NULL);
166 m_pmt_rate_distributions .resize (number_of_timeslice_types , NULL);
167 }
168
169 /*
170 * Initializes histograms for a given timeslice type.
171 *
172 * \param du_ids The list of ids for the DUs in the detector
173 * \param modules_per_string the number of modules in a string
174 * \param ts_type Index of the timeslice types on the JDAQTimesliceTypes_t typelist.
175 * \param ts_name The name of the timeslice type
176 */
177 void initialize(std::set<int> du_ids , int modules_per_string , int ts_type , std::string ts_name){
178
179 h_dom_mean_rates[ts_type] = h2_detector_map_factory(MAKE_STRING ("h_mean_dom_rates_" + ts_name).c_str(),
180 MAKE_STRING (ts_name + " ; DU number ; Floor number ; time slice averaged rate [Hz]").c_str(),
181
182 du_ids, modules_per_string);
183
184 m_mean_ToT[ts_type] = new JManager < string , TH2D > (new TH2D (MAKE_STRING ("%/h_mean_ToT_" + ts_name).c_str(),
185 MAKE_STRING (ts_name + " ; TDC channel ; Floor number ; mean ToT [ns] ").c_str(),
186 NUMBER_OF_PMTS , -0.5 , NUMBER_OF_PMTS - 0.5 ,
187 modules_per_string , 0.5 , 0.5 + modules_per_string));
188
189 m_mean_ToT_distribution[ts_type] = new JManager < string , TH1D > (new TH1D (MAKE_STRING ("%/h_mean_ToT_distribution" + ts_name).c_str(),
190 MAKE_STRING (ts_name + " ; ToT [ns] ; # PMTS ").c_str(),
191 256, -0.5, 255.5));
192
193 TH2D* h = new TH2D (MAKE_STRING ("%_" + ts_name + "_2SToT").c_str(),
194 MAKE_STRING (ts_name + " ; TDC channel ; ToT [ns] ; counts").c_str(),
195 NUMBER_OF_PMTS , -0.5 , NUMBER_OF_PMTS - 0.5 , 256, -0.5, 255.5);
196
197 h->Sumw2();
198
199 m_pmt_tot_distributions[ts_type] = new JManager < string , TH2D > (h);
200
201 int NBinsY=0;
202 double Ymax=0;
203
205 NBinsY=100;
206 Ymax=25;
208 NBinsY=100;
209 Ymax=2;
211 NBinsY=10;
212 Ymax=0.1;
213 }
214
215 TH2D* h_pmt_rate_distributions = new TH2D (MAKE_STRING ("%/h_pmt_rate_distributions_" + ts_name).c_str(),
216 MAKE_STRING (ts_name + " ; TDC channel ; rate [kHz] ; counts ").c_str(),
217 NUMBER_OF_PMTS , -0.5 , NUMBER_OF_PMTS - 0.5 ,
218 NBinsY, 0, Ymax);
219
220 m_pmt_rate_distributions[ts_type] = new JManager < string , TH2D > (h_pmt_rate_distributions);
221
222 }
223
224 /*
225 * Fills the mean ToT as a function of the PMT and floor number for a given DU
226 *
227 * \param table table with the ToT distributions for each PMT in a module, for every timeslice type
228 * \param string The string number
229 * \param floor The floor number
230 */
232
233 int i = 0 ;
234
235 for (typename vector < JManager < string , TH2D >* >::const_iterator it = m_pmt_tot_distributions.begin() ; it != m_pmt_tot_distributions.end() ; ++it , ++i){
236
237 if ((*it)){
238
239 for (typename JManager < string , TH2D >::const_iterator j = (*it) -> begin() ; j != (*it) -> end() ; ++j){
240
241 TString s (MAKE_STRING(j -> first).c_str());
242 TPRegexp r ("(\\w+)/(\\DU)(\\d+)/(F)(\\d+)");
243
244 TObjArray* o = r.MatchS(s);
245
246 int String = ((TObjString *)o->At(3))->GetString().Atoi();
247 int Floor = ((TObjString *)o->At(5))->GetString().Atoi();
248
249 for (int pmt = 1 ; pmt <= (j -> second) -> GetXaxis() -> GetNbins() ; pmt++){
250
251 (*m_mean_ToT[i])[MAKE_STRING("Detector/DU" + to_string(String))] -> Fill((j->second) -> GetXaxis() -> GetBinCenter(pmt) , Floor , (j -> second) -> ProjectionY ("" , pmt , pmt) -> GetMean () );
252 (*m_mean_ToT_distribution[i])[MAKE_STRING("Detector/DU" + to_string(String))] -> Fill((j -> second) -> ProjectionY ("" , pmt , pmt) -> GetMean () );
253 }
254 }
255 }
256 }
257 }
258
259};
260
261/*
262 * Structure to store histograms obtained from the JDAQEvent TTree
263 */
265
281
282 /* One histogram per DU*/
283 JManager < string , TH2D >* m_Snapshot_hits_per_pmt;
284
285 /*
286 * Constructor
287 */
306
307 /*
308 * Initializes the histograms.
309 * \param du_ids The list of du ids in the detector
310 * \param frame_index_range The range of frame indices
311 * \param modules_per_string The number of modules in a string.
312 */
313 void initialize(std::set<int> & du_ids , int modules_per_string){
314
315 h_Trigger_bit_event = new TH1D ("h_Trigger_bit_event",
316 "Number of events as a function of trigger bit in event ; Trigger Bit ; Counts",
317 NUMBER_OF_TRIGGER_BITS , -0.5, NUMBER_OF_TRIGGER_BITS - 0.5 );
318
319 for (int i = 0 ; i < (int) NUMBER_OF_TRIGGER_BITS ; i++){
320
321 if (getTriggerName(i)){
322
323 h_Trigger_bit_event -> GetXaxis() -> SetBinLabel(i+1 , getTriggerName(i) );
324 }else{
325
326 h_Trigger_bit_event -> GetXaxis() -> SetBinLabel(i+1 , "" );
327 }
328 }
329 h_Trigger_bit_event -> GetXaxis() -> LabelsOption("v");
330
331
332 h_Trigger_bit_hit = new TH1D ("h_Trigger_bit_hit",
333 "Number of hits per event as a function of trigger bit in hit ; Trigger Bit ; #Events",
334 NUMBER_OF_TRIGGER_BITS, -0.5, NUMBER_OF_TRIGGER_BITS - 0.5 );
335
336 for (int i = 0 ; i < (int) NUMBER_OF_TRIGGER_BITS ; i++){
337
338 if (getTriggerName(i)){
339
340 h_Trigger_bit_hit -> GetXaxis() -> SetBinLabel(i+1 , getTriggerName(i) );
341 }else{
342
343 h_Trigger_bit_hit -> GetXaxis() -> SetBinLabel(i+1 , "" );
344 }
345 }
346
347 h_Trigger_bit_hit -> GetXaxis() -> LabelsOption("v");
348
349 h_Snapshot_hits = new TH1D ("h_Snapshot_hits", " ; Number of snapshot hits; Events ", 50, 0, 4 );
351
352 h_Triggered_hits = new TH1D ("h_Triggered_hits", " ; Number of triggered hits; Events ", 50 , 0, 4 );
354
355 h_Triggered_hits_3dmuon = new TH1D ("h_Triggered_hits_3dmuon", " ; Number of triggered hits for JTRIGGER3DMUON; Events ", 50 , 0, 3 );
357
358 h_Triggered_hits_3dmuon_per_module = h2_detector_map_factory("h_Triggered_hits_3dmuon_per_module",
359 "Number of triggered hits for JTRIGGER3DMUON; String ; Floor ; Number of JTRIGGER3DMUON hits",
360 du_ids, modules_per_string);
361
362 h_Triggered_over_Snapshot_hits = new TH1D ("h_Triggered_over_Snapshot_hits", " ; Triggered/Snapshot hits; Events", 100 , 0, 0.5 );
363
364 h_event_duration = new TH1D ("h_event_duration", " ; Event Duration [ns]; Events", 60 , 1, 6 );
366
367 h_Number_of_overlays = new TH1D ("h_Number_of_overlays", " ; Number of overlays; Events ", 1000, -0.5, 1000 - 0.5 );
368 h_Snapshot_hits_per_module = h2_detector_map_factory("h_Snapshot_hits_per_module",
369 " ; String ; Floor ; Number of snapshot hits ",
370 du_ids, modules_per_string);
371
372 h_Triggered_hits_per_module = h2_detector_map_factory("h_Triggered_hits_per_module",
373 " ; String ; Floor ; Number of triggered hits",
374 du_ids, modules_per_string);
375
376 m_Snapshot_hits_per_pmt = new JManager < string , TH2D > ( new TH2D ("%/h_Snapshot_hits_per_pmt",
377 " ; TDC Channel ; Floor ; Number of snapshot hits",
378 NUMBER_OF_PMTS , -0.5 , NUMBER_OF_PMTS - 0.5 ,
379 modules_per_string , 0.5 , modules_per_string + 0.5 ) );
380
381 h_pmt_distribution_triggered_hits = new TH1D ("h_pmt_distribution_triggered_hits",
382 " ; TDC Channel ; Counts [a.u.]", NUMBER_OF_PMTS , -0.5 , NUMBER_OF_PMTS - 0.5);
383 h_tot_distribution_triggered_hits = new TH1D ("h_tot_distribution_triggered_hits", " ; ToT [ns] ; Counts [a.u.]", 256 , -0.5, 255.5);
384 h_pmt_distribution_snapshot_hits = new TH1D ("h_pmt_distribution_snapshot_hits", " ; TDC Channel ; Counts [a.u.]", NUMBER_OF_PMTS , -0.5 , NUMBER_OF_PMTS - 0.5);
385 h_tot_distribution_snapshot_hits = new TH1D ("h_tot_distribution_snapshot_hits", " ; ToT [ns] ; Counts [a.u.]", 256, -0.5, 255.5);
386 }
387};
388
389/*
390 * Class to manage the histograms produced by JRunAnalyzer.
391 */
393
394public:
395
399 std::set <int> du_ids;
401
403
412
413 /*
414 * Convert String number into histo string index
415 *
416 */
417 int string_to_index (int string){
418 return std::distance(du_ids.begin(), du_ids.find(string));
419 }
420
421
422 /*
423 * Initializes summary slice histograms.
424 *
425 * \param range The range of values for the summary slice indices in the run. The range of some histograms is based on the number of time slices, and their indices.
426 */
431
432 /*
433 * Initializes summary slice histograms.
434 *
435 * \param range The range of values for the summary slice indices in the run. The range of some histograms is based on the number of time slices, and their indices.
436 */
437 template <class T>
439
441 const string prefix = "KM3NETDAQ::JDAQ" ;
442 string ts_name = T::Class_Name();
443 string::size_type pos = ts_name.find(prefix);
444
445 if (pos != string::npos) ts_name.replace(ts_name.find(prefix) , prefix.length() , "");
446
447 h_timeslice.initialize(du_ids , modules_per_string , index , ts_name);
448 }
449
450 /*
451 * Initializes JDAQEvent histograms.
452 */
457
458 /*
459 * Checks if the histograms in a table have been initialized. If yes, it writes them on the indicated directory of a root file
460 *
461 * \param f The root file
462 * \param dirname The directory where the histograms will be written
463 * \param table A vector of histograms.
464 */
465 template <class T>
466 void Write_histogram_table_to_file(TFile & f , string dirname , vector < vector < T* > > table){
467
468 if(f.GetDirectory(dirname.c_str()) == 0) f.mkdir (dirname.c_str());
469
470 f.cd (dirname.c_str());
471
472 for (int i=0 ; i < (int)table.size(); i++){
473
474 for (int j=0 ; j< (int)table[i].size(); j++){
475
476 if (table[i][j]) table [i][j] -> Write();
477 }
478 }
479 }
480
481 /*
482 * Checks if the histograms in a table have been initialized. If yes, it writes them on the indicated directory of a root file
483 *
484 * \param f The root file
485 * \param dirname The directory where the histograms will be written
486 * \param table A vector of histograms.
487 */
488 template <class T>
489 void Write_histogram_table_to_file(TFile & f , string dirname , vector < T* > table){
490
491 if(f.GetDirectory(dirname.c_str()) == 0) f.mkdir (dirname.c_str());
492
493 f.cd (dirname.c_str());
494
495 for (int i=0 ; i < (int)table.size(); i++){
496
497 if (table[i]) table[i] -> Write();
498 }
499 }
500
501 /*
502 * Checks if the histograms in a table have been initialized. If yes, it writes them on the indicated directory of a root file
503 *
504 * \param f The root file
505 * \param dirname The directory where the histograms will be written
506 * \param table A vector of histograms.
507 */
508 template <class T , class V>
509 void Write_manager_to_file(TFile & f , string dirname , JManager < T , V >* table){
510
511 if(f.GetDirectory(dirname.c_str()) == 0) f.mkdir (dirname.c_str());
512
513 f.cd (dirname.c_str());
514
515 for (typename JManager < T , V >::const_iterator i = table -> begin() ; i != table -> end() ; ++i){
516
517 i -> second -> Write();
518 }
519 }
520
521 /*
522 * Replaces wildcard in manager objects titles by their keys.
523 * \param A manager.
524 * \param wc The wildcard
525 */
526 template <class T , class V>
527 void Replace_wildcard_in_name(JManager < T , V >* manager , char wc = '%'){
528
529 for (typename JManager < T , V >::const_iterator i = manager -> begin() ; i != manager -> end() ; ++i){
530
531 if (i -> second -> GetTitle()){
532
533 std::string buffer = i -> second -> GetTitle();
534 string::size_type ipos = buffer.find(wc);
535
536 if (ipos!=std::string::npos){
537
538 ostringstream os;
539
540 os << i -> first ;
541
542 buffer.replace(ipos, 1, os.str());
543
544 i -> second -> SetTitle(buffer.c_str());
545 }
546 }
547 }
548 }
549
550 /*
551 * Writes the contents of a JManager into a file. Each object in each will be written in a directory specified by its key.
552 * \param f The root file
553 * \param A manager.
554 */
555 template < class T , class V >
556 void Write_manager_in_key_dir(TFile & f ,JManager <T , V>* manager){
557
558 for (typename JManager < T , V >::const_iterator i = manager -> begin() ; i != manager -> end() ; ++i){
559
560 std::string fullpath = MAKE_STRING(i->second->GetName());
561
562 int pos = fullpath.rfind ('/');
563 std::string name = fullpath.substr (pos + 1);
564 std::string path = fullpath.substr (0 , pos);
565
566 if (f.GetDirectory(path.c_str()) == 0) f.mkdir (path.c_str());
567
568 f.cd(path.c_str());
569
570 i -> second -> SetName(name.c_str());
571 i -> second -> Write();
572 }
573 }
574
575 /*
576 * Writes the contents of a vector of JManager objects into a file. Each object in each will be written in a directory specified by its key.
577 * \param f The root file
578 * \param table the list of JManager objects.
579 */
580 template < class T , class V >
581 void Write_manager_table_in_key_dir(TFile & f , vector < JManager < T , V >* > table){
582
583 for (typename vector < JManager < T , V >* >::const_iterator i = table.begin() ; i != table.end() ; ++i){
584
585 if ((*i)){
586
587 for (typename JManager < T , V >::const_iterator j = (*i) -> begin() ; j != (*i) -> end() ; ++j){
588
589 std::string fullpath = MAKE_STRING(j-> second -> GetName());
590
591 int pos = fullpath.rfind ('/');
592 std::string name = fullpath.substr (pos + 1);
593 std::string path = fullpath.substr (0 , pos);
594
595 if (f.GetDirectory(path.c_str()) == 0) f.mkdir (path.c_str());
596
597 f.cd(path.c_str());
598
599 j -> second -> SetName(name.c_str());
600 j -> second -> Write();
601 }
602 }
603 }
604 }
605
606};
607
608#endif
KM3NeT DAQ constants, bit handling, etc.
Dynamic ROOT object management.
Direct access to module in detector data structure.
#define MAKE_CSTRING(A)
Make C-string.
Definition JPrint.hh:72
#define MAKE_STRING(A)
Make string.
Definition JPrint.hh:63
TH2D * h2_detector_map_factory(const char *name, const char *title, std::set< int > &du_ids, int modules_per_string)
Support methods.
ROOT TTree parameter settings of various packages.
Setting of trigger bits.
Detector data structure.
Definition JDetector.hh:96
TimesliceHistograms h_timeslice
void initialize_trigger_histograms()
void Write_manager_table_in_key_dir(TFile &f, vector< JManager< T, V > * > table)
void Write_histogram_table_to_file(TFile &f, string dirname, vector< T * > table)
int string_to_index(int string)
void Write_manager_to_file(TFile &f, string dirname, JManager< T, V > *table)
TriggerHistograms h_trigger
void Write_manager_in_key_dir(TFile &f, JManager< T, V > *manager)
void initialize_summary_histograms()
void initialize_timeslice_histograms()
JRA_Histograms(const JDetector &detector)
void Write_histogram_table_to_file(TFile &f, string dirname, vector< vector< T * > > table)
std::set< int > du_ids
SummaryHistograms h_summary
void Replace_wildcard_in_name(JManager< T, V > *manager, char wc='%')
T getUpperLimit() const
Get upper limit.
Definition JRange.hh:213
static int getN()
Get number of bins.
static const double * getData(const double factor=1.0)
Get abscissa values.
floor_range getRangeOfFloors(const JDetector &detector)
Get range of floors.
std::set< int > getStringIDs(const JDetector &detector)
Get list of strings identifiers.
void setLogarithmicX(TList *list)
Make x-axis of objects in list logarithmic (e.g. after using log10()).
void setLogarithmicY(TList *list)
Make y-axis of objects in list logarithmic (e.g. after using log10()).
Auxiliary classes and methods for language specific functionality.
std::string to_string(const T &value)
Convert value to string.
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
Support classes and methods for experiment specific I/O.
int j
Definition JPolint.hh:801
const char * getTriggerName(JTriggerbit_t bit)
Get trigger name.
KM3NeT DAQ data structures and auxiliaries.
Definition DataQueue.cc:39
static const unsigned int NUMBER_OF_TRIGGER_BITS
Number of trigger bits.
General purpose string class.
Definition JHead.hh:152
Detector file.
Definition JHead.hh:227
Indexing of data type in type list.
Definition JTypeList.hh:310
Length of type list.
Definition JTypeList.hh:176
JManager< string, TH1D > * m_mean_summary_rate_distribution
JManager< string, TH2D > * m_mean_summary_rate
void initialize(std::set< int > &du_ids, int modules_per_string)
JManager< string, TH2D > * m_summary_rate_distribution
vector< JManager< string, TH2D > * > m_pmt_tot_distributions
vector< JManager< string, TH2D > * > m_pmt_rate_distributions
vector< JManager< string, TH2D > * > m_mean_ToT
vector< TH2D * > h_dom_mean_rates
void initialize(std::set< int > du_ids, int modules_per_string, int ts_type, std::string ts_name)
vector< JManager< string, TH1D > * > m_mean_ToT_distribution
void initialize(std::set< int > &du_ids, int modules_per_string)
TH2D * h_Triggered_hits_3dmuon_per_module
TH2D * h_Triggered_hits_per_module
JManager< string, TH2D > * m_Snapshot_hits_per_pmt
TH1D * h_tot_distribution_snapshot_hits
TH1D * h_pmt_distribution_snapshot_hits
TH2D * h_Snapshot_hits_per_module
TH1D * h_pmt_distribution_triggered_hits
TH1D * h_tot_distribution_triggered_hits
TH1D * h_Triggered_over_Snapshot_hits