Jpp
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
RunAnalyzer.hh
Go to the documentation of this file.
1 #ifndef __NEWRUNANALYZER__
2 #define __NEWRUNANALYZER__
3 
4 #include "JSupport/JSupport.hh"
7 
10 
11 
12 #include "JGizmo/JManager.hh"
13 
14 
15 #include "JTools/JQuantile.hh"
16 #include "TFile.h"
17 #include "TH1D.h"
18 
19 #include "JRunHistograms.hh"
20 
21 using namespace std ;
22 
23 using namespace KM3NETDAQ ;
24 using namespace JDETECTOR ;
25 using namespace JLANG ;
26 using namespace JPP ;
27 
28 /**
29  * \author rgruiz
30  */
31 
32 
33 /**
34  * Class dedicated to the analysis of KM3NeT runs.
35  *
36  */
37 class RunAnalyzer {
38 
39  string inputfile ;
40 
42 
44 
46 
48 
50 
52 
53 public :
54 
55  /**
56  * Default constructor.
57  */
59 
60 
61 
62  /**
63  * Constructor.
64  *
65  * \param file Path to the file to be analyzed
66  * \param detector_file path to a detector file
67  */
68  RunAnalyzer (string file , string detector_file):
69 
70  inputfile (file),
71 
72  module_router (NULL)
73 
74  {
75 
76  try {
77 
78  load (detector_file, detector);
79 
80  }
81 
82  catch(const JException & error) {
83 
84  cerr << "FATAL ERROR. Could not open detector file '" << detector_file << "'." << endl ;
85 
86  exit(1) ;
87 
88  }
89 
90  histograms = JRA_Histograms (detector) ;
91 
92  module_router = new JModuleRouter(detector) ;
93 
94  DU_IDs = getStringIDs(detector) ;
95 
97 
98  }
99 
100 
101 
102  /**
103  * Destructor.
104  */
106 
107  }
108 
109 
110 
111  /*
112  * Function template used to read time ordered events from a tree.
113  *
114  * \param scanner A JTreeScanner
115  */
117 
118  while (scanner.hasNext()){
119 
120  JDAQEvent event = *(scanner.next()) ;
121 
122  int starting_frame_index = event.getFrameIndex () ;
123 
124  histograms.h_trigger.h_Snapshot_hits -> Fill((Double_t) event.size<JDAQSnapshotHit> ()) ;
125 
126  histograms.h_trigger.h_Triggered_hits -> Fill((Double_t) event.size<JDAQTriggeredHit> ()) ;
127 
128  for (unsigned int i = 0; i != NUMBER_OF_TRIGGER_BITS; ++i) {
129 
130  if (event.hasTriggerBit(i)) {
131 
132  histograms.h_trigger.h_Trigger_bit_event->Fill((Double_t) i);
133 
134  (*histograms.h_trigger.m_trigger_rates)[getTriggerName(i)] -> Fill (starting_frame_index - frame_index_range.first) ;
135 
136  }
137 
138  }
139 
140  histograms.h_trigger.h_n_triggered_hits_distribution -> Fill (event.size<JDAQTriggeredHit>()) ;
141 
143 
144  histograms.h_trigger.h_pmt_distribution_triggered_hits -> Fill (hit -> getPMT()) ;
145 
146  histograms.h_trigger.h_tot_distribution_triggered_hits -> Fill (hit -> getToT()) ;
147 
148  int String = module_router -> getModule (hit -> getModuleID()).getString() ;
149 
150  int Floor = module_router -> getModule (hit -> getModuleID()).getFloor() ;
151 
152  (*histograms.h_trigger.m_Trigger_map)[MAKE_STRING("Detector/DU" + to_string(String))] -> Fill (starting_frame_index , Floor) ;
153 
154  histograms.h_trigger.h_Triggered_hits_per_module -> Fill (String , Floor) ;
155 
156  for (unsigned int i = 0; i != NUMBER_OF_TRIGGER_BITS; ++i) {
157 
158  if (hit -> hasTriggerBit(i)) {
159 
160  histograms.h_trigger.h_Trigger_bit_hit->Fill((Double_t) i);
161 
162  }
163 
164  }
165 
166  }
167 
168  for (JDAQEvent::const_iterator<JDAQSnapshotHit> hit = event.begin<JDAQSnapshotHit>() ; hit != event.end<JDAQSnapshotHit>() ; ++hit){
169 
170  int String = module_router -> getModule (hit -> getModuleID()).getString() ;
171 
172  int Floor = module_router -> getModule (hit -> getModuleID()).getFloor() ;
173 
174  int pmt = hit -> getPMT() ;
175 
176  histograms.h_trigger.h_Snapshot_hits_per_module -> Fill (String , Floor) ;
177 
178  (*histograms.h_trigger.m_Snapshot_hits_per_pmt)[MAKE_STRING("Detector/DU" + to_string(String))] -> Fill (pmt , Floor) ;
179 
180  histograms.h_trigger.h_pmt_distribution_snapshot_hits -> Fill (hit -> getPMT()) ;
181 
182  histograms.h_trigger.h_tot_distribution_snapshot_hits -> Fill (hit -> getToT()) ;
183 
184  }
185 
186  }
187 
188  }
189 
190 
191  /*
192  * Function template to read and process time ordered summary slices.
193  *
194  * \param scanner A JTreeScanner from which time ordered summary slices are accessed.
195  * \param frame_index_range the range of the indices of the frames to be accessed.
196  */
198 
199  std::map < int , vector < vector < JQuantile > > > PMT_rate_quantiles ;
200 
201  std::map < int , vector <JQuantile> > DOM_rate_quantiles;
202 
203  while (scanner.hasNext()){
204 
205  JDAQSummaryslice slice = *(scanner.next()) ;
206 
207  double time_since_run_start = (slice.getFrameIndex() - frame_index_range.first) * getFrameTime() * 1.0e-9 ;
208 
209  JQuantile HRV_quantile ;
210 
211  JQuantile FIFO_quantile ;
212 
213  for (JDAQSummaryslice::const_iterator s_frame = slice.begin() ; s_frame != slice.end() ; ++ s_frame){
214 
215  PMT_rate_quantiles[module_router -> getModule (s_frame->getModuleID()).getString()].resize(modules_per_string , vector < JQuantile > (NUMBER_OF_PMTS)) ;
216 
217  JDAQFrameStatus status = s_frame -> getDAQFrameStatus() ;
218 
219  int nFIFOcount = status.countFIFOStatus() ;
220 
221  int nHRVcount = status.countHighRateVeto() ;
222 
223  HRV_quantile.put ( nHRVcount ) ;
224 
225  FIFO_quantile.put ( nFIFOcount ) ;
226 
227  int string = module_router -> getModule (s_frame->getModuleID()).getString() ;
228 
229  int floor = module_router -> getModule (s_frame->getModuleID()).getFloor() ;
230 
231  if ( nHRVcount > 0 ){
232 
233  histograms.h_summary.h_hrv_per_dom -> Fill (module_router -> getModule (s_frame->getModuleID()).getString() , module_router -> getModule (s_frame->getModuleID()).getFloor() ) ;
234 
235  }
236 
237  if (status.testDAQStatus() == false){
238 
239  histograms.h_summary.h_daq_status -> Fill(slice.getFrameIndex() - frame_index_range.first) ;
240 
241  histograms.h_summary.h_daq_status_per_dom -> Fill(module_router -> getModule (s_frame->getModuleID()).getString() , module_router -> getModule (s_frame->getModuleID()).getFloor()) ;
242 
243  }
244 
245  histograms.h_summary.h_fifo_per_dom -> Fill (module_router -> getModule (s_frame->getModuleID()).getString() , module_router -> getModule (s_frame->getModuleID()).getFloor() , nFIFOcount ) ;
246 
247  double rate = 0 ;
248 
249  for (int i = 0 ; i < NUMBER_OF_PMTS ; i++){
250 
251  (*histograms.h_summary.m_summary_rate_vs_time)[MAKE_STRING("Detector/DU" + to_string(string) + "/F" + to_string(floor))] -> Fill(time_since_run_start , i , s_frame -> getRate(i)) ;
252 
253  (*histograms.h_summary.m_summary_rate_distribution)[MAKE_STRING("Detector/DU" + to_string(string) + "/F" + to_string(floor))] -> Fill(s_frame -> getRate(i) * 1e-3 , i) ;
254 
255  (*histograms.h_summary.m_fifo_full)[MAKE_STRING("Detector/DU" + to_string(string) + "/F" + to_string(floor))] -> Fill(time_since_run_start , i , status.testFIFOStatus(i)) ;
256 
257  (*histograms.h_summary.m_hrv)[MAKE_STRING("Detector/DU" + to_string(string) + "/F" + to_string(floor))] -> Fill(time_since_run_start , i , status.testHighRateVeto(i)) ;
258 
259  rate += s_frame -> getRate(i) ;
260 
261  histograms.h_summary.h_pmt_rate_distribution -> Fill (s_frame -> getRate(i) * 1e-3) ;
262 
263  PMT_rate_quantiles[module_router -> getModule (s_frame->getModuleID()).getString()][module_router -> getModule (s_frame->getModuleID()).getFloor()-1][i].put(s_frame ->getRate(i) * 1e-3) ;
264 
265  }
266 
267  DOM_rate_quantiles[module_router -> getModule (s_frame->getModuleID()).getString()].resize(modules_per_string) ;
268 
269  DOM_rate_quantiles[module_router -> getModule (s_frame->getModuleID()).getString()][module_router -> getModule (s_frame->getModuleID()).getFloor()-1].put(rate * 1e-3) ;
270 
271  (*histograms.h_summary.m_module_rates_vs_time)[MAKE_STRING("Detector/DU" + to_string(string) + "/F" + to_string(floor))] -> Fill(time_since_run_start , rate ) ;
272 
273  }
274 
275  histograms.h_summary.h_fifo -> Fill ( slice.getFrameIndex() - frame_index_range.first , FIFO_quantile.getMean() ) ;
276 
277  histograms.h_summary.h_hrv -> Fill ( slice.getFrameIndex() - frame_index_range.first , HRV_quantile.getMean() ) ;
278 
279  }
280 
281  for (std::map<int , vector<JQuantile> >::const_iterator i = DOM_rate_quantiles.begin() ; i!= DOM_rate_quantiles.end() ; ++i){
282 
283  for (int j=1 ; j < modules_per_string + 1 ; j++){
284 
285  if (i -> second[j-1].getCount() > 0) histograms.h_summary.h_rate_summary ->Fill (i -> first , j , i -> second[j-1].getMean() ) ;
286 
287  }
288 
289  }
290 
291  for (std::map<int , vector < vector<JQuantile> > >::const_iterator i = PMT_rate_quantiles.begin() ; i!= PMT_rate_quantiles.end() ; ++i){
292 
293  for (int j = 0 ; j < modules_per_string ; j++){
294 
295  for (int k = 0 ; k < NUMBER_OF_PMTS ; k++){
296 
297  if (i -> second[j][k].getCount() > 0){
298 
299  (*histograms.h_summary.m_mean_summary_rate)[MAKE_STRING("Detector/DU" + to_string(i->first))] -> Fill(k , j+1 , i -> second[j][k].getMean()) ;
300 
301  }
302 
303  }
304 
305  }
306 
307  }
308 
309  }
310 
311 
312 
313  /*
314  * Checks through the use of a JTreeScanner, if summary slices are present in the run file.
315  * In case they are, the corresponding tree is iterated.
316  *
317  */
319 
320  JFrameIndexRange range ;
321 
323 
324  if( scanner.hasNext() ){
325 
326  range = JFrameIndexRange (scanner.begin()->getFrameIndex() , scanner.rbegin()->getFrameIndex()) ;
327 
328  histograms.initialize_summary_histograms (range) ;
329 
330  scanner.rewind () ;
331 
332  iterate_summaryslice_tree (scanner , range) ;
333 
334  }
335 
336  }
337 
338 
339 
340  /*
341  * Function template to read and process time ordered timeslices.
342  *
343  * \param scanner A JTreeScanner from which time ordered timeslices are accessed.
344  * \param frame_index_range the range of the indices of the frames to be accessed.
345  */
346  template <class T>
348 
350 
351  std::map < int , std::map <int , JQuantile> > DOM_rate_quantiles ;
352 
353  while (scanner.hasNext()){
354 
355  T slice = *(scanner.next()) ;
356 
357  double time_since_run_start = (slice.getFrameIndex() - frame_index_range.first) * getFrameTime() * 1.0e-9 ;
358 
359  JQuantile detector_rate_quantile ;
360 
361  std::map < int , JQuantile > active_DOM_DU_quantiles ;
362 
363  JQuantile active_DOM_quantile ;
364 
365  for(auto & s_frame : slice){
366 
367  JModule module = module_router -> getModule (s_frame.getModuleID()) ;
368 
369  double rate = s_frame.numberOfHits / (1.0e-9 * getFrameTime()) ;
370 
371  int string = module.getString() ;
372 
373  int floor = module.getFloor() ;
374 
375  DOM_rate_quantiles[string][floor].put(rate) ;
376 
377  detector_rate_quantile.put (rate) ;
378 
379  (*histograms.h_timeslice.m_module_rates_vs_time[ts_type])[MAKE_STRING("Detector/DU" + to_string(string) + "/F" + to_string(floor))] -> Fill(time_since_run_start , rate ) ;
380 
381  if ( rate > 0 ) {
382 
383  active_DOM_quantile.put ( 1 ) ;
384 
385  active_DOM_DU_quantiles [string].put(1) ;
386 
387  }
388 
389  std::map <int , vector < double > > hit_time_buffers ;
390 
391  std::map < int , JQuantile > tot_quantiles ;
392 
393  std::vector < JQuantile > pmt_rate_quantiles ;
394 
395  pmt_rate_quantiles.resize(NUMBER_OF_PMTS);
396 
397  for (JDAQSuperFrame::const_iterator hit = s_frame.begin() ; hit != s_frame.end() ; ++hit){
398 
399  (*histograms.h_timeslice.m_pmt_tot_distributions[ts_type])[MAKE_STRING("Detector/DU" + to_string(string) + "/F" + to_string(floor) + "/" + to_string(s_frame.getModuleID()))] -> Fill(hit->getPMT() , hit->getToT()) ;
400 
401  hit_time_buffers[hit->getPMT()].push_back (hit->getT()) ;
402 
403  tot_quantiles[hit->getPMT()].put (hit->getToT()) ;
404 
405  pmt_rate_quantiles[hit->getPMT()].put (1) ;
406 
407  if (hit->getToT() == 255) {
408 
409  (*histograms.h_timeslice.m_ToT_255[ts_type])[MAKE_STRING("Detector/DU" + to_string(string))] -> Fill(floor , hit->getPMT()) ;
410 
411  histograms.h_timeslice.h_ToT_255_vs_time[ts_type] -> Fill(hit -> getT()) ;
412 
413  histograms.h_timeslice.h_ToT_255_Floor_vs_time[ts_type] -> Fill(hit -> getT() , floor) ;
414 
415  histograms.h_timeslice.h_ToT_255_Floor_vs_time_2[ts_type] -> Fill(hit -> getT() , floor) ;
416 
417  }
418 
419  }
420 
421  for (std::vector<JQuantile>::const_iterator pmt = pmt_rate_quantiles.begin() ; pmt != pmt_rate_quantiles.end() ; pmt++){
422 
423  (*histograms.h_timeslice.m_pmt_rates_vs_time[ts_type])[MAKE_STRING("Detector/DU" + to_string(string) + "/F" + to_string(floor))] -> Fill(time_since_run_start , pmt - pmt_rate_quantiles.begin() , pmt->getCount()/getFrameTime()/1.0e-9 ) ;
424 
425  }
426 
427  for(std::map<int,JQuantile>::const_iterator pmt = tot_quantiles.begin() ; pmt != tot_quantiles.end() ; pmt++){
428 
429  if (pmt->second.getCount() > 0){
430 
431  (*histograms.h_timeslice.m_pmt_rate_distributions[ts_type])[MAKE_STRING("Detector/DU" + to_string(string) + "/F" + to_string(floor))] -> Fill(1e-3 * pmt->second.getCount()/getFrameTime()/1.0e-9 , pmt->first ) ;
432 
433  (*histograms.h_timeslice.m_pmt_tot_vs_time[ts_type])[MAKE_STRING("Detector/DU" + to_string(string) + "/F" + to_string(floor))] -> Fill (time_since_run_start , pmt->first , pmt->second.getMean() ) ;
434 
435  }
436 
437  for (int i=1 ; i < (int)hit_time_buffers[pmt->first].size() ; i++){
438 
439  (*histograms.h_timeslice.m_pmt_dt_consecutive_hits[ts_type])[MAKE_STRING("Detector/DU" + to_string(string) + "/F" + to_string(floor))] -> Fill ( log10(hit_time_buffers[pmt->first][i] - hit_time_buffers[pmt->first][i-1]) , pmt->first) ;
440 
441  }
442 
443  }
444 
445  }
446 
447  histograms.h_timeslice.h_rate[ts_type] -> Fill ( slice.getFrameIndex() - frame_index_range.first , detector_rate_quantile.getMean() ) ;
448 
449  histograms.h_timeslice.h_slice_start_time[ts_type] -> Fill ( slice.getFrameIndex() - frame_index_range.first , time_since_run_start ) ;
450 
451  histograms.h_timeslice.h_active_modules[ts_type] -> Fill ( slice.getFrameIndex() - frame_index_range.first , double(active_DOM_quantile.getCount()) / JDETECTOR::getNumberOfModules (detector) ) ;
452 
453  for ( std::map <int , JQuantile>::const_iterator i = active_DOM_DU_quantiles.begin() ; i != active_DOM_DU_quantiles.end() ; i++){
454 
455  histograms.h_timeslice.h_du_active_modules[ts_type] -> Fill(slice.getFrameIndex() - frame_index_range.first , i->first , i->second.getCount() ) ;
456 
457  }
458 
459  }
460 
461  for (std::map< int , std::map<int,JQuantile>>::const_iterator i = DOM_rate_quantiles.begin() ; i != DOM_rate_quantiles.end() ; i++){
462 
463  for (std::map<int , JQuantile>::const_iterator j = i->second.begin() ; j != i->second.end() ; j++){
464 
465  if (j->second.getCount() > 0) histograms.h_timeslice.h_dom_mean_rates[ts_type] -> Fill (i->first , j->first , j->second.getMean() ) ;
466 
467  }
468 
469  }
470 
471  }
472 
473  /*
474  * Function template that checks through the use of a JTreeScanner, if objects of different timeslice classes are present in the run file.
475  * In case they are, the corresponding tree is iterated.
476  *
477  */
478  template <class T>
480 
481  JFrameIndexRange range ;
482 
483  JTreeScanner < T , JDAQEvaluator > scanner(inputfile) ;
484 
485  if( scanner.hasNext() ){
486 
487  range = JFrameIndexRange (scanner.begin()->getFrameIndex() , scanner.rbegin()->getFrameIndex()) ;
488 
489  histograms.initialize_timeslice_histograms < T > (range) ;
490 
491  scanner.rewind () ;
492 
493  iterate_timeslice_tree (scanner , range) ;
494 
495  }
496 
497  }
498 
499 
500 
501  /*
502  * Function template that checks through the use of a JTreeScanner, if JDAQEvent objects are present in the run file.
503  * In case they are, the corresponding tree is read.
504  *
505  */
507 
508  JFrameIndexRange range ;
509 
510  JTreeScanner < JDAQEvent , JDAQEvaluator > scanner(inputfile) ;
511 
512  if( scanner.hasNext() ){
513 
514  range = JFrameIndexRange (scanner.begin()->getFrameIndex() , scanner.rbegin()->getFrameIndex()) ;
515 
516  histograms.initialize_trigger_histograms (range) ;
517 
518  scanner.rewind() ;
519 
520  iterate_daqevent_tree (scanner , range) ;
521 
522  }
523 
524  }
525 
526  /*
527  * Returns the histograms produced from the data in the file.
528  */
530 
531  return histograms ;
532 
533  }
534 
535 
536 };
537 
538 #endif
JRA_Histograms histograms
Definition: RunAnalyzer.hh:47
General exception.
Definition: JException.hh:40
JTOOLS::JRange< int > JFrameIndexRange
Type definition for frame index range.
DAQ key hit.
Definition: JDAQKeyHit.hh:24
void read_daqevents_from_file()
Definition: RunAnalyzer.hh:506
void read_summaryslices_from_file()
Definition: RunAnalyzer.hh:318
static const unsigned int NUMBER_OF_TRIGGER_BITS
Number of trigger bits.
int countFIFOStatus() const
Count FIFO status.
void iterate_timeslice_tree(JTreeScanner< T, JDAQEvaluator > &scanner, JFrameIndexRange frame_index_range)
Definition: RunAnalyzer.hh:347
Data structure for a composite optical module.
Definition: JModule.hh:47
int getNumberOfStrings(const JDetector &detector)
Get number of strings.
RunAnalyzer(string file, string detector_file)
Constructor.
Definition: RunAnalyzer.hh:68
double getCount(TH1D *hptr, int muon_threshold)
Detector data structure.
Definition: JDetector.hh:77
Template const_iterator.
Definition: JDAQEvent.hh:69
void put(const double x, const double w=1.0)
Put value.
Definition: JQuantile.hh:77
Router for direct addressing of module data in detector data structure.
std::set< int > getStringIDs(const JDetector &detector)
Get list of strings IDs.
Structure to store the ToT mean and standard deviation of the hits produced by a nanobeacon in a sour...
void read_timeslices_from_file()
Definition: RunAnalyzer.hh:479
Dynamic ROOT object management.
JModuleLocation default_module_location
Definition: RunAnalyzer.hh:45
Template definition for direct access of elements in ROOT TChain.
Definition: JTreeScanner.hh:89
double getMean(vector< double > &v)
get mean of vector content
Quantile calculator.
Definition: JQuantile.hh:34
#define MAKE_STRING(A)
Make string.
Definition: JPrint.hh:602
int getFrameIndex() const
Get frame index.
double getToT(const T &tot, const JCalibration &cal)
Get calibrated time-over-threshold of hit.
int getFloor() const
Get floor number.
int getString() const
Get string number.
JKey_t first
Definition: JPair.hh:128
Detector file.
Definition: JHead.hh:126
~RunAnalyzer()
Destructor.
Definition: RunAnalyzer.hh:105
JModuleRouter * module_router
Definition: RunAnalyzer.hh:43
Class dedicated to the analysis of KM3NeT runs.
Definition: RunAnalyzer.hh:37
Hit data structure.
Definition: JDAQHit.hh:36
double getFrameTime()
Get frame time duration.
Definition: JDAQClock.hh:162
JRA_Histograms getHistograms()
Definition: RunAnalyzer.hh:529
Support methods.
std::set< int > DU_IDs
Definition: RunAnalyzer.hh:49
void iterate_daqevent_tree(JTreeScanner< JDAQEvent, JDAQEvaluator > &scanner, JFrameIndexRange frame_index_range)
Definition: RunAnalyzer.hh:116
void iterate_summaryslice_tree(JTreeScanner< JDAQSummaryslice, JDAQEvaluator > &scanner, JFrameIndexRange frame_index_range)
Definition: RunAnalyzer.hh:197
void load(const JString &file_name, JDetector &detector)
Load detector from input file.
Scanning of objects from multiple files according a format that follows from the extension of each fi...
double getMean() const
Get mean value.
Definition: JQuantile.hh:131
Direct access to module in detector data structure.
bool testHighRateVeto() const
Test high-rate veto status.
int modules_per_string
Definition: RunAnalyzer.hh:51
int countHighRateVeto() const
Count high-rate veto status.
const char * getTriggerName(JTriggerbit_t bit)
Get trigger name.
long long int getCount() const
Get total count.
Definition: JQuantile.hh:97
General purpose string class.
Definition: JHead.hh:98
string inputfile
Definition: RunAnalyzer.hh:39
ROOT TTree parameter settings.
std::string to_string(const T &value)
Convert value to string.
JDetector detector
Definition: RunAnalyzer.hh:41
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.
int getNumberOfModules(const JDetector &detector)
Get number of modules.
RunAnalyzer()
Default constructor.
Definition: RunAnalyzer.hh:58
const JModule & getModule(const JDetector &detector, const JModuleLocation &location)
find module with a given string and floor number
Logical location of module.
bool testFIFOStatus() const
Test FIFO status.