Jpp 20.0.0-rc.2
the software that should make you happy
Loading...
Searching...
No Matches
JRunAnalyzer Class Reference

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

#include <JRunAnalyzer.hh>

Classes

struct  dom_type
 

Public Member Functions

 JRunAnalyzer (const JSingleFileScanner<> &file, const JDetector &detector, TFile *out, JLimit_t nTimeslices, JLimit_t nSummaryslices, JLimit_t nEvents, int analysislevel)
 Constructor.
 
 ~JRunAnalyzer ()
 Destructor.
 
void iterateEventTree (JTreeScanner< JDAQEvent > &scanner)
 
void iterateRecoEventTree (JParallelFileScanner< JTypeList< JDAQEvent, JFIT::JEvt > > &scanner, TFile &f)
 
void iterateSummarysliceTree (JTreeScanner< JDAQSummaryslice > &scanner)
 
template<class T >
void iterateTimesliceTree (JTreeScanner< T > &scanner)
 
void readSummaryData ()
 
template<class T >
void readTimesliceData ()
 
void readEvents ()
 
void readRecoEvents ()
 
void writeToFile (TFile *f, int analysis_level)
 

Private Attributes

JSingleFileScanner inputFile
 
JModuleRouter router
 
TFile * outputFile
 
JRA_Histograms histograms
 
JLimit_t numberOfTimeslices
 
JLimit_t numberOfSummaryslices
 
JLimit_t numberOfEvents
 
int analysis_level
 

Detailed Description

Class dedicated to the analysis of KM3NeT runs.

Definition at line 37 of file JRunAnalyzer.hh.

Constructor & Destructor Documentation

◆ JRunAnalyzer()

JRunAnalyzer::JRunAnalyzer ( const JSingleFileScanner<> & file,
const JDetector & detector,
TFile * out,
JLimit_t nTimeslices,
JLimit_t nSummaryslices,
JLimit_t nEvents,
int analysislevel )
inline

Constructor.

Parameters
filefile name
detectordetector
outpointer to output file
nTimeslicesnumber of timeslices to be read
nSummaryslicesnumber of summary slices to be read
nEventsnumber of events to be read
analysisleveloption for analysis of PMT or reco data

Definition at line 61 of file JRunAnalyzer.hh.

62 :
63 inputFile (file),
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
75 }
string outputFile
JSingleFileScanner inputFile
JRA_Histograms histograms
JLimit_t numberOfTimeslices
JLimit_t numberOfEvents
JModuleRouter router
JLimit_t numberOfSummaryslices
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
Detector file.
Definition JHead.hh:227

◆ ~JRunAnalyzer()

JRunAnalyzer::~JRunAnalyzer ( )
inline

Destructor.

Definition at line 80 of file JRunAnalyzer.hh.

80{}

Member Function Documentation

◆ iterateEventTree()

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 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
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
135 if(trigger_mask.hasTriggerBit(getTriggerBit<JTrigger3DMuon>())){
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
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 }
#define FATAL(A)
Definition JMessage.hh:67
#define MAKE_STRING(A)
Make string.
Definition JPrint.hh:63
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
int string_to_index(int string)
TriggerHistograms h_trigger
Range of values.
Definition JRange.hh:42
Template const_iterator.
Definition JDAQEvent.hh:68
Auxiliary class for trigger mask.
std::string to_string(const T &value)
Convert value to string.
const char * getTime()
Get current local time conform ISO-8601 standard.
JTriggerbit_t getTriggerBit()
Get the trigger bit.
static const unsigned int NUMBER_OF_TRIGGER_BITS
Number of trigger bits.
Type definition of range.
Definition JHead.hh:43
General purpose string class.
Definition JHead.hh:152
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

◆ iterateRecoEventTree()

void JRunAnalyzer::iterateRecoEventTree ( JParallelFileScanner< JTypeList< JDAQEvent, JFIT::JEvt > > & scanner,
TFile & f )
inline

Definition at line 191 of file JRunAnalyzer.hh.

191 {
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 }
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
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.
virtual bool hasNext() override
Check availability of next element.
Template L0 hit builder.
Definition JBuildL0.hh:38
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()).
Acoustic event fit.

◆ iterateSummarysliceTree()

void JRunAnalyzer::iterateSummarysliceTree ( JTreeScanner< JDAQSummaryslice > & scanner)
inline

Definition at line 262 of file JRunAnalyzer.hh.

262 {
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
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
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 }
std::vector< T >::difference_type distance(typename std::vector< T >::const_iterator first, typename PhysicsEvent::const_iterator< T > second)
Specialisation of STL distance.
SummaryHistograms h_summary
int countHighRateVeto() const
Count high-rate veto status.
int countFIFOStatus() const
Count FIFO status.
bool testDAQStatus() const
Test DAQ status of packets.
static const int NUMBER_OF_PMTS
Total number of PMTs in module.
Definition JDAQ.hh:26
JManager< string, TH1D > * m_mean_summary_rate_distribution
JManager< string, TH2D > * m_mean_summary_rate
JManager< string, TH2D > * m_summary_rate_distribution

◆ iterateTimesliceTree()

template<class T >
void JRunAnalyzer::iterateTimesliceTree ( JTreeScanner< T > & scanner)
inline

Definition at line 377 of file JRunAnalyzer.hh.

377 {
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 }
TimesliceHistograms h_timeslice
Hit data structure.
Definition JDAQHit.hh:35
int j
Definition JPolint.hh:801
double getFrameTime()
Get frame time duration.
Definition JDAQClock.hh:162
Indexing of data type in type list.
Definition JTypeList.hh:310
vector< JManager< string, TH2D > * > m_pmt_tot_distributions
vector< JManager< string, TH2D > * > m_pmt_rate_distributions
vector< TH2D * > h_dom_mean_rates

◆ readSummaryData()

void JRunAnalyzer::readSummaryData ( )
inline

Definition at line 441 of file JRunAnalyzer.hh.

441 {
442
443 JTreeScanner <JDAQSummaryslice> scanner(inputFile, numberOfSummaryslices);
444 if (scanner.hasNext()) {
446 scanner.rewind();
448 }
449 }
void initialize_summary_histograms()
void iterateSummarysliceTree(JTreeScanner< JDAQSummaryslice > &scanner)

◆ readTimesliceData()

template<class T >
void JRunAnalyzer::readTimesliceData ( )
inline

Definition at line 457 of file JRunAnalyzer.hh.

457 {
458
459 JTreeScanner <T> scanner(inputFile, numberOfTimeslices);
460 if(scanner.hasNext()) {
461
463 scanner.rewind();
464 iterateTimesliceTree(scanner);
465 }
466 }
void initialize_timeslice_histograms()
void iterateTimesliceTree(JTreeScanner< T > &scanner)

◆ readEvents()

void JRunAnalyzer::readEvents ( )
inline

Definition at line 473 of file JRunAnalyzer.hh.

473 {
474
475 JTreeScanner <JDAQEvent> scanner(inputFile, numberOfEvents);
476
477 if(scanner.hasNext()) {
478
480 scanner.rewind();
481 iterateEventTree(scanner);
482 }
483 }
void initialize_trigger_histograms()
void iterateEventTree(JTreeScanner< JDAQEvent > &scanner)

◆ readRecoEvents()

void JRunAnalyzer::readRecoEvents ( )
inline

Definition at line 489 of file JRunAnalyzer.hh.

489 {
490
491 JParallelFileScanner< JTypeList<JDAQEvent, JFIT::JEvt> > scanner(inputFile, numberOfEvents);
492
493 if(scanner.hasNext()) {
494
495 scanner.rewind();
497 }
498 }
void iterateRecoEventTree(JParallelFileScanner< JTypeList< JDAQEvent, JFIT::JEvt > > &scanner, TFile &f)

◆ writeToFile()

void JRunAnalyzer::writeToFile ( TFile * f,
int analysis_level )
inline

Definition at line 504 of file JRunAnalyzer.hh.

504 {
505
506 f->mkdir("Detector");
507 f->cd("Detector");
508
515
517
518 for (typename vector < JManager < string , TH2D >* >::const_iterator i = histograms.h_timeslice.m_pmt_tot_distributions.begin();
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
537
539
542
543 f->mkdir ( MAKE_STRING ("JDAQEvent").c_str());
544 f->cd ("JDAQEvent");
545
562
563 }
void Write_manager_table_in_key_dir(TFile &f, vector< JManager< T, V > * > table)
void Write_manager_in_key_dir(TFile &f, JManager< T, V > *manager)
void Write_histogram_table_to_file(TFile &f, string dirname, vector< vector< T * > > table)
vector< JManager< string, TH2D > * > m_mean_ToT
vector< JManager< string, TH1D > * > m_mean_ToT_distribution

Member Data Documentation

◆ inputFile

JSingleFileScanner JRunAnalyzer::inputFile
private

Definition at line 39 of file JRunAnalyzer.hh.

◆ router

JModuleRouter JRunAnalyzer::router
private

Definition at line 40 of file JRunAnalyzer.hh.

◆ outputFile

TFile* JRunAnalyzer::outputFile
private

Definition at line 41 of file JRunAnalyzer.hh.

◆ histograms

JRA_Histograms JRunAnalyzer::histograms
private

Definition at line 42 of file JRunAnalyzer.hh.

◆ numberOfTimeslices

JLimit_t JRunAnalyzer::numberOfTimeslices
private

Definition at line 43 of file JRunAnalyzer.hh.

◆ numberOfSummaryslices

JLimit_t JRunAnalyzer::numberOfSummaryslices
private

Definition at line 44 of file JRunAnalyzer.hh.

◆ numberOfEvents

JLimit_t JRunAnalyzer::numberOfEvents
private

Definition at line 45 of file JRunAnalyzer.hh.

◆ analysis_level

int JRunAnalyzer::analysis_level
private

Definition at line 46 of file JRunAnalyzer.hh.


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