Jpp  17.3.0-rc.1
the software that should make you happy
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Public Member Functions | Private Attributes | List of all members
RunAnalyzer Class Reference

Class dedicated to the analysis of KM3NeT runs. More...

#include <RunAnalyzer.hh>

Public Member Functions

 RunAnalyzer ()
 Default constructor. More...
 
 RunAnalyzer (string file, string detectorFile, JLimit_t nTimeslices, JLimit_t nSummaryslices, JLimit_t nEvents, bool pmtanalysis)
 Constructor. More...
 
 ~RunAnalyzer ()
 Destructor. More...
 
void iterateEventTree (JTreeScanner< JDAQEvent > &scanner)
 
void iterateSummarysliceTree (JTreeScanner< JDAQSummaryslice > &scanner)
 
template<class T >
void iterateTimesliceTree (JTreeScanner< T > &scanner)
 
void readSummaryData ()
 
template<class T >
void readTimesliceData ()
 
void readEvents ()
 
JRA_Histograms getHistograms ()
 

Private Attributes

string inputFile
 
JDetector detector
 
JModuleRouterrouter
 
JRA_Histograms histograms
 
int modulesPerString
 
JLimit_t numberOfTimeslices
 
JLimit_t numberOfSummaryslices
 
JLimit_t numberOfEvents
 
bool pmt_analysis
 

Detailed Description

Class dedicated to the analysis of KM3NeT runs.

Author
rgruiz, adomi

Definition at line 27 of file RunAnalyzer.hh.

Constructor & Destructor Documentation

RunAnalyzer::RunAnalyzer ( )
inline

Default constructor.

Definition at line 44 of file RunAnalyzer.hh.

44 {}
RunAnalyzer::RunAnalyzer ( string  file,
string  detectorFile,
JLimit_t  nTimeslices,
JLimit_t  nSummaryslices,
JLimit_t  nEvents,
bool  pmtanalysis 
)
inline

Constructor.

Parameters
filePath to the file to be analyzed
detectorFilepath to a detector file
nTimeslicesNumber of timeslices to be read
nSummaryslicesNumber of summary slices to be read
nEventsNumber of frames to be read
pmtanalysisPMT analysis option

Definition at line 56 of file RunAnalyzer.hh.

56  :
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 
75  }
JLimit_t numberOfSummaryslices
Definition: RunAnalyzer.hh:35
JRA_Histograms histograms
Definition: RunAnalyzer.hh:32
General exception.
Definition: JException.hh:23
JLimit_t numberOfTimeslices
Definition: RunAnalyzer.hh:34
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
Router for direct addressing of module data in detector data structure.
JLimit_t numberOfEvents
Definition: RunAnalyzer.hh:36
exit
Definition: JPizza.sh:36
string inputFile
Definition: RunAnalyzer.hh:29
Detector file.
Definition: JHead.hh:226
int modulesPerString
Definition: RunAnalyzer.hh:33
bool pmt_analysis
Definition: RunAnalyzer.hh:37
static const JStringCounter getNumberOfStrings
Function object to count unique strings.
void load(const std::string &file_name, JDetector &detector)
Load detector from input file.
JModuleRouter * router
Definition: RunAnalyzer.hh:31
int getNumberOfModules(const JDetector &detector)
Get number of modules.
RunAnalyzer::~RunAnalyzer ( )
inline

Destructor.

Definition at line 80 of file RunAnalyzer.hh.

80 {}

Member Function Documentation

void RunAnalyzer::iterateEventTree ( JTreeScanner< JDAQEvent > &  scanner)
inline

Definition at line 88 of file RunAnalyzer.hh.

88  {
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 
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());
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  }
TH1D * h_pmt_distribution_snapshot_hits
JRA_Histograms histograms
Definition: RunAnalyzer.hh:32
TH2D * h_Snapshot_hits_per_module
DAQ key hit.
Definition: JDAQKeyHit.hh:19
static const unsigned int NUMBER_OF_TRIGGER_BITS
Number of trigger bits.
TH1D * h_pmt_distribution_triggered_hits
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:68
TH2D * h_Triggered_hits_per_module
Template const_iterator.
Definition: JDAQEvent.hh:68
#define MAKE_STRING(A)
Make string.
Definition: JPrint.hh:127
TH1D * h_tot_distribution_snapshot_hits
TriggerHistograms h_trigger
bool pmt_analysis
Definition: RunAnalyzer.hh:37
#define FATAL(A)
Definition: JMessage.hh:67
int getString() const
Get string number.
Definition: JLocation.hh:134
General purpose string class.
Definition: JHead.hh:152
std::string to_string(const T &value)
Convert value to string.
bool hasModule(const JObjectID &id) const
Has module.
TH1D * h_tot_distribution_triggered_hits
JModuleRouter * router
Definition: RunAnalyzer.hh:31
JManager< string, TH2D > * m_Snapshot_hits_per_pmt
void RunAnalyzer::iterateSummarysliceTree ( JTreeScanner< JDAQSummaryslice > &  scanner)
inline

Definition at line 166 of file RunAnalyzer.hh.

166  {
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 
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  }
JRA_Histograms histograms
Definition: RunAnalyzer.hh:32
then fatal No hydrophone data file $HYDROPHONE_TXT fi sort gr k
int countFIFOStatus() const
Count FIFO status.
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:68
then set_variable singlesRate set_variable doublesRate set_variable numberOfSlices echo Generating random background echo Singles rate
JManager< string, TH2D > * m_mean_summary_rate
TH1D * h_dom_rate_distribution
#define MAKE_STRING(A)
Make string.
Definition: JPrint.hh:127
then awk string
int modulesPerString
Definition: RunAnalyzer.hh:33
TH1D * h_pmt_rate_distribution
bool pmt_analysis
Definition: RunAnalyzer.hh:37
#define FATAL(A)
Definition: JMessage.hh:67
SummaryHistograms h_summary
int getString() const
Get string number.
Definition: JLocation.hh:134
int countHighRateVeto() const
Count high-rate veto status.
JManager< string, TH2D > * m_summary_rate_distribution
std::string to_string(const T &value)
Convert value to string.
bool hasModule(const JObjectID &id) const
Has module.
int getCount(const T &hit)
Get hit count.
JModuleRouter * router
Definition: RunAnalyzer.hh:31
int j
Definition: JPolint.hh:703
static const int NUMBER_OF_PMTS
Total number of PMTs in module.
Definition: JDAQ.hh:26
bool testDAQStatus() const
Test DAQ status of packets.
template<class T >
void RunAnalyzer::iterateTimesliceTree ( JTreeScanner< T > &  scanner)
inline

Definition at line 279 of file RunAnalyzer.hh.

279  {
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  }
JRA_Histograms histograms
Definition: RunAnalyzer.hh:32
vector< JManager< string, TH2D > * > m_pmt_tot_distributions
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:68
vector< TH2D * > h_dom_mean_rates
then set_variable singlesRate set_variable doublesRate set_variable numberOfSlices echo Generating random background echo Singles rate
#define MAKE_STRING(A)
Make string.
Definition: JPrint.hh:127
TimesliceHistograms h_timeslice
Hit data structure.
Definition: JDAQHit.hh:34
double getFrameTime()
Get frame time duration.
Definition: JDAQClock.hh:162
do set_variable OUTPUT_DIRECTORY $WORKDIR T
then awk string
bool pmt_analysis
Definition: RunAnalyzer.hh:37
#define FATAL(A)
Definition: JMessage.hh:67
int getString() const
Get string number.
Definition: JLocation.hh:134
std::string to_string(const T &value)
Convert value to string.
vector< JManager< string, TH2D > * > m_pmt_rate_distributions
bool hasModule(const JObjectID &id) const
Has module.
JModuleRouter * router
Definition: RunAnalyzer.hh:31
int j
Definition: JPolint.hh:703
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
void RunAnalyzer::readSummaryData ( )
inline

Definition at line 345 of file RunAnalyzer.hh.

345  {
346 
348  if (scanner.hasNext()) {
350  scanner.rewind();
351  iterateSummarysliceTree(scanner);
352  }
353  }
JLimit_t numberOfSummaryslices
Definition: RunAnalyzer.hh:35
JRA_Histograms histograms
Definition: RunAnalyzer.hh:32
void iterateSummarysliceTree(JTreeScanner< JDAQSummaryslice > &scanner)
Definition: RunAnalyzer.hh:166
void initialize_summary_histograms()
Template definition for direct access of elements in ROOT TChain.
string inputFile
Definition: RunAnalyzer.hh:29
template<class T >
void RunAnalyzer::readTimesliceData ( )
inline

Definition at line 362 of file RunAnalyzer.hh.

362  {
363 
365  if(scanner.hasNext()) {
366 
368  scanner.rewind();
369  iterateTimesliceTree(scanner);
370  }
371  }
JRA_Histograms histograms
Definition: RunAnalyzer.hh:32
JLimit_t numberOfTimeslices
Definition: RunAnalyzer.hh:34
void iterateTimesliceTree(JTreeScanner< T > &scanner)
Definition: RunAnalyzer.hh:279
string inputFile
Definition: RunAnalyzer.hh:29
do set_variable OUTPUT_DIRECTORY $WORKDIR T
void initialize_timeslice_histograms()
void RunAnalyzer::readEvents ( )
inline

Definition at line 379 of file RunAnalyzer.hh.

379  {
380 
382 
383  if(scanner.hasNext()) {
384 
386  scanner.rewind();
387  iterateEventTree(scanner);
388  }
389  }
JRA_Histograms histograms
Definition: RunAnalyzer.hh:32
void iterateEventTree(JTreeScanner< JDAQEvent > &scanner)
Definition: RunAnalyzer.hh:88
JLimit_t numberOfEvents
Definition: RunAnalyzer.hh:36
Template definition for direct access of elements in ROOT TChain.
string inputFile
Definition: RunAnalyzer.hh:29
void initialize_trigger_histograms()
JRA_Histograms RunAnalyzer::getHistograms ( )
inline

Definition at line 394 of file RunAnalyzer.hh.

394  {
395 
396  return histograms ;
397  }
JRA_Histograms histograms
Definition: RunAnalyzer.hh:32

Member Data Documentation

string RunAnalyzer::inputFile
private

Definition at line 29 of file RunAnalyzer.hh.

JDetector RunAnalyzer::detector
private

Definition at line 30 of file RunAnalyzer.hh.

JModuleRouter* RunAnalyzer::router
private

Definition at line 31 of file RunAnalyzer.hh.

JRA_Histograms RunAnalyzer::histograms
private

Definition at line 32 of file RunAnalyzer.hh.

int RunAnalyzer::modulesPerString
private

Definition at line 33 of file RunAnalyzer.hh.

JLimit_t RunAnalyzer::numberOfTimeslices
private

Definition at line 34 of file RunAnalyzer.hh.

JLimit_t RunAnalyzer::numberOfSummaryslices
private

Definition at line 35 of file RunAnalyzer.hh.

JLimit_t RunAnalyzer::numberOfEvents
private

Definition at line 36 of file RunAnalyzer.hh.

bool RunAnalyzer::pmt_analysis
private

Definition at line 37 of file RunAnalyzer.hh.


The documentation for this class was generated from the following file: