Jpp  debug
the software that should make you happy
JPlotPMTParameters2D.cc
Go to the documentation of this file.
1 #include <string>
2 #include <iostream>
3 #include <iomanip>
4 #include <vector>
5 #include <map>
6 
7 #include "TROOT.h"
8 #include "TFile.h"
9 #include "TH2D.h"
10 
12 #include "JDetector/JDetector.hh"
18 #include "JROOT/JManager.hh"
19 #include "JGizmo/JGizmoToolkit.hh"
20 
21 #include "Jeep/JPrint.hh"
22 #include "Jeep/JParser.hh"
23 #include "Jeep/JMessage.hh"
24 
25 namespace {
26 
28 
29  /**
30  * Get efficiency.
31  *
32  * \param parameters PMT parameters
33  * \return efficiency
34  */
35  inline double getEfficiency(const JPMTParameters& parameters)
36  {
37  return getSurvivalProbability(parameters) * parameters.QE;
38  }
39 
40  static const char* const EFFICIENCY = "EFFICIENCY";
41 }
42 
43 
44 /**
45  * \file
46  * Auxiliary application to plot PMT parameters.
47  *
48  * \author mdejong, rgruiz
49  */
50 
51 int main(int argc, char **argv)
52 {
53  using namespace std;
54  using namespace JPP;
55  using namespace KM3NETDAQ;
56 
57  string detectorFile;
58  vector<string> inputFile;
59  string outputFile;
60  string regexp;
61  int labelInterval;
62  bool showPMTAddress;
63  bool relative;
64  int debug;
65 
66  try {
67 
68  JParser<> zap("Auxiliary application to plot PMT parameters.");
69 
70  zap['a'] = make_field(detectorFile, "detector file.");
71  zap['P'] = make_field(inputFile, "PMT calibration data file(s)");
72  zap['o'] = make_field(outputFile, "output file.") = "pmt_parameters.root";
73  zap['r'] = make_field(regexp, "regular expresion to extract bin labels for the x-axis") = "";
74  zap['A'] = make_field(showPMTAddress, "show PMT address on y-axis");
75  zap['L'] = make_field(labelInterval, "interval between x-axis bins for which labels are shown") = 1;
76  zap['R'] = make_field(relative, "monitor changes relative to first input");
77  zap['d'] = make_field(debug, "debug") = 0;
78 
79  zap(argc, argv);
80  }
81  catch(const exception &error) {
82  FATAL(error.what() << endl);
83  }
84 
86 
87  try {
88  load(detectorFile, detector);
89  }
90  catch(const JException& error) {
91  FATAL(error);
92  }
93 
94  if (detector.empty()) {
95  FATAL("Empty detector." << endl);
96  }
97 
98  vector<JPMTParametersMap> parameters;
99 
100  for (vector<string>::const_iterator i = inputFile.begin(); i != inputFile.end(); ++i) {
101  parameters.push_back(JPMTParametersMap(i->c_str()));
102  }
103 
104 
105  const int NUMBER_OF_FILES = parameters.size();
106 
107  JManager<string, TH2D> manager(new TH2D("%", NULL,
108  NUMBER_OF_FILES, -0.5, NUMBER_OF_FILES - 0.5,
109  NUMBER_OF_PMTS, -0.5, NUMBER_OF_PMTS - 0.5));
110 
111  manager->Sumw2(kFALSE);
112 
113  if (regexp != "") {
114 
115  const int n = (NUMBER_OF_FILES < labelInterval) ? 1 : labelInterval;
116 
117  const TPRegexp buffer(regexp);
118 
119  for (int i = 0; i != NUMBER_OF_FILES; ++i) {
120 
121  if(i%n == 0)
122  manager->GetXaxis()->SetBinLabel(i+1 , parse(buffer, TString(inputFile[i].c_str())));
123  else
124  manager->GetXaxis()->SetBinLabel(i+1 , " ");
125  }
126  }
127 
128  for (JDetector::iterator module = detector.begin(); module != detector.end(); ++module) {
129 
130  DEBUG("Module " << setw(10) << module->getID() << endl);
131 
132  for (int pmt = 0; pmt != NUMBER_OF_PMTS; ++pmt) {
133 
134  map<string, double> buffer; // buffer values of first input
135 
136  for (int i = 0; i != NUMBER_OF_FILES; ++i) {
137 
138  const JProperties properties = parameters[i].getPMTParameters(JPMTIdentifier(module->getID(), pmt)).getProperties();
139 
140  for (JProperties::const_iterator p = properties.begin(); p != properties.end(); ++p) {
141 
142  double value = 0.0;
143 
144  try { value = (p->second.getValue<const double>()); } catch(const exception& error) {}
145  try { value = (p->second.getValue<const bool>() ? 1.0 : 0.0); } catch(const exception& error) {}
146 
147  if (i == 0) { buffer[p->first] = value; }
148  if (relative) { value -= buffer[p->first]; }
149 
150  manager[MAKE_CSTRING(module->getID() << "." << p->first)]->SetBinContent(i + 1, pmt + 1, value);
151  }
152 
153  double value = getEfficiency(parameters[i].getPMTParameters(JPMTIdentifier(module->getID(), pmt)));
154 
155  if (i == 0) { buffer[EFFICIENCY] = value; }
156  if (relative) { value -= buffer[EFFICIENCY]; }
157 
158  manager[MAKE_CSTRING(module->getID() << "." << EFFICIENCY)]->SetBinContent(i + 1, pmt + 1, value);
159  }
160  }
161  }
162 
163  for (JManager<string, TH2D>::iterator i = manager.begin(); i != manager.end(); ++i) {
164  i->second->Sumw2(kFALSE);
165  }
166 
167  if (showPMTAddress) {
168 
169  const JDetectorAddressMap& demo = getDetectorAddressMap(detector.getID());
170 
171  for (JManager<string, TH2D>::iterator i = manager.begin(); i != manager.end(); ++i) {
172 
173  int id;
174 
175  istringstream(i->first) >> id;
176 
177  setAxisLabels(*i->second, "Y", demo.get(id));
178  }
179  }
180 
181  manager.Write(outputFile.c_str());
182 }
string outputFile
KM3NeT DAQ constants, bit handling, etc.
Detector support kit.
Data structure for detector geometry and calibration.
Dynamic ROOT object management.
General purpose messaging.
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:62
#define FATAL(A)
Definition: JMessage.hh:67
int debug
debug level
Definition: JSirene.cc:69
Utility class to parse command line options.
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:2158
int main(int argc, char **argv)
I/O formatting auxiliaries.
#define MAKE_CSTRING(A)
Make C-string.
Definition: JPrint.hh:136
Lookup table for PMT addresses in detector.
const JModuleAddressMap & get(const int id) const
Get module address map.
Detector data structure.
Definition: JDetector.hh:96
Auxiliary class for map of PMT parameters.
Data structure for PMT parameters.
Utility class to parse parameter values.
Definition: JProperties.hh:501
General exception.
Definition: JException.hh:24
Utility class to parse command line options.
Definition: JParser.hh:1714
void Write(TDirectory &out, const bool wm=false)
Write objects to file.
Definition: JManager.hh:304
void load(const std::string &file_name, JDetector &detector)
Load detector from input file.
double getSurvivalProbability(const JPMTParameters &parameters)
Get model dependent probability that a one photo-electron hit survives the simulation of the PMT assu...
JDetectorAddressMap & getDetectorAddressMap()
Get detector address map.
void setAxisLabels(TAxis *axis, const JModuleAddressMap &memo)
Set axis with PMT address labels.
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
const int n
Definition: JPolint.hh:786
KM3NeT DAQ data structures and auxiliaries.
Definition: DataQueue.cc:39
static const int NUMBER_OF_PMTS
Total number of PMTs in module.
Definition: JDAQ.hh:26
boost::property_tree::ptree parse(std::string str)
Definition: configure.hh:24
Definition: JSTDTypes.hh:14
Detector file.
Definition: JHead.hh:227