1#ifndef __JRUNHISTOGRAMS__
2#define __JRUNHISTOGRAMS__
26#include "TProfile2D.h"
34#include "TDirectory.h"
37#include "TObjString.h"
41using namespace JLANG ;
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);
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)));
104 " Fraction of wrong DAQ Status [%]; String ; Floor ; Fraction of slices with wrong DAQ status of packets [%]", du_ids, modules_per_string);
110 h_pmt_rate_distribution =
new TH1D (
"h_pmt_rate_distribution",
"PMT rate distribution from summary slices ; rate [kHz] ; Counts",
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));
116 TH1D* h =
new TH1D(
"%/h_mean_summary_rate_distribution",
" ; rate [kHz] ; # PMTs", 40 , 0 , log10(1000));
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));
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));
163 m_mean_ToT .resize (number_of_timeslice_types , NULL);
180 MAKE_STRING (ts_name +
" ; DU number ; Floor number ; time slice averaged rate [Hz]").c_str(),
182 du_ids, modules_per_string);
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));
190 MAKE_STRING (ts_name +
" ; ToT [ns] ; # PMTS ").c_str(),
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);
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 ,
239 for (
typename JManager < string , TH2D >::const_iterator
j = (*it) -> begin() ;
j != (*it) -> end() ; ++
j){
242 TPRegexp r (
"(\\w+)/(\\DU)(\\d+)/(F)(\\d+)");
244 TObjArray* o = r.MatchS(s);
246 int String = ((TObjString *)o->At(3))->GetString().Atoi();
247 int Floor = ((TObjString *)o->At(5))->GetString().Atoi();
249 for (
int pmt = 1 ; pmt <= (
j -> second) -> GetXaxis() -> GetNbins() ; pmt++){
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 () );
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 );
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 );
349 h_Snapshot_hits =
new TH1D (
"h_Snapshot_hits",
" ; Number of snapshot hits; Events ", 50, 0, 4 );
352 h_Triggered_hits =
new TH1D (
"h_Triggered_hits",
" ; Number of triggered hits; Events ", 50 , 0, 4 );
355 h_Triggered_hits_3dmuon =
new TH1D (
"h_Triggered_hits_3dmuon",
" ; Number of triggered hits for JTRIGGER3DMUON; Events ", 50 , 0, 3 );
359 "Number of triggered hits for JTRIGGER3DMUON; String ; Floor ; Number of JTRIGGER3DMUON hits",
360 du_ids, modules_per_string);
364 h_event_duration =
new TH1D (
"h_event_duration",
" ; Event Duration [ns]; Events", 60 , 1, 6 );
367 h_Number_of_overlays =
new TH1D (
"h_Number_of_overlays",
" ; Number of overlays; Events ", 1000, -0.5, 1000 - 0.5 );
369 " ; String ; Floor ; Number of snapshot hits ",
370 du_ids, modules_per_string);
373 " ; String ; Floor ; Number of triggered hits",
374 du_ids, modules_per_string);
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 ) );
382 " ; TDC Channel ; Counts [a.u.]", NUMBER_OF_PMTS , -0.5 , NUMBER_OF_PMTS - 0.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);
418 return std::distance(
du_ids.begin(),
du_ids.find(
string));
441 const string prefix =
"KM3NETDAQ::JDAQ" ;
442 string ts_name = T::Class_Name();
443 string::size_type pos = ts_name.find(prefix);
445 if (pos != string::npos) ts_name.replace(ts_name.find(prefix) , prefix.length() ,
"");
468 if(f.GetDirectory(dirname.c_str()) == 0) f.mkdir (dirname.c_str());
470 f.cd (dirname.c_str());
472 for (
int i=0 ; i < (int)table.size(); i++){
474 for (
int j=0 ;
j< (int)table[i].size();
j++){
476 if (table[i][
j]) table [i][
j] -> Write();
491 if(f.GetDirectory(dirname.c_str()) == 0) f.mkdir (dirname.c_str());
493 f.cd (dirname.c_str());
495 for (
int i=0 ; i < (int)table.size(); i++){
497 if (table[i]) table[i] -> Write();
508 template <
class T ,
class V>
511 if(f.GetDirectory(dirname.c_str()) == 0) f.mkdir (dirname.c_str());
513 f.cd (dirname.c_str());
515 for (
typename JManager < T , V >::const_iterator i = table -> begin() ; i != table -> end() ; ++i){
517 i -> second -> Write();
526 template <
class T ,
class V>
529 for (
typename JManager < T , V >::const_iterator i = manager -> begin() ; i != manager -> end() ; ++i){
531 if (i -> second -> GetTitle()){
533 std::string buffer = i -> second -> GetTitle();
534 string::size_type ipos = buffer.find(wc);
536 if (ipos!=std::string::npos){
542 buffer.replace(ipos, 1, os.str());
544 i -> second -> SetTitle(buffer.c_str());
555 template <
class T ,
class V >
558 for (
typename JManager < T , V >::const_iterator i = manager -> begin() ; i != manager -> end() ; ++i){
560 std::string fullpath =
MAKE_STRING(i->second->GetName());
562 int pos = fullpath.rfind (
'/');
563 std::string name = fullpath.substr (pos + 1);
564 std::string path = fullpath.substr (0 , pos);
566 if (f.GetDirectory(path.c_str()) == 0) f.mkdir (path.c_str());
570 i -> second -> SetName(name.c_str());
571 i -> second -> Write();
580 template <
class T ,
class V >
583 for (
typename vector < JManager < T , V >* >::const_iterator i = table.begin() ; i != table.end() ; ++i){
587 for (
typename JManager < T , V >::const_iterator
j = (*i) -> begin() ;
j != (*i) -> end() ; ++
j){
589 std::string fullpath =
MAKE_STRING(
j-> second -> GetName());
591 int pos = fullpath.rfind (
'/');
592 std::string name = fullpath.substr (pos + 1);
593 std::string path = fullpath.substr (0 , pos);
595 if (f.GetDirectory(path.c_str()) == 0) f.mkdir (path.c_str());
599 j -> second -> SetName(name.c_str());
600 j -> second -> Write();
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.
#define MAKE_STRING(A)
Make string.
TH2D * h2_detector_map_factory(const char *name, const char *title, std::set< int > &du_ids, int modules_per_string)
ROOT TTree parameter settings of various packages.
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)
SummaryHistograms h_summary
void Replace_wildcard_in_name(JManager< T, V > *manager, char wc='%')
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.
const char * getTriggerName(JTriggerbit_t bit)
Get trigger name.
KM3NeT DAQ data structures and auxiliaries.
static const unsigned int NUMBER_OF_TRIGGER_BITS
Number of trigger bits.
General purpose string class.
Indexing of data type in type list.
JManager< string, TH1D > * m_mean_summary_rate_distribution
JManager< string, TH2D > * m_mean_summary_rate
TH1D * h_pmt_rate_distribution
void initialize(std::set< int > &du_ids, int modules_per_string)
TH2D * h_daq_status_per_dom
JManager< string, TH2D > * m_summary_rate_distribution
TH1D * h_dom_rate_distribution
vector< JManager< string, TH2D > * > m_pmt_tot_distributions
void Fill_mean_ToT_histograms()
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_Trigger_bit_event
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_Number_of_overlays
TH1D * h_Triggered_over_Snapshot_hits
TH1D * h_Triggered_hits_3dmuon