#include #include #include #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" /** * \file * * Example program to plot acoustic fit results. * \author mdejong */ int main(int argc, char **argv) { using namespace std; using namespace JPP; typedef JTYPELIST::typelist typelist; JMultipleFileScanner inputFile; JLimit_t& numberOfEvents = inputFile.getLimit(); string outputFile; double Z; int debug; try { JParser<> zap("Example program to plot acoustic fit results."); zap['f'] = make_field(inputFile, "input file (output of JKatoomba[.sh]/JFremantle[.sh])"); zap['n'] = make_field(numberOfEvents) = JLimit::max(); zap['o'] = make_field(outputFile) = "parramatta.root"; zap['Z'] = make_field(Z, "detector height (for 2nd order tilt correction)") = 0.0; zap['d'] = make_field(debug) = 2; zap(argc, argv); } catch(const exception &error) { FATAL(error.what() << endl); } Z *= 0.5; // half height const JFormat_t format(4, 0, std::ios_base::fmtflags(), '0'); TH1D h1("chi2/NDF", NULL, 500, 0.0, 10.0); TH1D h2("amplitude", NULL, 500, 0.0, 100.0); JGraph_t g0; JGraph_t g1; map GO; map GA; map GS; JManager H1(new TProfile("stretching[%]", NULL, 50, 0.0, 1.0e2)); JManager H2(new TH2D("[%].tiltdeviation", NULL, 300, -4.0, +4.0, 300, -4.0, +4.0), '%', format); JObjectMultiplexer in(inputFile); for (counter_type counter = 0; in.hasNext() && counter != inputFile.getLimit(); ++counter) { STATUS("event: " << setw(10) << counter << '\r'); DEBUG(endl); const JEvt* evt = in.next(); h1.Fill(evt->chi2 / evt->ndf); const double t1 = 0.5 * (evt->UNIXTimeStart + evt->UNIXTimeStop); g0.put(t1, evt->nhit - evt->npar); g1.put(t1, evt->chi2 / evt->ndf); double Tx = 0.0; double Ty = 0.0; double Ts = 0.0; double Vs = 0.0; for (JEvt::const_iterator i = evt->begin(); i != evt->end(); ++i) { const int id = i->id; const double tx = (i->tx + i->tx2 * Z) * 1.0e3; // [mrad] const double ty = (i->ty + i->ty2 * Z) * 1.0e3; // [mrad] const double vs = i->vs * 1.0e2; // [%] const double ts = sqrt(tx*tx + ty*ty); h2.Fill(ts); Tx += tx; Ty += ty; Vs += vs; Ts += ts; H1 ->Fill(ts, vs); H1[id]->Fill(ts, vs); GO[id].put(t1, atan2(ty, tx)); GA[id].put(t1, ts); GS[id].put(t1, vs); } Tx /= evt->size(); Ty /= evt->size(); Vs /= evt->size(); Ts /= evt->size(); for (JEvt::const_iterator i = evt->begin(); i != evt->end(); ++i) { const double tx = (i->tx + i->tx2 * Z) * 1.0e3; // [mrad] const double ty = (i->ty + i->ty2 * Z) * 1.0e3; // [mrad] H2 ->Fill(tx - Tx, ty - Ty); H2[i->id]->Fill(tx - Tx, ty - Ty); GO[-1].put(t1, atan2(Ty, Tx)); GA[-1].put(t1, Ts); GS[-1].put(t1, Vs); } } STATUS(endl); TH1D hx("hx", NULL, H2.size(), -0.5, H2.size() + 0.5); TH1D hy("hy", NULL, H2.size(), -0.5, H2.size() + 0.5); JQuantile Qx, Qy; for (JManager::const_iterator i = H2.begin(); i != H2.end(); ++i) { const int ix = distance(H2.cbegin(), i) + 1; hx.GetXaxis()->SetBinLabel(ix, MAKE_CSTRING(i->first)); hy.GetXaxis()->SetBinLabel(ix, MAKE_CSTRING(i->first)); hx.SetBinContent(ix, i->second->GetMean(1)); hy.SetBinContent(ix, i->second->GetMean(2)); hx.SetBinError (ix, i->second->GetStdDev(1)); hy.SetBinError (ix, i->second->GetStdDev(2)); Qx.put(i->second->GetMean(1)); Qy.put(i->second->GetMean(2)); } if (Qx.getCount() > 1 && Qy.getCount() > 1) { cout << "deviation: " << FIXED(7,3) << Qx.getSTDev() << ' ' << FIXED(7,3) << Qy.getSTDev() << endl; } TFile out(outputFile.c_str(), "recreate"); out << h1 << h2 << hx << hy << JGraph(g0, "g0") << JGraph(g1, "g1"); out << H1 << *H1 << H2 << *H2; for (map::const_iterator i = GO.begin(); i != GO.end(); ++i) { out << JGraph(i->second, MAKE_CSTRING("G[" << FILL(4,'0') << i->first << "].orientation")); } for (map::const_iterator i = GA.begin(); i != GA.end(); ++i) { out << JGraph(i->second, MAKE_CSTRING("G[" << FILL(4,'0') << i->first << "].amplitude")); } for (map::const_iterator i = GS.begin(); i != GS.end(); ++i) { out << JGraph(i->second, MAKE_CSTRING("G[" << FILL(4,'0') << i->first << "].stretching")); } out.Write(); out.Close(); }