Jpp 20.0.0-rc.9-29-gccc23c492-D
the software that should make you happy
Loading...
Searching...
No Matches
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/JStats.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
45 JMultipleFileScanner<typelist> inputFile;
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 JGraph_t g2;
77
78 map<int, JGraph_t> GO;
79 map<int, JGraph_t> GA;
80 map<int, JGraph_t> GS;
81
82 JManager<int, TProfile> H1(new TProfile("stretching[%]", NULL, 50, 0.0, 1.0e2));
83 JManager<int, TH2D> H2(new TH2D("[%].tiltdeviation", NULL, 300, -4.0, +4.0, 300, -4.0, +4.0), '%', format);
84
85 JObjectMultiplexer<typelist, JEvt> in(inputFile);
86
87 for (counter_type counter = 0; in.hasNext() && counter != inputFile.getLimit(); ++counter) {
88
89 STATUS("event: " << setw(10) << counter << '\r'); DEBUG(endl);
90
91 const JEvt* evt = in.next();
92
93 h1.Fill(evt->chi2 / evt->ndf);
94
95 const double t1 = 0.5 * (evt->UNIXTimeStart + evt->UNIXTimeStop);
96
97 g0.put(t1, evt->nhit - evt->npar);
98 g1.put(t1, evt->chi2 / evt->ndf);
99 g2.put(t1, evt->numberOfIterations);
100
101 double Tx = 0.0;
102 double Ty = 0.0;
103 double Ts = 0.0;
104 double Vs = 0.0;
105
106 for (JEvt::const_iterator i = evt->begin(); i != evt->end(); ++i) {
107
108 const int id = i->id;
109 const double tx = (i->tx + i->tx2 * Z) * 1.0e3; // [mrad]
110 const double ty = (i->ty + i->ty2 * Z) * 1.0e3; // [mrad]
111 const double vs = i->vs * 1.0e2; // [%]
112 const double ts = sqrt(tx*tx + ty*ty);
113
114 h2.Fill(ts);
115
116 Tx += tx;
117 Ty += ty;
118 Vs += vs;
119 Ts += ts;
120
121 H1 ->Fill(ts, vs);
122 H1[id]->Fill(ts, vs);
123
124 GO[id].put(t1, atan2(ty, tx));
125 GA[id].put(t1, ts);
126 GS[id].put(t1, vs);
127 }
128
129 Tx /= evt->size();
130 Ty /= evt->size();
131 Vs /= evt->size();
132 Ts /= evt->size();
133
134 for (JEvt::const_iterator i = evt->begin(); i != evt->end(); ++i) {
135
136 const double tx = (i->tx + i->tx2 * Z) * 1.0e3; // [mrad]
137 const double ty = (i->ty + i->ty2 * Z) * 1.0e3; // [mrad]
138
139 H2 ->Fill(tx - Tx, ty - Ty);
140 H2[i->id]->Fill(tx - Tx, ty - Ty);
141
142 GO[-1].put(t1, atan2(Ty, Tx));
143 GA[-1].put(t1, Ts);
144 GS[-1].put(t1, Vs);
145 }
146 }
147 STATUS(endl);
148
149 TH1D hx("hx", NULL, H2.size(), -0.5, H2.size() + 0.5);
150 TH1D hy("hy", NULL, H2.size(), -0.5, H2.size() + 0.5);
151
152 JStats Qx, Qy;
153
154 for (JManager<int, TH2D>::const_iterator i = H2.begin(); i != H2.end(); ++i) {
155
156 const int ix = distance(H2.cbegin(), i) + 1;
157
158 hx.GetXaxis()->SetBinLabel(ix, MAKE_CSTRING(i->first));
159 hy.GetXaxis()->SetBinLabel(ix, MAKE_CSTRING(i->first));
160
161 hx.SetBinContent(ix, i->second->GetMean(1));
162 hy.SetBinContent(ix, i->second->GetMean(2));
163 hx.SetBinError (ix, i->second->GetStdDev(1));
164 hy.SetBinError (ix, i->second->GetStdDev(2));
165
166 Qx.put(i->second->GetMean(1));
167 Qy.put(i->second->GetMean(2));
168 }
169
170 if (Qx.getCount() > 1 && Qy.getCount() > 1) {
171 cout << "deviation: " << FIXED(7,3) << Qx.getSTDev() << ' ' << FIXED(7,3) << Qy.getSTDev() << endl;
172 }
173
174 TFile out(outputFile.c_str(), "recreate");
175
176 out << h1 << h2 << hx << hy
177 << JGraph(g0, "g0")
178 << JGraph(g1, "g1")
179 << JGraph(g2, "g2");
180
181 out << H1 << *H1 << H2 << *H2;
182
183 for (map<int, JGraph_t>::const_iterator i = GO.begin(); i != GO.end(); ++i) {
184 out << JGraph(i->second, MAKE_CSTRING("G[" << FILL(4,'0') << i->first << "].orientation"));
185 }
186
187 for (map<int, JGraph_t>::const_iterator i = GA.begin(); i != GA.end(); ++i) {
188 out << JGraph(i->second, MAKE_CSTRING("G[" << FILL(4,'0') << i->first << "].amplitude"));
189 }
190
191 for (map<int, JGraph_t>::const_iterator i = GS.begin(); i != GS.end(); ++i) {
192 out << JGraph(i->second, MAKE_CSTRING("G[" << FILL(4,'0') << i->first << "].stretching"));
193 }
194
195 out.Write();
196 out.Close();
197}
Acoustic event fit.
Definition JEvt.hh:310
int nhit
number of hits
Definition JEvt.hh:263
double UNIXTimeStop
stop time
Definition JEvt.hh:262
double ndf
weighed number of degrees of freedom
Definition JEvt.hh:266
double chi2
chi2
Definition JEvt.hh:267
int npar
number of fit parameters
Definition JEvt.hh:265
double UNIXTimeStart
start time
Definition JEvt.hh:261
int numberOfIterations
number of iterations
Definition JEvt.hh:268