Jpp  16.0.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
JRunAnalyzer Class Reference

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

#include <JRunAnalyzer.hh>

Public Member Functions

 JRunAnalyzer ()
 Default constructor. More...
 
 JRunAnalyzer (string file, string detectorFile, JLimit_t nTimeslices, JLimit_t nSummaryslices, JLimit_t nEvents, bool pmtanalysis)
 Constructor. More...
 
 ~JRunAnalyzer ()
 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 JRunAnalyzer.hh.

Constructor & Destructor Documentation

JRunAnalyzer::JRunAnalyzer ( )
inline

Default constructor.

Definition at line 44 of file JRunAnalyzer.hh.

44 {}
JRunAnalyzer::JRunAnalyzer ( 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
pmtanalysisOption for analysis of PMT data

Definition at line 56 of file JRunAnalyzer.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  }
General exception.
Definition: JException.hh:23
then set_variable PMT_FILE set_variable DAQ_FILE set_variable OUTPUT_FILE set_variable DETECTOR else fatal Wrong number of arguments fi set_variable RUNBYRUN file
int modulesPerString
Definition: JRunAnalyzer.hh:33
JRA_Histograms histograms
Definition: JRunAnalyzer.hh:32
Router for direct addressing of module data in detector data structure.
exit
Definition: JPizza.sh:36
string inputFile
Definition: JRunAnalyzer.hh:29
JLimit_t numberOfTimeslices
Definition: JRunAnalyzer.hh:34
JModuleRouter * router
Definition: JRunAnalyzer.hh:31
Detector file.
Definition: JHead.hh:224
static const JStringCounter getNumberOfStrings
Function object to count unique strings.
void load(const std::string &file_name, JDetector &detector)
Load detector from input file.
JLimit_t numberOfEvents
Definition: JRunAnalyzer.hh:36
JLimit_t numberOfSummaryslices
Definition: JRunAnalyzer.hh:35
int getNumberOfModules(const JDetector &detector)
Get number of modules.
JRunAnalyzer::~JRunAnalyzer ( )
inline

Destructor.

Definition at line 80 of file JRunAnalyzer.hh.

80 {}

Member Function Documentation

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

Definition at line 88 of file JRunAnalyzer.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_Triggered_over_Snapshot_hits -> Fill((Double_t) event.size<JDAQTriggeredHit> () / event.size<JDAQSnapshotHit > ());
97  histograms.h_trigger.h_Number_of_overlays -> Fill(event.getOverlays());
98 
99  for (unsigned int i = 0; i != NUMBER_OF_TRIGGER_BITS; ++i) {
100 
101  if (event.hasTriggerBit(i)) {
102 
103  histograms.h_trigger.h_Trigger_bit_event -> Fill((Double_t) i);
104  }
105  }
106 
107  int counter_3dmuon = 0;
108 
109  for (JDAQEvent::const_iterator<JDAQTriggeredHit> hit = event.begin<JDAQTriggeredHit>(); hit != event.end<JDAQTriggeredHit>(); ++hit) {
110 
111  JDAQTriggerMask trigger_mask(event.getTriggerMask(*hit));
112 
113  if(trigger_mask.hasTriggerBit(getTriggerBit<JTrigger3DMuon>())) counter_3dmuon++;
114 
115  if (router->hasModule(hit->getModuleID())) {
116 
117  const JModule& module = router->getModule (hit->getModuleID());
118 
121 
122  int String = module.getString();
123  int Floor = module.getFloor();
124 
125  histograms.h_trigger.h_Triggered_hits_per_module -> Fill(String , Floor);
126 
127  if(trigger_mask.hasTriggerBit(getTriggerBit<JTrigger3DMuon>())){
129  }
130 
131  for (unsigned int i = 0; i != NUMBER_OF_TRIGGER_BITS; ++i) {
132 
133  if (hit -> hasTriggerBit(i)) {
134 
135  histograms.h_trigger.h_Trigger_bit_hit->Fill((Double_t) i);
136  }
137  }
138  }else{
139  FATAL("JModuleRouter trying to access non existing identifier: "<< hit->getModuleID());
140  }
141  }
142 
143  if(counter_3dmuon != 0) histograms.h_trigger.h_Triggered_hits_3dmuon->Fill(counter_3dmuon);
144 
145  for (JDAQEvent::const_iterator<JDAQSnapshotHit> hit = event.begin<JDAQSnapshotHit>() ; hit != event.end<JDAQSnapshotHit>() ; ++hit){
146 
147  if (router->hasModule(hit->getModuleID())) {
148 
149  const JModule& module = router->getModule (hit->getModuleID());
150 
151  int String = module.getString();
152  int Floor = module.getFloor();
153  int pmt = hit-> getPMT();
154 
155  if(pmt_analysis == true){
156 
157  (*histograms.h_trigger.m_Snapshot_hits_per_pmt)[MAKE_STRING("Detector/DU" + to_string(String))] -> Fill(pmt, Floor);
158 
159  }
160 
161  histograms.h_trigger.h_pmt_distribution_snapshot_hits -> Fill(hit->getPMT());
163  histograms.h_trigger.h_Snapshot_hits_per_module -> Fill(String, Floor);
164 
165 
166  }else{
167  FATAL("JModuleRouter trying to access non existing identifier: "<< hit->getModuleID());
168  }
169  }
170  }
171  }
TH1D * h_pmt_distribution_snapshot_hits
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
Auxiliary class for trigger mask.
JRA_Histograms histograms
Definition: JRunAnalyzer.hh:32
TH2D * h_Triggered_hits_per_module
Template const_iterator.
Definition: JDAQEvent.hh:68
JModuleRouter * router
Definition: JRunAnalyzer.hh:31
#define MAKE_STRING(A)
Make string.
Definition: JPrint.hh:142
TH1D * h_tot_distribution_snapshot_hits
TriggerHistograms h_trigger
#define FATAL(A)
Definition: JMessage.hh:67
int getString() const
Get string number.
Definition: JLocation.hh:134
General purpose string class.
Definition: JHead.hh:150
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
TH1D * h_Triggered_hits_3dmuon
TH2D * h_Triggered_hits_3dmuon_per_module
JManager< string, TH2D > * m_Snapshot_hits_per_pmt
TH1D * h_Triggered_over_Snapshot_hits
void JRunAnalyzer::iterateSummarysliceTree ( JTreeScanner< JDAQSummaryslice > &  scanner)
inline

Definition at line 179 of file JRunAnalyzer.hh.

179  {
180 
181  std::map <int,vector<vector <JQuantile>>> PMT_rate_quantiles;
182  std::map <int,vector<JQuantile>> DOM_rate_quantiles;
183 
184  while (scanner.hasNext()){
185 
186  JDAQSummaryslice slice = *(scanner.next());
187 
188  for (JDAQSummaryslice::const_iterator frame = slice.begin() ; frame != slice.end() ; ++ frame) {
189 
190  if (router->hasModule(frame->getModuleID())) {
191 
192  const JModule& module = router->getModule (frame->getModuleID());
193 
194  int string = module.getString();
195  int floor = module.getFloor ();
196 
197  PMT_rate_quantiles[string].resize(modulesPerString , vector <JQuantile>(NUMBER_OF_PMTS));
198 
199  JDAQFrameStatus status = frame -> getDAQFrameStatus();
200  int nFIFOcount = status.countFIFOStatus();
201  int nHRVcount = status.countHighRateVeto();
202 
203  if (nHRVcount > 0) {
204  histograms.h_summary.h_hrv_per_dom->Fill(string , floor, nHRVcount);
205  }
206 
207  if (status.testDAQStatus() == false) {
208  histograms.h_summary.h_daq_status_per_dom->Fill(string , floor);
209  }
210 
211  histograms.h_summary.h_fifo_per_dom->Fill(string , floor , nFIFOcount);
212 
213  if(pmt_analysis == true){
214 
215  TH2D* h2 = (*histograms.h_summary.m_summary_rate_distribution)[MAKE_STRING("Detector/DU" + to_string(string) + "/F" + to_string(floor))];
216 
217  double rate = 0;
218  const double factor = 1.0e-3; // [kHz]
219 
220  for (int i = 0 ; i < NUMBER_OF_PMTS ; i++){
221 
222  rate += frame->getRate(i, factor);
223 
224  h2->Fill(i , frame->getRate(i));
225 
226  PMT_rate_quantiles[string][floor-1][i].put(frame->getRate(i, factor));
227 
228  histograms.h_summary.h_pmt_rate_distribution->Fill(frame->getRate(i, factor), frame->getWeight(i, factor));
229 
230  }
231 
233 
234  DOM_rate_quantiles[string].resize(modulesPerString);
235 
236  DOM_rate_quantiles[string][floor-1].put(rate * 1e-3);
237 
238  } else {
239 
240  double rate = 0;
241  const double factor = 1.0e-3; // [kHz]
242 
243  for (int i = 0 ; i < NUMBER_OF_PMTS ; i++){
244 
245  rate += frame->getRate(i, factor);
246 
247  histograms.h_summary.h_pmt_rate_distribution->Fill(frame->getRate(i, factor), frame->getWeight(i, factor));
248  }
249 
251 
252  DOM_rate_quantiles[string].resize(modulesPerString);
253 
254  DOM_rate_quantiles[string][floor-1].put(rate);
255 
256  }
257 
258  }else{
259  FATAL("JModuleRouter trying to access non existing identifier: "<< frame->getModuleID());
260  }
261  }
262  }
263 
264 
265  for (std::map<int , vector<JQuantile> >::const_iterator i = DOM_rate_quantiles.begin() ; i!= DOM_rate_quantiles.end() ; ++i) {
266 
267  for (int j=1 ; j < modulesPerString + 1 ; j++) {
268 
269  if (i->second[j-1].getCount() > 0) histograms.h_summary.h_rate_summary->Fill(i->first, j, i->second[j-1].getMean());
270  }
271  }
272 
273 
274  if(pmt_analysis == true){
275 
276  for (std::map<int , vector < vector<JQuantile> > >::const_iterator i = PMT_rate_quantiles.begin() ; i!= PMT_rate_quantiles.end() ; ++i) {
277 
278  for (int j=0 ; j < modulesPerString ; j++){
279 
280  for (int k=0 ; k < NUMBER_OF_PMTS ; k++){
281 
282  if (i -> second[j][k].getCount() > 0){
283 
284  (*histograms.h_summary.m_mean_summary_rate) [MAKE_STRING("Detector/DU" + to_string(i->first))]->Fill(k, j+1, i->second[j][k].getMean());
285  (*histograms.h_summary.m_mean_summary_rate_distribution)[MAKE_STRING("Detector/DU" + to_string(i->first))]->Fill( i->second[j][k].getMean());
286  }
287  }
288  }
289  }
290  }
291  }
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
int modulesPerString
Definition: JRunAnalyzer.hh:33
JRA_Histograms histograms
Definition: JRunAnalyzer.hh:32
then set_variable singlesRate set_variable doublesRate set_variable numberOfSlices echo Generating random background echo Singles rate
JModuleRouter * router
Definition: JRunAnalyzer.hh:31
JManager< string, TH2D > * m_mean_summary_rate
TH1D * h_dom_rate_distribution
#define MAKE_STRING(A)
Make string.
Definition: JPrint.hh:142
TH1D * h_pmt_rate_distribution
#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.
int j
Definition: JPolint.hh:682
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 JRunAnalyzer::iterateTimesliceTree ( JTreeScanner< T > &  scanner)
inline

Definition at line 300 of file JRunAnalyzer.hh.

300  {
301 
303 
304  std::map < int , std::map <int , JQuantile> > DOM_rate_quantiles;
305 
306  double inverseFrameTimeSec = 1. / (1.0e-9 * getFrameTime());
307 
308  while (scanner.hasNext()){
309 
310  T slice = *(scanner.next());
311 
312  for(auto & frame : slice) {
313  if (router->hasModule(frame.getModuleID())) {
314 
315  const JModule& module = router->getModule (frame.getModuleID());
316  double rate = frame.numberOfHits * inverseFrameTimeSec;
317  int string = module.getString();
318  int floor = module.getFloor ();
319 
320  DOM_rate_quantiles[string][floor].put(rate);
321 
322  vector <int> pmt_hit_count (NUMBER_OF_PMTS , 0) ;
323 
324  if(pmt_analysis == true){
325 
326  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()))];
327 
328  for (JDAQSuperFrame::const_iterator hit = frame.begin() ; hit != frame.end() ; ++hit){
329 
330  h2 -> Fill(hit->getPMT() , hit->getToT()) ;
331 
332  pmt_hit_count[hit->getPMT()]++;
333 
334  }
335 
336  for (int pmt = 0 ; pmt != NUMBER_OF_PMTS ; ++pmt) {
337 
338  (*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);
339 
340  }
341 
342  }
343 
344  }else{
345  FATAL("JModuleRouter trying to access non existing identifier: "<< frame.getModuleID());
346  }
347  }
348  }
349 
350 
351  for (std::map< int , std::map<int,JQuantile>>::const_iterator i = DOM_rate_quantiles.begin() ; i != DOM_rate_quantiles.end() ; i++){
352 
353  for (std::map<int , JQuantile>::const_iterator j = i->second.begin() ; j != i->second.end() ; j++){
354 
355  if (j->second.getCount() > 0) histograms.h_timeslice.h_dom_mean_rates[ts_type] -> Fill (i->first , j->first , j->second.getMean() ) ;
356  }
357  }
358  }
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
JRA_Histograms histograms
Definition: JRunAnalyzer.hh:32
then set_variable singlesRate set_variable doublesRate set_variable numberOfSlices echo Generating random background echo Singles rate
JModuleRouter * router
Definition: JRunAnalyzer.hh:31
#define MAKE_STRING(A)
Make string.
Definition: JPrint.hh:142
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
#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.
int j
Definition: JPolint.hh:682
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 JRunAnalyzer::readSummaryData ( )
inline

Definition at line 366 of file JRunAnalyzer.hh.

366  {
367 
369  if (scanner.hasNext()) {
371  scanner.rewind();
372  iterateSummarysliceTree(scanner);
373  }
374  }
JRA_Histograms histograms
Definition: JRunAnalyzer.hh:32
string inputFile
Definition: JRunAnalyzer.hh:29
void initialize_summary_histograms()
void iterateSummarysliceTree(JTreeScanner< JDAQSummaryslice > &scanner)
Template definition for direct access of elements in ROOT TChain.
Definition: JTreeScanner.hh:91
JLimit_t numberOfSummaryslices
Definition: JRunAnalyzer.hh:35
template<class T >
void JRunAnalyzer::readTimesliceData ( )
inline

Definition at line 383 of file JRunAnalyzer.hh.

383  {
384 
386  if(scanner.hasNext()) {
387 
389  scanner.rewind();
390  iterateTimesliceTree(scanner);
391  }
392  }
void iterateTimesliceTree(JTreeScanner< T > &scanner)
JRA_Histograms histograms
Definition: JRunAnalyzer.hh:32
string inputFile
Definition: JRunAnalyzer.hh:29
JLimit_t numberOfTimeslices
Definition: JRunAnalyzer.hh:34
Template definition for direct access of elements in ROOT TChain.
Definition: JTreeScanner.hh:91
do set_variable OUTPUT_DIRECTORY $WORKDIR T
void initialize_timeslice_histograms()
void JRunAnalyzer::readEvents ( )
inline

Definition at line 400 of file JRunAnalyzer.hh.

400  {
401 
403 
404  if(scanner.hasNext()) {
405 
407  scanner.rewind();
408  iterateEventTree(scanner);
409  }
410  }
JRA_Histograms histograms
Definition: JRunAnalyzer.hh:32
string inputFile
Definition: JRunAnalyzer.hh:29
Template definition for direct access of elements in ROOT TChain.
Definition: JTreeScanner.hh:91
void iterateEventTree(JTreeScanner< JDAQEvent > &scanner)
Definition: JRunAnalyzer.hh:88
JLimit_t numberOfEvents
Definition: JRunAnalyzer.hh:36
void initialize_trigger_histograms()
JRA_Histograms JRunAnalyzer::getHistograms ( )
inline

Definition at line 415 of file JRunAnalyzer.hh.

415  {
416 
417  return histograms ;
418  }
JRA_Histograms histograms
Definition: JRunAnalyzer.hh:32

Member Data Documentation

string JRunAnalyzer::inputFile
private

Definition at line 29 of file JRunAnalyzer.hh.

JDetector JRunAnalyzer::detector
private

Definition at line 30 of file JRunAnalyzer.hh.

JModuleRouter* JRunAnalyzer::router
private

Definition at line 31 of file JRunAnalyzer.hh.

JRA_Histograms JRunAnalyzer::histograms
private

Definition at line 32 of file JRunAnalyzer.hh.

int JRunAnalyzer::modulesPerString
private

Definition at line 33 of file JRunAnalyzer.hh.

JLimit_t JRunAnalyzer::numberOfTimeslices
private

Definition at line 34 of file JRunAnalyzer.hh.

JLimit_t JRunAnalyzer::numberOfSummaryslices
private

Definition at line 35 of file JRunAnalyzer.hh.

JLimit_t JRunAnalyzer::numberOfEvents
private

Definition at line 36 of file JRunAnalyzer.hh.

bool JRunAnalyzer::pmt_analysis
private

Definition at line 37 of file JRunAnalyzer.hh.


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