Jpp  17.3.0-rc.1
the software that should make you happy
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
JParramatta.cc
Go to the documentation of this file.
1 #include <iostream>
2 #include <iomanip>
3 
4 #include "TROOT.h"
5 #include "TFile.h"
6 #include "TH1D.h"
7 #include "TH2D.h"
8 #include "TProfile.h"
9 #include "TGraph.h"
10 
12 #include "JSupport/JTreeScanner.hh"
13 
14 #include "JROOT/JRootToolkit.hh"
15 #include "JROOT/JGraph.hh"
16 #include "JROOT/JManager.hh"
17 
18 #include "JAcoustics/JEvt.hh"
19 #include "JAcoustics/JSupport.hh"
20 
21 #include "Jeep/JPrint.hh"
22 #include "Jeep/JParser.hh"
23 #include "Jeep/JMessage.hh"
24 
25 
26 /**
27  * \file
28  *
29  * Example program to plot acoustic fit results.
30  * \author mdejong
31  */
32 int main(int argc, char **argv)
33 {
34  using namespace std;
35  using namespace JPP;
36 
38  JLimit_t& numberOfEvents = inputFile.getLimit();
39  string outputFile;
40  int debug;
41 
42  try {
43 
44  JParser<> zap("Example program to plot acoustic fit results.");
45 
46  zap['f'] = make_field(inputFile, "input file (output of JKatoomba[.sh])");
47  zap['n'] = make_field(numberOfEvents) = JLimit::max();
48  zap['o'] = make_field(outputFile) = "parramatta.root";
49  zap['d'] = make_field(debug) = 2;
50 
51  zap(argc, argv);
52  }
53  catch(const exception &error) {
54  FATAL(error.what() << endl);
55  }
56 
57  const JFormat_t format(4, 0, std::ios_base::fmtflags(), '0');
58 
59  TH1D h1("chi2/NDF", NULL, 500, 0.0, 10.0);
60 
61  JGraph_t g0;
62  JGraph_t g1;
63 
67 
68  JManager<int, TProfile> H1(new TProfile("stretching[%]", NULL, 50, 0.0, 1.0e2));
69  JManager<int, TH2D> H2(new TH2D ("[%].tiltdeviation", NULL, 300, -4.0, +4.0, 300, -4.0, +4.0), '%', format);
70 
71  while (inputFile.hasNext()) {
72 
73  STATUS("event: " << setw(10) << inputFile.getCounter() << '\r'); DEBUG(endl);
74 
75  const JEvt* evt = inputFile.next();
76 
77  h1.Fill(evt->chi2 / evt->ndf);
78 
79  const double t1 = 0.5 * (evt->UNIXTimeStart + evt->UNIXTimeStop);
80 
81  g0.put(t1, evt->nhit - evt->npar);
82  g1.put(t1, evt->chi2 / evt->ndf);
83 
84  double sum_tx = 0;
85  double sum_ty = 0;
86 
87  for (JEvt::const_iterator i = evt->begin(); i != evt->end(); ++i) {
88 
89  const int id = i->id;
90  const double tx = i->tx * 1.0e3; // [mrad]
91  const double ty = i->ty * 1.0e3; // [mrad]
92  const double vs = i->vs * 1.0e2; // [%]
93  const double ts = sqrt(tx*tx + ty*ty);
94 
95  sum_tx += tx; // [mrad]
96  sum_ty += ty; // [mrad]
97 
98  H1 ->Fill(ts, vs);
99  H1[id]->Fill(ts, vs);
100 
101  GO[id].put(t1, atan2(ty, tx));
102  GA[id].put(t1, ts);
103  GS[id].put(t1, vs);
104 
105  }
106 
107  sum_tx /= evt->size();
108  sum_ty /= evt->size();
109 
110  for (JEvt::const_iterator i = evt->begin(); i != evt->end(); ++i) {
111 
112  H2[i->id]->Fill(i->tx * 1.0e3 - sum_tx, i->ty * 1.0e3 - sum_ty); // [mrad]
113  }
114  }
115  STATUS(endl);
116 
117  TFile out(outputFile.c_str(), "recreate");
118 
119  out << h1
120  << JGraph(g0, "g0")
121  << JGraph(g1, "g1");
122 
123  out << H1 << *H1 << H2;
124 
125  for (map<int, JGraph_t>::const_iterator i = GO.begin(); i != GO.end(); ++i) {
126  out << JGraph(i->second, MAKE_CSTRING("G[" << FILL(4,'0') << i->first << "].orientation"));
127  }
128 
129  for (map<int, JGraph_t>::const_iterator i = GA.begin(); i != GA.end(); ++i) {
130  out << JGraph(i->second, MAKE_CSTRING("G[" << FILL(4,'0') << i->first << "].amplitude"));
131  }
132 
133  for (map<int, JGraph_t>::const_iterator i = GS.begin(); i != GS.end(); ++i) {
134  out << JGraph(i->second, MAKE_CSTRING("G[" << FILL(4,'0') << i->first << "].stretching"));
135  }
136 
137  out.Write();
138  out.Close();
139 }
Utility class to parse command line options.
Definition: JParser.hh:1517
int npar
number of fit parameters
int main(int argc, char *argv[])
Definition: Main.cc:15
Data structure for graph data.
Definition: JGraph.hh:21
#define STATUS(A)
Definition: JMessage.hh:63
ROOT TTree parameter settings.
#define MAKE_CSTRING(A)
Make C-string.
Definition: JPrint.hh:136
Dynamic ROOT object management.
Auxiliary data structure to build TGraph.
Definition: JGraph.hh:42
string outputFile
double UNIXTimeStop
stop time
Auxiliary class for defining the range of iterations of objects.
Definition: JLimit.hh:41
Auxiliary class to manage set of compatible ROOT objects (e.g. histograms) using unique keys...
Definition: JManager.hh:43
I/O formatting auxiliaries.
std::vector< double > vs
Acoustic event fit.
double ndf
weighed number of degrees of freedom
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1993
General purpose messaging.
Auxiliary data structure for sequence of same character.
Definition: JManip.hh:328
#define FATAL(A)
Definition: JMessage.hh:67
Scanning of objects from multiple files according a format that follows from the extension of each fi...
int nhit
number of hits
double UNIXTimeStart
start time
General purpose class for object reading from a list of file names.
Utility class to parse command line options.
const JLimit & getLimit() const
Get limit.
Definition: JLimit.hh:73
Acoustic event fit.
Data structure for format specifications.
Definition: JManip.hh:522
int debug
debug level
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:62
Double_t g1(const Double_t x)
Function.
Definition: JQuantiles.cc:25