Jpp  18.0.0-rc.2
the software that should make you happy
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
JTriggerEfficiency2D.cc
Go to the documentation of this file.
1 
2 #include <string>
3 #include <iostream>
4 
5 #include "TROOT.h"
6 #include "TFile.h"
7 #include "TH2D.h"
8 
11 
12 #include "JDAQ/JDAQEventIO.hh"
13 
14 #include "JAAnet/JAAnetToolkit.hh"
15 
17 
20 #include "JSupport/JSupport.hh"
21 
22 #include "Jeep/JParser.hh"
23 #include "Jeep/JMessage.hh"
24 
25 
26 /**
27  * \file
28  *
29  * Example program to histogram trigger efficiency.
30  * \author mdejong, bjjung
31  */
32 int main(int argc, char **argv)
33 {
34  using namespace std;
35  using namespace JPP;
36  using namespace KM3NETDAQ;
37 
38  typedef JAbstractHistogram<double> JAbstractHistogram_t;
39 
40  JTriggeredFileScanner<> inputFile;
41  string outputFile;
42  JAbstractHistogram_t X;
43  JAbstractHistogram_t Y;
44  int debug;
45 
46  try {
47 
48  JParser<> zap("Example program to histogram trigger efficiency.");
49 
50  zap['f'] = make_field(inputFile);
51  zap['o'] = make_field(outputFile) = "efficiency.root";
52  zap['X'] = make_field(X, "Energy binning") = JAbstractHistogram_t(20, 1e-3, 100.0);
53  zap['Y'] = make_field(Y, "Zenith-angle binning") = JAbstractHistogram_t(20, -1.0, 1.0);
54  zap['d'] = make_field(debug) = 1;
55 
56  zap(argc, argv);
57  }
58  catch(const exception &error) {
59  FATAL(error.what() << endl);
60  }
61 
62 
63  // output
64 
65  TFile out(outputFile.c_str(), "recreate");
66 
67  TH2D h0("h0", NULL,
68  X.getNumberOfBins(), X.getLowerLimit(), X.getUpperLimit(),
69  Y.getNumberOfBins(), Y.getLowerLimit(), Y.getUpperLimit());
70  TH2D h1("h1", NULL,
71  X.getNumberOfBins(), X.getLowerLimit(), X.getUpperLimit(),
72  Y.getNumberOfBins(), Y.getLowerLimit(), Y.getUpperLimit());
73  TH2D h2("h2", NULL,
74  X.getNumberOfBins(), X.getLowerLimit(), X.getUpperLimit(),
75  Y.getNumberOfBins(), Y.getLowerLimit(), Y.getUpperLimit());
76 
77  h2.Sumw2();
78 
79 
80  while (inputFile.hasNext()) {
81 
82  STATUS("event: " << setw(10) << inputFile.getCounter() << '\r'); DEBUG(endl);
83 
84  JTriggeredFileScanner<>::multi_pointer_type ps = inputFile.next();
85 
86  const Evt* event = ps;
87 
88  vector<Hit>::const_iterator hit = find_if_not(event->mc_hits.cbegin(), event->mc_hits.cend(), &is_noise);
89 
90  vector<Trk>::const_iterator primary = find_if(event->mc_trks.cbegin(), event->mc_trks.cend(), &is_initialstate);
91 
92  if (hit != event->mc_hits.cend() && primary != event->mc_trks.cend()) {
93  h1.Fill(primary->E, primary->dir.z, 1.0);
94  }
95  }
96  STATUS(endl);
97 
98 
99  for (JMultipleFileScanner<Evt> in(inputFile); in.hasNext(); ) {
100 
101  STATUS("event: " << setw(10) << in.getCounter() << '\r'); DEBUG(endl);
102 
103  const Evt* event = in.next();
104 
105  vector<Hit>::const_iterator hit = find_if_not(event->mc_hits.cbegin(), event->mc_hits.cend(), &is_noise);
106 
107  vector<Trk>::const_iterator primary = find_if(event->mc_trks.cbegin(), event->mc_trks.cend(), &is_initialstate);
108 
109  if (hit != event->mc_hits.cend() && primary != event->mc_trks.cend()) {
110  h0.Fill(primary->E, primary->dir.z, 1.0);
111  }
112  }
113  STATUS(endl);
114 
115 
116  for (Int_t i = 1; i <= h1.GetNbinsX(); ++i) {
117  for (Int_t j = 1; j <= h1.GetNbinsY(); ++j) {
118 
119  const Double_t z1 = h1.GetBinContent(i,j);
120  const Double_t z0 = h0.GetBinContent(i,j);
121 
122  if (z0 != 0.0) {
123 
124  const Double_t z3 = z1 / z0;
125  const Double_t w3 = sqrt(z1 * (z0 - z1) / (z0*z0*z0));
126 
127  h2.SetBinContent(i, j, z3);
128  h2.SetBinError (i, j, w3);
129 
130  } else {
131  ERROR("Bin " << h0.GetName() << "[" << i << "," << j << "] empty." << endl);
132  }
133  }
134  }
135 
136  out.Write();
137  out.Close();
138 }
139 
140 
Utility class to parse command line options.
Definition: JParser.hh:1514
int main(int argc, char *argv[])
Definition: Main.cc:15
ROOT TTree parameter settings of various packages.
Synchronously read DAQ events and Monte Carlo events (and optionally other events).
#define STATUS(A)
Definition: JMessage.hh:63
Auxiliary class to synchronously read DAQ events and Monte Carlo events (and optionally other events)...
then fatal Wrong number of arguments fi set_variable STRING $argv[1] set_variable DETECTORXY_TXT $WORKDIR $DETECTORXY_TXT tail read X Y CHI2 RMS printf optimum n $X $Y $CHI2 $RMS awk v Y
Simple data structure for histogram binning.
string outputFile
bool is_noise(const Hit &hit)
Verify hit origin.
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1989
Definition of hit and track types and auxiliary methods for handling Monte Carlo data.
#define ERROR(A)
Definition: JMessage.hh:66
bool is_initialstate(const Trk &track)
Test whether given track corresponds to an initial state particle.
General purpose messaging.
#define FATAL(A)
Definition: JMessage.hh:67
General purpose class for object reading from a list of file names.
Utility class to parse command line options.
Primary particle.
Definition: JHead.hh:1161
no fit printf nominal n $STRING awk v X
int j
Definition: JPolint.hh:703
General purpose class for multiple pointers.
then fatal Wrong number of arguments fi set_variable DETECTOR $argv[1] set_variable INPUT_FILE $argv[2] eval JPrintDetector a $DETECTOR O IDENTIFIER eval JPrintDetector a $DETECTOR O SUMMARY JAcoustics sh $DETECTOR_ID source JAcousticsToolkit sh CHECK_EXIT_CODE typeset A EMITTERS get_tripods $WORKDIR tripod txt EMITTERS get_transmitters $WORKDIR transmitter txt EMITTERS for EMITTER in
Definition: JCanberra.sh:46
int debug
debug level
The Evt class respresent a Monte Carlo (MC) event as well as an offline event.
Definition: Evt.hh:20
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:62