Jpp  18.0.1-rc.2
the software that should make you happy
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
RunAnalyzer.hh
Go to the documentation of this file.
1 #ifndef __RUNANALYZER__
2 #define __RUNANALYZER__
3 
4 #include "JSupport/JSupport.hh"
9 #include "JROOT/JManager.hh"
10 #include "JTools/JQuantile.hh"
11 #include "TFile.h"
12 #include "TH1D.h"
13 #include "JRunHistograms.hh"
14 
15 using namespace std ;
16 using namespace KM3NETDAQ ;
17 using namespace JPP ;
18 
19 /**
20  * \author rgruiz, adomi
21  */
22 
23 /**
24  * Class dedicated to the analysis of %KM3NeT runs.
25  *
26  */
27 class RunAnalyzer {
28 
29  string inputFile;
38 
39 public :
40 
41  /**
42  * Default constructor.
43  */
45 
46  /**
47  * Constructor.
48  *
49  * \param file Path to the file to be analyzed
50  * \param detectorFile path to a detector 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 frames to be read
54  * \param pmtanalysis PMT analysis option
55  */
56  RunAnalyzer (string file , string detectorFile , JLimit_t nTimeslices , JLimit_t nSummaryslices , JLimit_t nEvents, bool pmtanalysis):
57  inputFile (file),
58  router (NULL),
59  numberOfTimeslices (nTimeslices),
60  numberOfSummaryslices (nSummaryslices),
61  numberOfEvents (nEvents),
62  pmt_analysis (pmtanalysis)
63  {
64  try {
65  load (detectorFile, detector);
66  }
67  catch(const JException & error) {
68  cerr << "FATAL ERROR. Could not open detector file '" << detectorFile << "'." << endl;
69  exit(1);
70  }
71 
72  histograms = JRA_Histograms(detector);
73  router = new JModuleRouter(detector);
75  }
76 
77  /**
78  * Destructor.
79  */
81 
82  /*
83  * Function template used to read time ordered events from a tree.
84  *
85  * \param scanner A JTreeScanner
86  * \param
87  */
89 
90  while (scanner.hasNext()) {
91 
92  JDAQEvent event = *(scanner.next());
93 
94  histograms.h_trigger.h_Snapshot_hits -> Fill((Double_t) event.size<JDAQSnapshotHit > ());
95  histograms.h_trigger.h_Triggered_hits -> Fill((Double_t) event.size<JDAQTriggeredHit> ());
96  histograms.h_trigger.h_Number_of_overlays -> Fill(event.getOverlays());
97 
98  for (unsigned int i = 0; i != NUMBER_OF_TRIGGER_BITS; ++i) {
99 
100  if (event.hasTriggerBit(i)) {
101 
102  histograms.h_trigger.h_Trigger_bit_event -> Fill((Double_t) i);
103  }
104  }
105 
106  for (JDAQEvent::const_iterator<JDAQTriggeredHit> hit = event.begin<JDAQTriggeredHit>(); hit != event.end<JDAQTriggeredHit>(); ++hit) {
107 
108  if (router->hasModule(hit->getModuleID())) {
109 
110  const JModule& module = router->getModule (hit->getModuleID());
111 
112  histograms.h_trigger.h_pmt_distribution_triggered_hits->Fill(hit->getPMT());
113  histograms.h_trigger.h_tot_distribution_triggered_hits->Fill(hit->getToT());
114 
115  int String = module.getString();
116  int Floor = module.getFloor();
117 
118  histograms.h_trigger.h_Triggered_hits_per_module -> Fill(String , Floor);
119 
120  for (unsigned int i = 0; i != NUMBER_OF_TRIGGER_BITS; ++i) {
121 
122  if (hit -> hasTriggerBit(i)) {
123 
124  histograms.h_trigger.h_Trigger_bit_hit->Fill((Double_t) i);
125  }
126  }
127  }else{
128  FATAL("JModuleRouter trying to access non existing identifier: "<< hit->getModuleID());
129  }
130  }
131 
132  for (JDAQEvent::const_iterator<JDAQSnapshotHit> hit = event.begin<JDAQSnapshotHit>() ; hit != event.end<JDAQSnapshotHit>() ; ++hit){
133 
134  if (router->hasModule(hit->getModuleID())) {
135 
136  const JModule& module = router->getModule (hit->getModuleID());
137 
138  int String = module.getString();
139  int Floor = module.getFloor();
140  int pmt = hit-> getPMT();
141 
142  if(pmt_analysis == true){
143 
144  (*histograms.h_trigger.m_Snapshot_hits_per_pmt)[MAKE_STRING("Detector/DU" + to_string(String))] -> Fill(pmt, Floor);
145 
146  }
147 
148  histograms.h_trigger.h_pmt_distribution_snapshot_hits -> Fill(hit->getPMT());
149  histograms.h_trigger.h_tot_distribution_snapshot_hits -> Fill(hit->getToT());
150  histograms.h_trigger.h_Snapshot_hits_per_module -> Fill(String, Floor);
151 
152 
153  }else{
154  FATAL("JModuleRouter trying to access non existing identifier: "<< hit->getModuleID());
155  }
156  }
157  }
158  }
159 
160 
161  /*
162  * Function template to read and process time ordered summary slices.
163  *
164  * \param scanner A JTreeScanner from which time ordered summary slices are accessed.
165  */
167 
168  std::map <int,vector<vector <JQuantile>>> PMT_rate_quantiles;
169  std::map <int,vector<JQuantile>> DOM_rate_quantiles;
170 
171  while (scanner.hasNext()){
172 
173  JDAQSummaryslice slice = *(scanner.next());
174 
175  for (JDAQSummaryslice::const_iterator frame = slice.begin() ; frame != slice.end() ; ++ frame) {
176 
177  if (router->hasModule(frame->getModuleID())) {
178 
179  const JModule& module = router->getModule (frame->getModuleID());
180 
181  int string = module.getString();
182  int floor = module.getFloor ();
183 
184  PMT_rate_quantiles[string].resize(modulesPerString , vector <JQuantile>(NUMBER_OF_PMTS));
185 
186  JDAQFrameStatus status = frame -> getDAQFrameStatus();
187  int nFIFOcount = status.countFIFOStatus();
188  int nHRVcount = status.countHighRateVeto();
189 
190  if (nHRVcount > 0) {
191  histograms.h_summary.h_hrv_per_dom->Fill(string , floor);
192  }
193 
194  if (status.testDAQStatus() == false) {
195  histograms.h_summary.h_daq_status_per_dom->Fill(string , floor);
196  }
197 
198  histograms.h_summary.h_fifo_per_dom->Fill(string , floor , nFIFOcount);
199 
200  if(pmt_analysis == true){
201 
202  TH2D* h2 = (*histograms.h_summary.m_summary_rate_distribution)[MAKE_STRING("Detector/DU" + to_string(string) + "/F" + to_string(floor))];
203 
204  double rate = 0;
205 
206  for (int i = 0 ; i < NUMBER_OF_PMTS ; i++){
207 
208  rate += frame->getRate(i);
209 
210  h2->Fill(i , frame->getRate(i) * 1e-3);
211 
212  PMT_rate_quantiles[string][floor-1][i].put(frame->getRate(i)*1e-3);
213 
214  }
215 
216  DOM_rate_quantiles[string].resize(modulesPerString);
217 
218  DOM_rate_quantiles[string][floor-1].put(rate * 1e-3);
219 
220  } else {
221 
222  double rate = 0;
223 
224  for (int i = 0 ; i < NUMBER_OF_PMTS ; i++){
225  rate += frame->getRate(i);
226  histograms.h_summary.h_pmt_rate_distribution->Fill(frame->getRate(i)*1e-3);
227  }
228 
229  histograms.h_summary.h_dom_rate_distribution->Fill(rate*1e-3);
230 
231  DOM_rate_quantiles[string].resize(modulesPerString);
232 
233  DOM_rate_quantiles[string][floor-1].put(rate * 1e-3);
234 
235  }
236 
237  }else{
238  FATAL("JModuleRouter trying to access non existing identifier: "<< frame->getModuleID());
239  }
240  }
241  }
242 
243 
244  for (std::map<int , vector<JQuantile> >::const_iterator i = DOM_rate_quantiles.begin() ; i!= DOM_rate_quantiles.end() ; ++i) {
245 
246  for (int j=1 ; j < modulesPerString + 1 ; j++) {
247 
248  if (i->second[j-1].getCount() > 0) histograms.h_summary.h_rate_summary->Fill(i->first, j, i->second[j-1].getMean());
249  }
250  }
251 
252 
253  if(pmt_analysis == true){
254 
255  for (std::map<int , vector < vector<JQuantile> > >::const_iterator i = PMT_rate_quantiles.begin() ; i!= PMT_rate_quantiles.end() ; ++i) {
256 
257  for (int j=0 ; j < modulesPerString ; j++){
258 
259  for (int k=0 ; k < NUMBER_OF_PMTS ; k++){
260 
261  if (i -> second[j][k].getCount() > 0){
262 
263  (*histograms.h_summary.m_mean_summary_rate) [MAKE_STRING("Detector/DU" + to_string(i->first))]->Fill(k, j+1, i->second[j][k].getMean());
264  (*histograms.h_summary.m_mean_summary_rate_distribution)[MAKE_STRING("Detector/DU" + to_string(i->first))]->Fill( i->second[j][k].getMean());
265  }
266  }
267  }
268  }
269  }
270  }
271 
272  /*
273  * Function template to read and process time ordered timeslices.
274  *
275  * \param scanner A JTreeScanner from which time ordered timeslices are accessed.
276  * \param frame_index_range the range of the indices of the frames to be accessed.
277  */
278  template <class T>
280 
282 
283  std::map < int , std::map <int , JQuantile> > DOM_rate_quantiles;
284 
285  double inverseFrameTimeSec = 1. / (1.0e-9 * getFrameTime());
286 
287  while (scanner.hasNext()){
288 
289  T slice = *(scanner.next());
290 
291  for(auto & frame : slice) {
292  if (router->hasModule(frame.getModuleID())) {
293 
294  const JModule& module = router->getModule (frame.getModuleID());
295  double rate = frame.numberOfHits * inverseFrameTimeSec;
296  int string = module.getString();
297  int floor = module.getFloor ();
298 
299  DOM_rate_quantiles[string][floor].put(rate);
300 
301  vector <int> pmt_hit_count (NUMBER_OF_PMTS , 0) ;
302 
303  if(pmt_analysis == true){
304 
305  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()))];
306 
307  for (JDAQSuperFrame::const_iterator hit = frame.begin() ; hit != frame.end() ; ++hit){
308 
309  h2 -> Fill(hit->getPMT() , hit->getToT()) ;
310 
311  pmt_hit_count[hit->getPMT()]++;
312 
313  }
314 
315  for (int pmt = 0 ; pmt != NUMBER_OF_PMTS ; ++pmt) {
316 
317  (*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);
318 
319  }
320 
321  }
322 
323  }else{
324  FATAL("JModuleRouter trying to access non existing identifier: "<< frame.getModuleID());
325  }
326  }
327  }
328 
329 
330  for (std::map< int , std::map<int,JQuantile>>::const_iterator i = DOM_rate_quantiles.begin() ; i != DOM_rate_quantiles.end() ; i++){
331 
332  for (std::map<int , JQuantile>::const_iterator j = i->second.begin() ; j != i->second.end() ; j++){
333 
334  if (j->second.getCount() > 0) histograms.h_timeslice.h_dom_mean_rates[ts_type] -> Fill (i->first , j->first , j->second.getMean() ) ;
335  }
336  }
337  }
338 
339 
340  /*
341  * Checks through the use of a JTreeScanner, if summary slices are present in the run file.
342  * In case they are, the corresponding tree is iterated.
343  *
344  */
346 
347  JTreeScanner <JDAQSummaryslice> scanner(inputFile, numberOfSummaryslices);
348  if (scanner.hasNext()) {
349  histograms.initialize_summary_histograms();
350  scanner.rewind();
351  iterateSummarysliceTree(scanner);
352  }
353  }
354 
355 
356  /*
357  * Function template that checks through the use of a JTreeScanner, if objects of different timeslice classes are present in the run file.
358  * In case they are, the corresponding tree is iterated.
359  *
360  */
361  template <class T>
363 
364  JTreeScanner <T> scanner(inputFile, numberOfTimeslices);
365  if(scanner.hasNext()) {
366 
367  histograms.initialize_timeslice_histograms <T>();
368  scanner.rewind();
369  iterateTimesliceTree(scanner);
370  }
371  }
372 
373 
374  /*
375  * Function template that checks through the use of a JTreeScanner, if JDAQEvent objects are present in the run file.
376  * In case they are, the corresponding tree is read.
377  *
378  */
379  void readEvents() {
380 
381  JTreeScanner <JDAQEvent> scanner(inputFile, numberOfEvents);
382 
383  if(scanner.hasNext()) {
384 
385  histograms.initialize_trigger_histograms();
386  scanner.rewind();
387  iterateEventTree(scanner);
388  }
389  }
390 
391  /*
392  * Returns the histograms produced from the data in the file.
393  */
395 
396  return histograms ;
397  }
398 };
399 
400 #endif
JLimit_t numberOfSummaryslices
Definition: RunAnalyzer.hh:35
JRA_Histograms histograms
Definition: RunAnalyzer.hh:32
General exception.
Definition: JException.hh:23
RunAnalyzer(string file, string detectorFile, JLimit_t nTimeslices, JLimit_t nSummaryslices, JLimit_t nEvents, bool pmtanalysis)
Constructor.
Definition: RunAnalyzer.hh:56
void readEvents()
Definition: RunAnalyzer.hh:379
DAQ key hit.
Definition: JDAQKeyHit.hh:19
then fatal No hydrophone data file $HYDROPHONE_TXT fi sort gr k
void iterateEventTree(JTreeScanner< JDAQEvent > &scanner)
Definition: RunAnalyzer.hh:88
ROOT TTree parameter settings of various packages.
static const unsigned int NUMBER_OF_TRIGGER_BITS
Number of trigger bits.
void iterateSummarysliceTree(JTreeScanner< JDAQSummaryslice > &scanner)
Definition: RunAnalyzer.hh:166
int countFIFOStatus() const
Count FIFO status.
JLimit_t numberOfTimeslices
Definition: RunAnalyzer.hh:34
int getFloor() const
Get floor number.
Definition: JLocation.hh:145
Data structure for a composite optical module.
Definition: JModule.hh:68
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
Detector data structure.
Definition: JDetector.hh:89
Template const_iterator.
Definition: JDAQEvent.hh:68
Router for direct addressing of module data in detector data structure.
JLimit_t numberOfEvents
Definition: RunAnalyzer.hh:36
then set_variable singlesRate set_variable doublesRate set_variable numberOfSlices echo Generating random background echo Singles rate
exit
Definition: JPizza.sh:36
Dynamic ROOT object management.
void iterateTimesliceTree(JTreeScanner< T > &scanner)
Definition: RunAnalyzer.hh:279
void readTimesliceData()
Definition: RunAnalyzer.hh:362
Template definition for direct access of elements in ROOT TChain.
#define MAKE_STRING(A)
Make string.
Definition: JPrint.hh:127
string inputFile
Definition: RunAnalyzer.hh:29
Auxiliary class for defining the range of iterations of objects.
Definition: JLimit.hh:41
Detector file.
Definition: JHead.hh:226
~RunAnalyzer()
Destructor.
Definition: RunAnalyzer.hh:80
Class dedicated to the analysis of KM3NeT runs.
Definition: RunAnalyzer.hh:27
Hit data structure.
Definition: JDAQHit.hh:34
double getFrameTime()
Get frame time duration.
Definition: JDAQClock.hh:162
JRA_Histograms getHistograms()
Definition: RunAnalyzer.hh:394
do set_variable OUTPUT_DIRECTORY $WORKDIR T
then awk string
int modulesPerString
Definition: RunAnalyzer.hh:33
Support methods.
bool pmt_analysis
Definition: RunAnalyzer.hh:37
static const JStringCounter getNumberOfStrings
Function object to count unique strings.
#define FATAL(A)
Definition: JMessage.hh:67
Scanning of objects from multiple files according a format that follows from the extension of each fi...
Direct access to module in detector data structure.
int getString() const
Get string number.
Definition: JLocation.hh:134
void load(const std::string &file_name, JDetector &detector)
Load detector from input file.
int countHighRateVeto() const
Count high-rate veto status.
void readSummaryData()
Definition: RunAnalyzer.hh:345
General purpose string class.
Definition: JHead.hh:152
std::string to_string(const T &value)
Convert value to string.
int getCount(const T &hit)
Get hit count.
JModuleRouter * router
Definition: RunAnalyzer.hh:31
int j
Definition: JPolint.hh:703
JDetector detector
Definition: RunAnalyzer.hh:30
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:44