Jpp 20.0.0-rc.1
the software that should make you happy
Loading...
Searching...
No Matches
JRunAnalyzer.hh
Go to the documentation of this file.
1#ifndef __JRUNANALYZER__
2#define __JRUNANALYZER__
3
7
13#include "JROOT/JManager.hh"
14#include "JTools/JQuantile.hh"
15#include "JTools/JRange.hh"
19#include "JMath/JConstants.hh"
20#include "JTrigger/JHitL0.hh"
21#include "JTrigger/JBuildL0.hh"
23#include "TFile.h"
24#include "TH1D.h"
25#include "JRunHistograms.hh"
26
27/**
28 * \author rgruiz, adomi
29 */
30
32
33/**
34 * Class dedicated to the analysis of %KM3NeT runs.
35 *
36 */
38
41 TFile *outputFile;
47
48public :
49
50 /**
51 * Constructor.
52 *
53 * \param file file name
54 * \param detector detector
55 * \param out pointer to output file
56 * \param nTimeslices number of timeslices to be read
57 * \param nSummaryslices number of summary slices to be read
58 * \param nEvents number of events to be read
59 * \param analysislevel option for analysis of PMT or reco data
60 */
61 JRunAnalyzer (const JSingleFileScanner<>& file, const JDetector& detector, TFile *out,
62 JLimit_t nTimeslices, JLimit_t nSummaryslices, JLimit_t nEvents, int analysislevel):
63 inputFile (file),
64 router (detector),
65 outputFile (out),
66 numberOfTimeslices (nTimeslices),
67 numberOfSummaryslices (nSummaryslices),
68 numberOfEvents (nEvents),
69 analysis_level (analysislevel)
70 {
71 using namespace JPP;
72 using namespace std;
73
74 histograms = JRA_Histograms(detector);
75 }
76
77 /**
78 * Destructor.
79 */
81
82 /*
83 * Function template used to read events from a tree.
84 *
85 * \param scanner A JTreeScanner
86 * \param
87 */
88 void iterateEventTree (JTreeScanner < JDAQEvent > & scanner) {
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 typedef JRange<double> JRange_t;
109 JRange_t range(JRange_t::DEFAULT_RANGE());
110
111 for (JDAQEvent::const_iterator<JDAQTriggeredHit> hit = event.begin<JDAQTriggeredHit>(); hit != event.end<JDAQTriggeredHit>(); ++hit) {
112
113 const JModule& module = router.getModule(hit->getModuleID());
114 const double t1 = getTime(hit->getT(), module.getPMT(hit->getPMT()));
115
116 range.include(t1);
117
118 JDAQTriggerMask trigger_mask(event.getTriggerMask(*hit));
119
120 if(trigger_mask.hasTriggerBit(getTriggerBit<JTrigger3DMuon>())) counter_3dmuon++;
121
122 if (router.hasModule(hit->getModuleID())) {
123
124 const JModule& module = router.getModule (hit->getModuleID());
125
126 histograms.h_trigger.h_pmt_distribution_triggered_hits->Fill(hit->getPMT());
127 histograms.h_trigger.h_tot_distribution_triggered_hits->Fill(hit->getToT());
128
129 int String = module.getString();
130 int StringIndex = histograms.string_to_index(String);
131 int Floor = module.getFloor();
132
133 histograms.h_trigger.h_Triggered_hits_per_module -> Fill(StringIndex , Floor);
134
136 histograms.h_trigger.h_Triggered_hits_3dmuon_per_module -> Fill(StringIndex , Floor);
137 }
138
139 for (unsigned int i = 0; i != NUMBER_OF_TRIGGER_BITS; ++i) {
140
141 if (hit -> hasTriggerBit(i)) {
142
143 histograms.h_trigger.h_Trigger_bit_hit->Fill((Double_t) i);
144 }
145 }
146 }else{
147 FATAL("JModuleRouter trying to access non existing identifier: "<< hit->getModuleID());
148 }
149 }
150
151 histograms.h_trigger.h_event_duration->Fill(range.getLength());
152
153 histograms.h_trigger.h_pmt_distribution_triggered_hits->Scale(1. / (Double_t) event.size<JDAQTriggeredHit> ());
154
155 if(counter_3dmuon != 0) histograms.h_trigger.h_Triggered_hits_3dmuon->Fill(counter_3dmuon);
156
157 for (JDAQEvent::const_iterator<JDAQSnapshotHit> hit = event.begin<JDAQSnapshotHit>() ; hit != event.end<JDAQSnapshotHit>() ; ++hit){
158
159 if (router.hasModule(hit->getModuleID())) {
160
161 const JModule& module = router.getModule (hit->getModuleID());
162
163 int String = module.getString();
164 int StringIndex = histograms.string_to_index(String);
165 int Floor = module.getFloor();
166 int pmt = hit-> getPMT();
167
168 if(analysis_level == 1){
169
170 (*histograms.h_trigger.m_Snapshot_hits_per_pmt)[MAKE_STRING("Detector/DU" + to_string(String))] -> Fill(pmt, Floor);
171
172 }
173
174 histograms.h_trigger.h_pmt_distribution_snapshot_hits -> Fill(hit->getPMT());
175 histograms.h_trigger.h_tot_distribution_snapshot_hits -> Fill(hit->getToT());
176 histograms.h_trigger.h_Snapshot_hits_per_module -> Fill(StringIndex, Floor);
177
178 }else{
179 FATAL("JModuleRouter trying to access non existing identifier: "<< hit->getModuleID());
180 }
181 }
182 }
183 }
184
185 /*
186 * Function template used to read reconstructed events from a tree.
187 *
188 * \param scanner A ParallelFileScanner for JDAQEvent and JEvt
189 * \param
190 */
192
193 TH1D* h_tres = new TH1D("h_tres", ";Time residuals [ns]; Entries", 100, -50, 150);
194 TH1D* h_likelihood = new TH1D ("h_likelihood", " ; Likelihood; Reconstructed Events ", 100, -1000, 1000);
195 TH1D* h_beta0 = new TH1D("h_beta0", "; beta0; Reconstructed Events", 20, 0, 1);
196 TH1D* h_energy = new TH1D ("h_energy", " ; Energy [GeV]; Reconstructed Events ", 65, 0, 9);
197 setLogarithmicX(h_energy);
198 TH1D* h_zenith = new TH1D("h_zenith", "; cosZenith; Reconstructed Events", 20, -1, 1);
199 TH1D* h_azimuth = new TH1D("h_azimuth", "; cosAzimuth; Reconstructed Events", 20, -1, 1);
200 TH1D* h_radial_position = new TH1D ("h_radial_position", "; Radial Position [m]; Reconstructed Events", 60, 0, 4.2);
201 setLogarithmicX(h_radial_position);
202 TH1D* h_z_position = new TH1D ("h_z_position", "; Z Position [m] ; Reconstructed Events", 50, 0, 3.2);
203 setLogarithmicX(h_z_position);
204
205 while (scanner.hasNext()) {
206
207 JParallelFileScanner< JTypeList<JDAQEvent, JFIT::JEvt> >::multi_pointer_type ps = scanner.next();
208
209 JEvt* evt = ps;
210 JDAQEvent* event = ps;
211
212 if (!evt->empty()) {
213
214 JEvt::iterator best = evt->begin();
215 h_energy -> Fill(best->getE());
216 h_likelihood -> Fill(best->getQ());
217 h_z_position -> Fill(best->getZ());
218 h_radial_position -> Fill(sqrt(best->getX()*best->getX() + best->getY()*best->getY()));
219
220 JTrack3D track = getTrack(*best);
221 h_zenith -> Fill(cos(track.getDirection().getTheta()));
222 h_azimuth -> Fill(cos(track.getDirection().getPhi()));
223 h_beta0 -> Fill(best->getW(JGANDALF_BETA0_RAD));
224
225 std::vector<JHitL0> dataL0;
226 const JBuildL0<JHitL0> buildL0;
227 buildL0(JDAQTimeslice(*event, false), router, back_inserter(dataL0));
228
229 for (std::vector<JHitL0>::const_iterator hit = dataL0.begin(); hit != dataL0.end(); ++hit) {
230
231 h_tres -> Fill(hit->getT() - track.getT(hit->getPosition()));
232 }
233 }
234 }
235
236 f.mkdir("Reco");
237 f.cd("Reco");
238
239 h_energy->Write();
240 h_likelihood -> Write();
241 h_z_position -> Write();
242 h_radial_position -> Write();
243 h_zenith -> Write();
244 h_azimuth -> Write();
245 h_beta0 -> Write();
246 h_tres -> Write();
247 }
248
249 struct dom_type :
250 public std::vector< JQuantile>
251 {
253 std::vector<JQuantile>(NUMBER_OF_PMTS)
254 {}
255 };
256
257 /*
258 * Function template to read and process time ordered summary slices.
259 *
260 * \param scanner A JTreeScanner from which time ordered summary slices are accessed.
261 */
262 void iterateSummarysliceTree (JTreeScanner < JDAQSummaryslice > & scanner) {
263
264 using namespace std;
265 using namespace JPP;
266
267 map<int, map<int, dom_type> > PMT_rate_quantiles;
268 map<int, map<int, JQuantile> > DOM_rate_quantiles;
269
270 while (scanner.hasNext()){
271
272 JDAQSummaryslice slice = *(scanner.next());
273
274 for (JDAQSummaryslice::const_iterator frame = slice.begin() ; frame != slice.end() ; ++ frame) {
275
276 if (router.hasModule(frame->getModuleID())) {
277
278 const JModule& module = router.getModule (frame->getModuleID());
279
280 int string = module.getString();
281 int stringIndex = histograms.string_to_index(string);
282 int floor = module.getFloor ();
283
284 JDAQFrameStatus status = frame -> getDAQFrameStatus();
285 int nFIFOcount = status.countFIFOStatus();
286 int nHRVcount = status.countHighRateVeto();
287
288 if (nHRVcount > 0) {
289 histograms.h_summary.h_hrv_per_dom->Fill(stringIndex , floor, nHRVcount);
290 }
291
292 if (status.testDAQStatus() == false) {
293 histograms.h_summary.h_daq_status_per_dom->Fill(stringIndex , floor, (1./distance(scanner.begin(), scanner.end()))*100);
294 }
295
296 histograms.h_summary.h_fifo_per_dom->Fill(stringIndex , floor , nFIFOcount);
297
298 if(analysis_level == 1){
299
300 TH2D* h2 = (*histograms.h_summary.m_summary_rate_distribution)[MAKE_STRING("Detector/DU" + to_string(string) + "/F" + to_string(floor))];
301
302 double rate = 0;
303 const double factor = 1.0e-3; // [kHz]
304
305 for (int i = 0 ; i < NUMBER_OF_PMTS ; i++){
306
307 rate += frame->getRate(i, factor);
308
309 h2->Fill(i , frame->getRate(i, factor));
310
311 PMT_rate_quantiles[string][floor][i].put(frame->getRate(i, factor));
312
313 histograms.h_summary.h_pmt_rate_distribution->Fill(frame->getRate(i, factor), frame->getWeight(i, factor));
314
315 }
316
317 histograms.h_summary.h_dom_rate_distribution->Fill(rate);
318
319 DOM_rate_quantiles[string][floor].put(rate);
320
321 } else {
322
323 double rate = 0;
324 const double factor = 1.0e-3; // [kHz]
325
326 for (int i = 0 ; i < NUMBER_OF_PMTS ; i++){
327
328 rate += frame->getRate(i, factor);
329
330 histograms.h_summary.h_pmt_rate_distribution->Fill(frame->getRate(i, factor), frame->getWeight(i, factor));
331 }
332
333 histograms.h_summary.h_dom_rate_distribution->Fill(rate);
334
335 DOM_rate_quantiles[string][floor].put(rate);
336
337 }
338
339 }else{
340 FATAL("JModuleRouter trying to access non existing identifier: "<< frame->getModuleID());
341 }
342 }
343 }
344
345 for (const auto& i1 : DOM_rate_quantiles) {
346 for (const auto& i2 : i1.second) {
347 if (i2.second.getCount() > 0) {
348 histograms.h_summary.h_rate_summary->Fill(histograms.string_to_index(i1.first), i2.first, i2.second.getMean());
349 }
350 }
351 }
352
353 if (analysis_level == 1){
354 for (const auto& i1 : PMT_rate_quantiles){
355
356 TH2D* h2 = (*histograms.h_summary.m_mean_summary_rate) [MAKE_STRING("Detector/DU" + to_string(i1.first))];
357 TH1D* h1 = (*histograms.h_summary.m_mean_summary_rate_distribution)[MAKE_STRING("Detector/DU" + to_string(i1.first))];
358
359 for (const auto& i2 : i1.second) {
360 for (int i3 = 0 ; i3 != NUMBER_OF_PMTS ; ++i3) {
361 if (i2.second[i3].getCount() > 0){
362 h2->Fill(i3, i2.first, i2.second[i3].getMean());
363 h1->Fill(i2.second[i3].getMean());
364 }
365 }
366 }
367 }
368 }
369 }
370
371 /*
372 * Function template to read and process time ordered timeslices.
373 *
374 * \param scanner A JTreeScanner from which time ordered timeslices are accessed.
375 */
376 template <class T>
377 void iterateTimesliceTree (JTreeScanner < T > & scanner) {
378
380
381 std::map < int , std::map <int , JQuantile> > DOM_rate_quantiles;
382
383 double inverseFrameTimeSec = 1. / (1.0e-9 * getFrameTime());
384
385 while (scanner.hasNext()){
386
387 T slice = *(scanner.next());
388
389 for(auto & frame : slice) {
390 if (router.hasModule(frame.getModuleID())) {
391
392 const JModule& module = router.getModule (frame.getModuleID());
393 double rate = frame.numberOfHits * inverseFrameTimeSec;
394 int string = module.getString();
395 int floor = module.getFloor ();
396
397 DOM_rate_quantiles[string][floor].put(rate);
398
399 vector <int> pmt_hit_count (NUMBER_OF_PMTS , 0) ;
400
401 if(analysis_level == 1){
402
403 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()))];
404
405 for (JDAQSuperFrame::const_iterator hit = frame.begin() ; hit != frame.end() ; ++hit){
406
407 h2 -> Fill(hit->getPMT() , hit->getToT()) ;
408
409 pmt_hit_count[hit->getPMT()]++;
410
411 }
412
413 for (int pmt = 0 ; pmt != NUMBER_OF_PMTS ; ++pmt) {
414
415 (*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);
416
417 }
418
419 }
420
421 }else{
422 FATAL("JModuleRouter trying to access non existing identifier: "<< frame.getModuleID());
423 }
424 }
425 }
426
427 for (std::map< int , std::map<int,JQuantile>>::const_iterator i = DOM_rate_quantiles.begin() ; i != DOM_rate_quantiles.end() ; i++){
428
429 for (std::map<int , JQuantile>::const_iterator j = i->second.begin() ; j != i->second.end() ; j++){
430
431 if (j->second.getCount() > 0) histograms.h_timeslice.h_dom_mean_rates[ts_type] -> Fill (histograms.string_to_index(i->first) , j->first , j->second.getMean() ) ;
432 }
433 }
434 }
435
436 /*
437 * Checks through the use of a JTreeScanner, if summary slices are present in the run file.
438 * In case they are, the corresponding tree is iterated.
439 *
440 */
442
443 JTreeScanner <JDAQSummaryslice> scanner(inputFile, numberOfSummaryslices);
444 if (scanner.hasNext()) {
446 scanner.rewind();
448 }
449 }
450
451 /*
452 * Function template that checks through the use of a JTreeScanner, if objects of different timeslice classes are present in the run file.
453 * In case they are, the corresponding tree is iterated.
454 *
455 */
456 template <class T>
458
459 JTreeScanner <T> scanner(inputFile, numberOfTimeslices);
460 if(scanner.hasNext()) {
461
462 histograms.initialize_timeslice_histograms <T>();
463 scanner.rewind();
464 iterateTimesliceTree(scanner);
465 }
466 }
467
468 /*
469 * Function template that checks through the use of a JTreeScanner, if JDAQEvent objects are present in the run file.
470 * In case they are, the corresponding tree is read.
471 *
472 */
473 void readEvents() {
474
475 JTreeScanner <JDAQEvent> scanner(inputFile, numberOfEvents);
476
477 if(scanner.hasNext()) {
478
480 scanner.rewind();
481 iterateEventTree(scanner);
482 }
483 }
484
485 /*
486 * Function template that checks through the use of a JTreeScanner, if JEvt objects are present in the run file.
487 * In case they are, the corresponding tree is read.
488 */
490
491 JParallelFileScanner< JTypeList<JDAQEvent, JFIT::JEvt> > scanner(inputFile, numberOfEvents);
492
493 if(scanner.hasNext()) {
494
495 scanner.rewind();
497 }
498 }
499
500 /*
501 * Writes the histograms to a root file.
502 * \param f The root file.
503 */
504 void writeToFile(TFile *f,int analysis_level){
505
506 f->mkdir("Detector");
507 f->cd("Detector");
508
509 if (histograms.h_summary.h_fifo_per_dom) histograms.h_summary.h_fifo_per_dom -> Write();
510 if (histograms.h_summary.h_daq_status_per_dom) histograms.h_summary.h_daq_status_per_dom -> Write();
511 if (histograms.h_summary.h_hrv_per_dom) histograms.h_summary.h_hrv_per_dom -> Write();
512 if (histograms.h_summary.h_rate_summary) histograms.h_summary.h_rate_summary -> Write();
513 if (histograms.h_summary.h_pmt_rate_distribution) histograms.h_summary.h_pmt_rate_distribution -> Write();
514 if (histograms.h_summary.h_dom_rate_distribution) histograms.h_summary.h_dom_rate_distribution -> Write();
515
516 histograms.Write_histogram_table_to_file(*f , MAKE_STRING("Detector"), histograms.h_timeslice.h_dom_mean_rates);
517
518 for (typename vector < JManager < string , TH2D >* >::const_iterator i = histograms.h_timeslice.m_pmt_tot_distributions.begin();
519 i != histograms.h_timeslice.m_pmt_tot_distributions.end() ; ++i){
520
521 if ((*i)){
522
523 for (typename JManager < string , TH2D >::const_iterator j = (*i) -> begin() ; j != (*i) -> end() ; ++j){
524
525 j->second->Scale(1., "width");
526 }
527 }
528 }
529
532
534
535 if (histograms.h_summary.m_mean_summary_rate) histograms.Write_manager_in_key_dir(*f , histograms.h_summary.m_mean_summary_rate);
537
539
540 histograms.Write_manager_table_in_key_dir (*f , histograms.h_timeslice.m_mean_ToT);
542
543 f->mkdir ( MAKE_STRING ("JDAQEvent").c_str());
544 f->cd ("JDAQEvent");
545
546 if (histograms.h_trigger.h_Trigger_bit_event) histograms.h_trigger.h_Trigger_bit_event -> Write();
547 if (histograms.h_trigger.h_Trigger_bit_hit) histograms.h_trigger.h_Trigger_bit_hit -> Write();
548 if (histograms.h_trigger.h_Triggered_hits) histograms.h_trigger.h_Triggered_hits -> Write();
549 if (histograms.h_trigger.h_Triggered_hits_3dmuon) histograms.h_trigger.h_Triggered_hits_3dmuon -> Write();
551 if (histograms.h_trigger.h_Snapshot_hits) histograms.h_trigger.h_Snapshot_hits -> Write();
553 if (histograms.h_trigger.h_event_duration) histograms.h_trigger.h_event_duration -> Write();
554 if (histograms.h_trigger.h_Number_of_overlays) histograms.h_trigger.h_Number_of_overlays -> Write();
558 if (histograms.h_trigger.h_tot_distribution_snapshot_hits) {histograms.h_trigger.h_tot_distribution_snapshot_hits->Scale(1, "width") ; histograms.h_trigger.h_tot_distribution_snapshot_hits -> Write();}
559 if (histograms.h_trigger.h_Triggered_hits_per_module) { histograms.h_trigger.h_Triggered_hits_per_module-> Write();}
560 if (histograms.h_trigger.h_Snapshot_hits_per_module) { histograms.h_trigger.h_Snapshot_hits_per_module -> Write();}
562
563 }
564};
565
566#endif
string outputFile
Basic data structure for L0 hit.
Dynamic ROOT object management.
Mathematical constants.
#define FATAL(A)
Definition JMessage.hh:67
Direct access to module in detector data structure.
Scanning of objects from multiple files according a format that follows from the extension of each fi...
Parallel scanning of objects from a single file or multiple files according a format that follows fro...
#define MAKE_STRING(A)
Make string.
Definition JPrint.hh:63
Auxiliary class to define a range between two values.
Support methods.
ROOT TTree parameter settings of various packages.
std::vector< T >::difference_type distance(typename std::vector< T >::const_iterator first, typename PhysicsEvent::const_iterator< T > second)
Specialisation of STL distance.
Detector data structure.
Definition JDetector.hh:96
Router for direct addressing of module data in detector data structure.
bool hasModule(const JObjectID &id) const
Has module.
const JModule & getModule(const JObjectID &id) const
Get module parameters.
Data structure for a composite optical module.
Definition JModule.hh:75
const JPMT & getPMT(const int index) const
Get PMT.
Definition JModule.hh:172
double getT(const JVector3D &pos) const
Get arrival time of Cherenkov light at given position.
Definition JTrack3D.hh:87
JVersor3D getDirection(const JVector3D &pos) const
Get photon direction of Cherenkov light on PMT.
Definition JTrack3D.hh:108
double getTheta() const
Get theta angle.
Definition JVersor3D.hh:128
double getPhi() const
Get phi angle.
Definition JVersor3D.hh:144
TimesliceHistograms h_timeslice
void initialize_trigger_histograms()
void Write_manager_table_in_key_dir(TFile &f, vector< JManager< T, V > * > table)
int string_to_index(int string)
TriggerHistograms h_trigger
void Write_manager_in_key_dir(TFile &f, JManager< T, V > *manager)
void initialize_summary_histograms()
void initialize_timeslice_histograms()
void Write_histogram_table_to_file(TFile &f, string dirname, vector< vector< T * > > table)
SummaryHistograms h_summary
Class dedicated to the analysis of KM3NeT runs.
void iterateSummarysliceTree(JTreeScanner< JDAQSummaryslice > &scanner)
JSingleFileScanner inputFile
JRA_Histograms histograms
JRunAnalyzer(const JSingleFileScanner<> &file, const JDetector &detector, TFile *out, JLimit_t nTimeslices, JLimit_t nSummaryslices, JLimit_t nEvents, int analysislevel)
Constructor.
void iterateEventTree(JTreeScanner< JDAQEvent > &scanner)
void iterateRecoEventTree(JParallelFileScanner< JTypeList< JDAQEvent, JFIT::JEvt > > &scanner, TFile &f)
JLimit_t numberOfTimeslices
TFile * outputFile
void readSummaryData()
void iterateTimesliceTree(JTreeScanner< T > &scanner)
void writeToFile(TFile *f, int analysis_level)
~JRunAnalyzer()
Destructor.
void readTimesliceData()
JLimit_t numberOfEvents
void readRecoEvents()
JModuleRouter router
JLimit_t numberOfSummaryslices
General purpose class for parallel reading of objects from a single file or multiple files.
virtual const multi_pointer_type & next() override
Get next element.
Object reading from a list of files.
virtual void rewind() override
Rewind.
virtual bool hasNext() override
Check availability of next element.
Range of values.
Definition JRange.hh:42
T getLength() const
Get length (difference between upper and lower limit).
Definition JRange.hh:289
range_type & include(argument_type x)
Include given value to range.
Definition JRange.hh:397
Template L0 hit builder.
Definition JBuildL0.hh:38
Template const_iterator.
Definition JDAQEvent.hh:68
int countHighRateVeto() const
Count high-rate veto status.
int countFIFOStatus() const
Count FIFO status.
bool testDAQStatus() const
Test DAQ status of packets.
Hit data structure.
Definition JDAQHit.hh:35
Auxiliary class for trigger mask.
bool hasTriggerBit(const unsigned int bit) const
Check trigger bit.
static const int JGANDALF_BETA0_RAD
KM3NeT Data Definitions v3.6.1-2-g905a24d https://git.km3net.de/common/km3net-dataformat.
JTrack3E getTrack(const Trk &track)
Get track.
void setLogarithmicX(TList *list)
Make x-axis of objects in list logarithmic (e.g. after using log10()).
std::string to_string(const T &value)
Convert value to string.
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
const char * getTime()
Get current local time conform ISO-8601 standard.
int j
Definition JPolint.hh:801
JTriggerbit_t getTriggerBit()
Get the trigger bit.
double getFrameTime()
Get frame time duration.
Definition JDAQClock.hh:162
static const unsigned int NUMBER_OF_TRIGGER_BITS
Number of trigger bits.
static const int NUMBER_OF_PMTS
Total number of PMTs in module.
Definition JDAQ.hh:26
Type definition of range.
Definition JHead.hh:43
General purpose string class.
Definition JHead.hh:152
Detector file.
Definition JHead.hh:227
Acoustic event fit.
Indexing of data type in type list.
Definition JTypeList.hh:310
Type list.
Definition JTypeList.hh:23
Auxiliary class for defining the range of iterations of objects.
Definition JLimit.hh:45
Auxiliary data structure for running average, standard deviation and quantiles.
Definition JQuantile.hh:41
JManager< string, TH1D > * m_mean_summary_rate_distribution
JManager< string, TH2D > * m_mean_summary_rate
JManager< string, TH2D > * m_summary_rate_distribution
vector< JManager< string, TH2D > * > m_pmt_tot_distributions
vector< JManager< string, TH2D > * > m_pmt_rate_distributions
vector< JManager< string, TH2D > * > m_mean_ToT
vector< TH2D * > h_dom_mean_rates
vector< JManager< string, TH1D > * > m_mean_ToT_distribution
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_pmt_distribution_snapshot_hits
TH2D * h_Snapshot_hits_per_module
TH1D * h_pmt_distribution_triggered_hits
TH1D * h_tot_distribution_triggered_hits
TH1D * h_Triggered_over_Snapshot_hits