Jpp
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
JMuonPostfit.cc
Go to the documentation of this file.
1 #include <string>
2 #include <iostream>
3 #include <iomanip>
4 #include <vector>
5 #include <set>
6 
7 #include "TROOT.h"
8 #include "TFile.h"
9 #include "TH1D.h"
10 #include "TH2D.h"
11 #include "TProfile.h"
12 #include "TMath.h"
13 
19 #include "JAAnet/JHead.hh"
20 #include "JAAnet/JHeadToolkit.hh"
21 #include "JDAQ/JDAQEventIO.hh"
22 #include "JTools/JConstants.hh"
23 #include "JTools/JQuantile.hh"
27 #include "JSupport/JSupport.hh"
28 #include "JReconstruction/JEvt.hh"
30 #include "JGizmo/JVolume.hh"
31 #include "JLang/JPredicate.hh"
32 
33 #include "Jeep/JParser.hh"
34 #include "Jeep/JMessage.hh"
35 
36 
37 namespace {
38 
40 
41 
42  /**
43  * Count number of unique identifiers in event.
44  *
45  * \param __begin begin of data
46  * \param __end end of data
47  * \return number of unique identifiers
48  */
49  template<class JHit_t>
50  inline int getCount(JDAQEvent::const_iterator<JHit_t> __begin,
52  {
53  using namespace std;
54  using namespace KM3NETDAQ;
55 
56  typedef JDAQModuleIdentifier JType_t;
57 
58  set<JType_t> buffer;
59 
60  for (JDAQEvent::const_iterator<JHit_t> i = __begin; i != __end; ++i) {
61  buffer.insert(JType_t(*i));
62  }
63 
64  return buffer.size();
65  }
66 
67 
68  const char* const muon_t = "muon";
69  const char* const neutrino_t = "neutrino";
70 }
71 
72 /**
73  * \file
74  *
75  * Example program to histogram fit results for the muon reconstruction chain.
76  * \author mdejong, bofearraigh.
77  */
78 int main(int argc, char **argv)
79 {
80  using namespace std;
81  using namespace JPP;
82  using namespace KM3NETDAQ;
83 
84  typedef JTriggeredFileScanner<JEvt> JTriggeredFileScanner_t;
85  typedef JTriggeredFileScanner_t::multi_pointer_type multi_pointer_type;
86 
87  JTriggeredFileScanner_t inputFile;
88  JLimit_t& numberOfEvents = inputFile.getLimit();
89  string outputFile;
90  size_t numberOfPrefits;
92  JAtmosphericMuon atmosphere;
93  double Emin_GeV;
94  double NPE;
95  int application;
96  string option;
97  string primary;
98  int debug;
99 
100  try {
101 
102  JParser<> zap("Example program to histogram fit results.");
103 
104  zap['f'] = make_field(inputFile);
105  zap['o'] = make_field(outputFile) = "postfit.root";
106  zap['n'] = make_field(numberOfEvents) = JLimit::max();
107  zap['N'] = make_field(numberOfPrefits) = 1;
109  zap['E'] = make_field(Emin_GeV) = 0.0;
110  zap['M'] = make_field(NPE) = 0.0;
112  zap['a'] = make_field(atmosphere) = JAtmosphericMuon(90.0, 90.0);
113  zap['O'] = make_field(option) = "E", "N", "LINE", "LOGE";
114  zap['p'] = make_field(primary) = muon_t, neutrino_t;
115  zap['d'] = make_field(debug) = 2;
116 
117  zap(argc, argv);
118  }
119  catch(const exception& error) {
120  FATAL(error.what() << endl);
121  }
122 
123  JHead head;
124 
125  try {
126  head = getHeader(inputFile);
127  }
128  catch(const JException& error) {
129  FATAL(error);
130  }
131 
132  const JVolume volume(head, option != "LINE");
133  const JPosition3D center = get<JPosition3D>(head);
134  JCylinder3D cylinder = get<JCylinder3D>(head);
135 
136  cylinder.add(center);
137 
138  NOTICE("Center [m]: " << center << endl);
139  NOTICE("Reposition can [m]: " << cylinder << endl);
140 
141 
142  TFile out(outputFile.c_str(), "recreate");
143 
144  TH1D job("job", NULL, 100, 0.5, 100.5);
145 
146  TH1D hn("hn", NULL, 250, -0.5, 249.5); // NDF
147  TH1D hq("hq", NULL, 300, 0.0, 600.0); // quality parameter
148  TH1D hi("hi", NULL, 101, -0.5, 100.5); // index of the best event in order of descending angular resolution
149  TH1D hv("hv", NULL, 200, -6.0, 0.0); // log of Poisson distribution of (#hits,number photo-electrons)
150  TH1D h1("h1", NULL, 200, -2.0, +2.0); // number of photo-electrons along the whole track
151  TH1D hc("hc", NULL, 200, -1.0, +1.0); // dz (z-slope)
152  TH1D hu("hu", NULL, 400, -1.0e3, 1.0e3); // quality difference between best up-going and down-going track (check of atmospheric muon)
153 
154  TH1D hx("hx", NULL, 70, -3.0, +2.3); // angular deviation [log(deg)]
155  TH1D hd("hd", NULL, 100, 0.0, 10.0); // distance between true and best event at intersection point [m]
156  TH1D hz("hz", NULL, 100, -200.0, 200.0); // intersection of track with neutrino vertex [m]
157  TH1D ht("ht", NULL, 100, -100.0, 100.0); // time between true and best track [ns]
158 
159  TH1D E_0("E_0", NULL, 100, volume.getXmin(), volume.getXmax()); // true muon energy
160  TH1D E_1("E_1", NULL, 100, volume.getXmin(), volume.getXmax()); // uncorrected energy of best track
161  TH1D E_2("E_2", NULL, 100, volume.getXmin(), volume.getXmax()); // corrected energy of best track
162  TH1D E_E("E_E", NULL, 100, -5.0, +5.0); // ratio of reconstructed energy and true energy [log(E/E)]
163  TH2D ExE("ExE", NULL,
164  40, volume.getXmin(), volume.getXmax(),
165  40, volume.getXmin(), volume.getXmax()); // (Etrue,Ereco)
166 
167  const Int_t ny = 28800;
168  const Double_t ymin = 0.0; // [deg]
169  const Double_t ymax = 180.0; // [deg]
170 
171  vector<double> X;
172 
173  if (option.find('E') != string::npos) {
174 
175  for (double x = volume.getXmin(); x <= volume.getXmax(); x += (volume.getXmax() - volume.getXmin()) / 20) {
176  X.push_back(x);
177  }
178 
179  } else {
180 
181  double x = -0.5;
182 
183  for ( ; x <= 15.5; x += 1.0) { X.push_back(x); }
184  for ( ; x <= 25.5; x += 2.0) { X.push_back(x); }
185  for ( ; x <= 50.5; x += 5.0) { X.push_back(x); }
186  for ( ; x <= 100.5; x += 10.0) { X.push_back(x); }
187  for ( ; x <= 250.5; x += 20.0) { X.push_back(x); }
188  }
189 
190  TH2D h2("h2", NULL, X.size() - 1, X.data(), ny, ymin, ymax);
191 
192  TH2D Va("Va", NULL,
193  100, 0.0, cylinder.getRadius()*cylinder.getRadius(),
194  100, cylinder.getZmin() - 50.0, cylinder.getZmax() + 50.0);
195  TH2D Vb("Vb", NULL,
196  100, 0.0, cylinder.getRadius()*cylinder.getRadius(),
197  100, cylinder.getZmin() - 50.0, cylinder.getZmax() + 50.0);
198 
199  JQuantile Q("Angle", true);
200  JQuantile O("Omega", true);
201 
202 
203  while (inputFile.hasNext()) {
204 
205  STATUS("event: " << setw(10) << inputFile.getCounter() << '\r'); DEBUG(endl);
206 
207  multi_pointer_type ps = inputFile.next();
208 
209  JDAQEvent* tev = ps;
210  JEvt* evt = ps;
211  Evt* event = ps;
212 
213  const JTimeConverter converter(*event, *tev);
214 
215  job.Fill(1.0);
216 
217 
218  double Enu = 0.0;
219  double Emu = 0.0;
220  double Emax = 0.0;
221 
222  if (has_neutrino(*event)) {
223  Enu = get_neutrino(*event).E;
224  }
225 
226  vector<Trk>::const_iterator muon = event->mc_trks.end();
227 
228  for (vector<Trk>::const_iterator track = event->mc_trks.begin(); track != event->mc_trks.end(); ++track) {
229 
230  if (is_muon(*track)) {
231 
232  Emu += track->E;
233 
234  if (track->E > Emax) {
235  muon = track;
236  Emax = track->E;
237  }
238  }
239  }
240 
241  if (muon == event->mc_trks.end()) {
242  continue;
243  }
244 
245  job.Fill(3.0);
246 
247 
248  // abscissa
249 
250  Double_t x = 0.0;
251 
252  if (option.find('E') != string::npos)
253  x = volume.getX(event->mc_trks[0].E);
254  else
255  x = getCount(tev->begin<JDAQTriggeredHit>(), tev->end<JDAQTriggeredHit>());
256 
257 
258  if (!evt->empty()) {
259 
260  JEvt::iterator __end = partition(evt->begin(), evt->end(), JHistory::is_event(application));
261 
262  if (evt->begin() == __end) {
263  continue;
264  }
265 
266  job.Fill(4.0);
267 
268  if (numberOfPrefits > 0) {
269 
270  JEvt::iterator __q = __end;
271 
272  advance(__end = evt->begin(), min(numberOfPrefits, (size_t) distance(evt->begin(), __q)));
273 
274  partial_sort(evt->begin(), __end, __q, quality_sorter);
275 
276  } else {
277 
278  sort(evt->begin(), __end, quality_sorter);
279  }
280 
281  double E = 0.0;
282  JPointing pointing;
283 
284  if (primary == muon_t) {
285  E = Emu;
286  pointing = JPointing(getDirection(*muon));
287  } else if (primary == neutrino_t) {
288  E = Enu;
289  pointing = JPointing(getDirection(get_neutrino(*event)));
290  }
291 
292  JEvt::iterator best = pointing(evt->begin(), __end);
293  const Double_t beta = pointing.getAngle(*best);
294  const double Efit = best->getE();
295  const double Eraw = best->getW(JENERGY_ENERGY, numeric_limits<double>::min());
296  const double mip = best->getW(JSTART_NPE_MIP, numeric_limits<double>::max());
297 
298  // selection of fit result
299 
300  bool ok = (Efit >= Emin_GeV && mip >= NPE);
301 
302  if (ok) {
303 
304  job.Fill(5.0);
305 
306  hn.Fill((Double_t) best->getNDF());
307  hq.Fill(best->getQ());
308  hi.Fill((Double_t) distance(evt->begin(), best));
309  hc.Fill(best->getDZ());
310 
311  // atmospheric muon hypothesis test.
312 
313  if (( has_neutrino(*event) && get_neutrino(*event).dir.z >= atmosphere.dot2) || // difference in quality between best upwards and best downgoing track.
314  (!has_neutrino(*event) && best->getDZ() >= atmosphere.dot2)) { // negative (positive) values imply a down(up) going track.
315  hu.Fill(atmosphere(evt->begin(), __end));
316  }
317 
318  hx.Fill(max(log10(beta), hx.GetXaxis()->GetXmin()));
319 
320  Q.put(beta);
321 
322  JTrack3E ta = getTrack(*muon);
323  JTrack3E tb = getTrack(*best);
324 
325  ta.add(center);
326  tb.sub(converter.putTime());
327 
328  static_cast<JTrack3D&>(ta).move(ta.getIntersection(tb), getSpeedOfLight()); // move to intersection point of both tracks
329  static_cast<JTrack3D&>(tb).move(tb.getIntersection(ta), getSpeedOfLight());
330 
331  hd.Fill((tb.getPosition() - ta.getPosition()).getLength());
332 
333  if (has_neutrino(*event)) {
334 
335  job.Fill(6.0);
336 
337  const Trk neutrino = get_neutrino(*event);
338 
339  JPosition3D vertex = getPosition(neutrino);
340 
341  vertex.add(center);
342 
343  if (cylinder.is_inside(vertex)) {
344 
345  job.Fill(7.0);
346 
347  JTrack3E tc = getTrack(*best);
348 
349  Va.Fill(best->getX()*best->getX() + best->getY()*best->getY(), best->getZ()); //R^{2} [m^2], z [m]
350  hz.Fill(tc.getIntersection(vertex));
351  }
352  }
353 
354  Vb.Fill(best->getX()*best->getX() + best->getY()*best->getY(), best->getZ()); //R^{2} [m^2], z [m]
355  ht.Fill(tb.getT() - ta.getT());
356 
357  E_0.Fill(volume.getX(E, true));
358  E_1.Fill(volume.getX(Eraw, true));
359  E_2.Fill(volume.getX(Efit, true));
360  E_E.Fill(volume.getX(Efit) - volume.getX(E));
361  ExE.Fill(volume.getX(E), volume.getX(Efit));
362 
363  h2.Fill(x, beta);
364 
365  if (best->hasW(JSTART_NPE_MIP)) {
366  h1.Fill(log10(best->getW(JSTART_NPE_MIP)));
367  }
368 
369  if (best->hasW(JVETO_NPE) && best->hasW(JVETO_NUMBER_OF_HITS)) {
370 
371  const double npe = best->getW(JVETO_NPE); // number of photo-electrons
372  const int count = best->getW(JVETO_NUMBER_OF_HITS); // number of hits
373  const double pv = TMath::PoissonI(count, npe); // Poisson distribution function for (#hits,npe)
374 
375  hv.Fill(max(log10(pv), hv.GetXaxis()->GetXmin()));
376  }
377  }
378  }
379  }
380  STATUS(endl);
381 
382  NOTICE("Number of events input " << setw(8) << right << job.GetBinContent(1) << endl);
383  NOTICE("Number of events with muon " << setw(8) << right << job.GetBinContent(3) << endl);
384  NOTICE("Number of events with fit " << setw(8) << right << job.GetBinContent(4) << endl);
385  NOTICE("Number of events selected " << setw(8) << right << job.GetBinContent(5) << endl);
386  NOTICE("Number of events with neutrino " << setw(8) << right << job.GetBinContent(6) << endl);
387  NOTICE("Number of events contained " << setw(8) << right << job.GetBinContent(7) << endl);
388 
389  if (Q.getCount() != 0) {
390  NOTICE("Median space angle [deg] " << FIXED (6,3) << Q.getQuantile(0.5) << endl);
391  }
392 
393  out.Write();
394  out.Close();
395 }
static const int JMUONSTART
Utility class to parse command line options.
Definition: JParser.hh:1493
General exception.
Definition: JException.hh:23
ROOT TTree parameter settings.
static const int JMUONPATH
static const int JVETO_NPE
number of photo-electrons from JVeto.cc
Auxiliary class to compare fit results with respect to a reference direction (e.g.
Definition: JEvtToolkit.hh:586
JTrack3E getTrack(const Trk &track)
Get track.
static const int JENERGY_ENERGY
uncorrected energy [GeV] from JEnergy.cc
double z
Definition: Vec.hh:14
Synchronously read DAQ events and Monte Carlo events (and optionally other events).
double putTime() const
Get Monte Carlo minus DAQ/trigger hit time.
std::vector< T >::difference_type distance(typename std::vector< T >::const_iterator first, typename PhysicsEvent::const_iterator< T > second)
Specialisation of STL distance.
JTime & sub(const JTime &value)
Subtraction operator.
const char *const neutrino_t
Definition: io_ascii.hh:28
Quantile calculator.
Definition: JQuantile.hh:83
const double getSpeedOfLight()
Number of bytes in a gigabyte.
Definition: JConstants.hh:89
#define STATUS(A)
Definition: JMessage.hh:63
bool has_neutrino(const Evt &evt)
Test whether given event has an incoming neutrino.
Template const_iterator.
Definition: JDAQEvent.hh:68
double getIntersection(const JVector3D &pos) const
Get longitudinal position along axis of position of closest approach with given position.
Definition: JAxis3D.hh:157
bool is_muon(const Trk &track)
Test whether given track is a (anti-)muon.
Vec dir
track direction
Definition: Trk.hh:16
General purpose sorter of fit results.
Definition: JEvtToolkit.hh:262
Auxiliary class to synchronously read DAQ events and Monte Carlo events (and optionally other events)...
then for HISTOGRAM in h0 h1
Definition: JMatrixNZ.sh:69
Empty structure for specification of parser element that is initialised (i.e. does not require input)...
Definition: JParser.hh:63
Double_t getX(const Double_t E, double constrain=false) const
Get abscissa value.
Definition: JVolume.hh:115
JCylinder3D & add(const JVector3D &pos)
Add position.
Definition: JCylinder3D.hh:146
Auxiliary data structure for floating point format specification.
Definition: JPrint.hh:481
string outputFile
double E
Energy [GeV] (either MC truth or reconstructed)
Definition: Trk.hh:18
3D track with energy.
Definition: JTrack3E.hh:29
double getQuantile(const double Q, const bool reverse=false) const
Get quantile.
Definition: JQuantile.hh:402
Double_t getXmin() const
Get minimal abscissa value.
Definition: JVolume.hh:91
static const int JSTART_NPE_MIP
number of photo-electrons up to the barycentre from JStart.cc
JTime & add(const JTime &value)
Addition operator.
Head getHeader(const JMultipleFileScanner_t &file_list)
Get Monte Carlo header.
Auxiliary class for defining the range of iterations of objects.
Definition: JLimit.hh:41
const_iterator< T > end() const
Get end of data.
Constants.
static const int JMUONPREFIT
Cylinder object.
Definition: JCylinder3D.hh:37
double getAngle(const JFit &fit) const
Get angle between reference direction and fit result.
Definition: JEvtToolkit.hh:635
const_iterator< T > begin() const
Get begin of data.
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1954
Auxiliary class to test history.
Definition: JHistory.hh:104
static const int JMUONGANDALF
#define NOTICE(A)
Definition: JMessage.hh:64
void put(const double x, const double w=1.0)
Put value.
Definition: JQuantile.hh:198
static const int JVETO_NUMBER_OF_HITS
number of hits from JVeto.cc
int debug
debug level
Definition: JSirene.cc:61
const JPosition3D & getPosition() const
Get position.
Definition: JPosition3D.hh:129
long long int getCount() const
Get total count.
Definition: JQuantile.hh:276
Reconstruction type dependent comparison of track quality.
Auxiliary class to evaluate atmospheric muon hypothesis.
Definition: JEvtToolkit.hh:772
double getT(const JVector3D &pos) const
Get arrival time of Cherenkov light at given position.
Definition: JTrack3D.hh:126
General purpose messaging.
counter_type advance(counter_type &counter, const counter_type value, const counter_type limit=std::numeric_limits< counter_type >::max())
Advance counter.
Monte Carlo run header.
Definition: JHead.hh:836
#define FATAL(A)
Definition: JMessage.hh:67
Data structure for set of track fit results.
Definition: JEvt.hh:294
std::vector< int > count
Definition: JAlgorithm.hh:184
static const int JMUONSIMPLEX
Utility class to parse command line options.
int getCount(const T &hit)
Get hit count.
Primary particle.
Definition: JHead.hh:776
JDirection3D getDirection(const Vec &v)
Get direction.
Auxiliary class to convert DAQ/trigger hit time to/from Monte Carlo hit time.
Data structure for position in three dimensions.
Definition: JPosition3D.hh:35
const JLimit & getLimit() const
Get limit.
Definition: JLimit.hh:73
const Trk & get_neutrino(const Evt &evt)
Get incoming neutrino.
The Trk class represents a Monte Carlo (MC) particle as well as a reconstructed track/shower.
Definition: Trk.hh:12
Double_t getXmax() const
Get maximal abscissa value.
Definition: JVolume.hh:102
Auxiliary class for histogramming of effective volume.
Definition: JVolume.hh:26
JVector3D & add(const JVector3D &vector)
Add vector.
Definition: JVector3D.hh:141
static const int JMUONENERGY
The Evt class respresent a Monte Carlo (MC) event as well as an offline event.
Definition: Evt.hh:19
then usage $script[input file[working directory[option]]] nWhere option can be E
Definition: JMuonPostfit.sh:37
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:62
JPosition3D getPosition(const Vec &v)
Get position.
int main(int argc, char *argv[])
Definition: Main.cpp:15