Jpp  19.1.0-rc.1
the software that should make you happy
Functions
JParramatta.cc File Reference

Example program to plot acoustic fit results. More...

#include <iostream>
#include <iomanip>
#include <map>
#include "TROOT.h"
#include "TFile.h"
#include "TH1D.h"
#include "TH2D.h"
#include "TProfile.h"
#include "TGraph.h"
#include "JSupport/JMultipleFileScanner.hh"
#include "JSupport/JTreeScanner.hh"
#include "JROOT/JRootToolkit.hh"
#include "JROOT/JGraph.hh"
#include "JROOT/JManager.hh"
#include "JTools/JQuantile.hh"
#include "JAcoustics/JEvt.hh"
#include "JAcoustics/JSuperEvt.hh"
#include "JAcoustics/JSupport.hh"
#include "JLang/JObjectMultiplexer.hh"
#include "Jeep/JPrint.hh"
#include "Jeep/JParser.hh"
#include "Jeep/JMessage.hh"

Go to the source code of this file.

Functions

int main (int argc, char **argv)
 

Detailed Description

Example program to plot acoustic fit results.

Author
mdejong

Definition in file JParramatta.cc.

Function Documentation

◆ main()

int main ( int  argc,
char **  argv 
)

Definition at line 38 of file JParramatta.cc.

39 {
40  using namespace std;
41  using namespace JPP;
42 
43  typedef JTYPELIST<JEvt, JSuperEvt>::typelist typelist;
44 
46  JLimit_t& numberOfEvents = inputFile.getLimit();
47  string outputFile;
48  double Z;
49  int debug;
50 
51  try {
52 
53  JParser<> zap("Example program to plot acoustic fit results.");
54 
55  zap['f'] = make_field(inputFile, "input file (output of JKatoomba[.sh]/JFremantle[.sh])");
56  zap['n'] = make_field(numberOfEvents) = JLimit::max();
57  zap['o'] = make_field(outputFile) = "parramatta.root";
58  zap['Z'] = make_field(Z, "detector height (for 2nd order tilt correction)") = 0.0;
59  zap['d'] = make_field(debug) = 2;
60 
61  zap(argc, argv);
62  }
63  catch(const exception &error) {
64  FATAL(error.what() << endl);
65  }
66 
67  Z *= 0.5; // half height
68 
69  const JFormat_t format(4, 0, std::ios_base::fmtflags(), '0');
70 
71  TH1D h1("chi2/NDF", NULL, 500, 0.0, 10.0);
72  TH1D h2("amplitude", NULL, 500, 0.0, 100.0);
73 
74  JGraph_t g0;
75  JGraph_t g1;
76 
80 
81  JManager<int, TProfile> H1(new TProfile("stretching[%]", NULL, 50, 0.0, 1.0e2));
82  JManager<int, TH2D> H2(new TH2D("[%].tiltdeviation", NULL, 300, -4.0, +4.0, 300, -4.0, +4.0), '%', format);
83 
85 
86  for (counter_type counter = 0; in.hasNext() && counter != inputFile.getLimit(); ++counter) {
87 
88  STATUS("event: " << setw(10) << counter << '\r'); DEBUG(endl);
89 
90  const JEvt* evt = in.next();
91 
92  h1.Fill(evt->chi2 / evt->ndf);
93 
94  const double t1 = 0.5 * (evt->UNIXTimeStart + evt->UNIXTimeStop);
95 
96  g0.put(t1, evt->nhit - evt->npar);
97  g1.put(t1, evt->chi2 / evt->ndf);
98 
99  double Tx = 0.0;
100  double Ty = 0.0;
101  double Ts = 0.0;
102  double Vs = 0.0;
103 
104  for (JEvt::const_iterator i = evt->begin(); i != evt->end(); ++i) {
105 
106  const int id = i->id;
107  const double tx = (i->tx + i->tx2 * Z) * 1.0e3; // [mrad]
108  const double ty = (i->ty + i->ty2 * Z) * 1.0e3; // [mrad]
109  const double vs = i->vs * 1.0e2; // [%]
110  const double ts = sqrt(tx*tx + ty*ty);
111 
112  h2.Fill(ts);
113 
114  Tx += tx;
115  Ty += ty;
116  Vs += vs;
117  Ts += ts;
118 
119  H1 ->Fill(ts, vs);
120  H1[id]->Fill(ts, vs);
121 
122  GO[id].put(t1, atan2(ty, tx));
123  GA[id].put(t1, ts);
124  GS[id].put(t1, vs);
125  }
126 
127  Tx /= evt->size();
128  Ty /= evt->size();
129  Vs /= evt->size();
130  Ts /= evt->size();
131 
132  for (JEvt::const_iterator i = evt->begin(); i != evt->end(); ++i) {
133 
134  const double tx = (i->tx + i->tx2 * Z) * 1.0e3; // [mrad]
135  const double ty = (i->ty + i->ty2 * Z) * 1.0e3; // [mrad]
136 
137  H2 ->Fill(tx - Tx, ty - Ty);
138  H2[i->id]->Fill(tx - Tx, ty - Ty);
139 
140  GO[-1].put(t1, atan2(Ty, Tx));
141  GA[-1].put(t1, Ts);
142  GS[-1].put(t1, Vs);
143  }
144  }
145  STATUS(endl);
146 
147  TH1D hx("hx", NULL, H2.size(), -0.5, H2.size() + 0.5);
148  TH1D hy("hy", NULL, H2.size(), -0.5, H2.size() + 0.5);
149 
150  JQuantile Qx, Qy;
151 
152  for (JManager<int, TH2D>::const_iterator i = H2.begin(); i != H2.end(); ++i) {
153 
154  const int ix = distance(H2.cbegin(), i) + 1;
155 
156  hx.GetXaxis()->SetBinLabel(ix, MAKE_CSTRING(i->first));
157  hy.GetXaxis()->SetBinLabel(ix, MAKE_CSTRING(i->first));
158 
159  hx.SetBinContent(ix, i->second->GetMean(1));
160  hy.SetBinContent(ix, i->second->GetMean(2));
161  hx.SetBinError (ix, i->second->GetStdDev(1));
162  hy.SetBinError (ix, i->second->GetStdDev(2));
163 
164  Qx.put(i->second->GetMean(1));
165  Qy.put(i->second->GetMean(2));
166  }
167 
168  if (Qx.getCount() > 1 && Qy.getCount() > 1) {
169  cout << "deviation: " << FIXED(7,3) << Qx.getSTDev() << ' ' << FIXED(7,3) << Qy.getSTDev() << endl;
170  }
171 
172  TFile out(outputFile.c_str(), "recreate");
173 
174  out << h1 << h2 << hx << hy
175  << JGraph(g0, "g0")
176  << JGraph(g1, "g1");
177 
178  out << H1 << *H1 << H2 << *H2;
179 
180  for (map<int, JGraph_t>::const_iterator i = GO.begin(); i != GO.end(); ++i) {
181  out << JGraph(i->second, MAKE_CSTRING("G[" << FILL(4,'0') << i->first << "].orientation"));
182  }
183 
184  for (map<int, JGraph_t>::const_iterator i = GA.begin(); i != GA.end(); ++i) {
185  out << JGraph(i->second, MAKE_CSTRING("G[" << FILL(4,'0') << i->first << "].amplitude"));
186  }
187 
188  for (map<int, JGraph_t>::const_iterator i = GS.begin(); i != GS.end(); ++i) {
189  out << JGraph(i->second, MAKE_CSTRING("G[" << FILL(4,'0') << i->first << "].stretching"));
190  }
191 
192  out.Write();
193  out.Close();
194 }
string outputFile
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:62
#define STATUS(A)
Definition: JMessage.hh:63
#define FATAL(A)
Definition: JMessage.hh:67
int debug
debug level
Definition: JSirene.cc:69
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:2158
#define MAKE_CSTRING(A)
Make C-string.
Definition: JPrint.hh:136
Double_t g1(const Double_t x)
Function.
Definition: JQuantiles.cc:25
std::vector< T >::difference_type distance(typename std::vector< T >::const_iterator first, typename PhysicsEvent::const_iterator< T > second)
Specialisation of STL distance.
Auxiliary class for multiplexing object iterators.
Utility class to parse command line options.
Definition: JParser.hh:1714
Auxiliary class to manage set of compatible ROOT objects (e.g. histograms) using unique keys.
Definition: JManager.hh:47
General purpose class for object reading from a list of file names.
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
Long64_t counter_type
Type definition for counter.
std::vector< double > vs
Definition: JSTDTypes.hh:14
Auxiliary data structure for sequence of same character.
Definition: JManip.hh:330
Auxiliary data structure for floating point format specification.
Definition: JManip.hh:448
Acoustic event fit.
int nhit
number of hits
double UNIXTimeStop
stop time
double ndf
weighed number of degrees of freedom
int npar
number of fit parameters
double UNIXTimeStart
start time
Data structure for format specifications.
Definition: JManip.hh:524
Type list.
Definition: JTypeList.hh:23
Data structure for graph data.
Definition: JGraph.hh:21
void put(const Double_t x, const Double_t y)
Put data.
Definition: JGraph.hh:28
Auxiliary data structure to build TGraph.
Definition: JGraph.hh:44
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:46
double getSTDev() const
Get standard deviation.
Definition: JQuantile.hh:281
void put(const double x, const double w=1.0)
Put value.
Definition: JQuantile.hh:133
long long int getCount() const
Get total count.
Definition: JQuantile.hh:186