Jpp  17.0.0-rc.1
the software that should make you happy
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Functions
JParramatta.cc File Reference

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

#include <iostream>
#include <iomanip>
#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 "JAcoustics/JEvt.hh"
#include "JAcoustics/JSupport.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

int main ( int  argc,
char **  argv 
)

Definition at line 32 of file JParramatta.cc.

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 
58  TH1D h1("chi2/NDF", NULL, 500, 0.0, 10.0);
59 
60  JGraph_t g0;
61  JGraph_t g1;
62 
66 
67  JManager<int, TProfile> H1(new TProfile("stretching[%]", NULL, 50, 0.0, 1.0e2));
68 
69  while (inputFile.hasNext()) {
70 
71  STATUS("event: " << setw(10) << inputFile.getCounter() << '\r'); DEBUG(endl);
72 
73  const JEvt* evt = inputFile.next();
74 
75  h1.Fill(evt->chi2 / evt->ndf);
76 
77  const double t1 = 0.5 * (evt->UNIXTimeStart + evt->UNIXTimeStop);
78 
79  g0.put(t1, evt->nhit - evt->npar);
80  g1.put(t1, evt->chi2 / evt->ndf);
81 
82  for (JEvt::const_iterator i = evt->begin(); i != evt->end(); ++i) {
83 
84  const int id = i->id;
85  const double tx = i->tx * 1.0e3; // [mrad]
86  const double ty = i->ty * 1.0e3; // [mrad]
87  const double vs = i->vs * 1.0e2; // [%]
88  const double ts = sqrt(tx*tx + ty*ty);
89 
90  H1 ->Fill(ts, vs);
91  H1[id]->Fill(ts, vs);
92 
93  GO[id].put(t1, atan2(ty, tx));
94  GA[id].put(t1, ts);
95  GS[id].put(t1, vs);
96  }
97  }
98  STATUS(endl);
99 
100 
101  TFile out(outputFile.c_str(), "recreate");
102 
103  out << h1
104  << JGraph(g0, "g0")
105  << JGraph(g1, "g1");
106 
107  out << H1 << *H1;
108 
109  for (map<int, JGraph_t>::const_iterator i = GO.begin(); i != GO.end(); ++i) {
110  out << JGraph(i->second, MAKE_CSTRING("G[" << FILL(4,'0') << i->first << "].orientation"));
111  }
112 
113  for (map<int, JGraph_t>::const_iterator i = GA.begin(); i != GA.end(); ++i) {
114  out << JGraph(i->second, MAKE_CSTRING("G[" << FILL(4,'0') << i->first << "].amplitude"));
115  }
116 
117  for (map<int, JGraph_t>::const_iterator i = GS.begin(); i != GS.end(); ++i) {
118  out << JGraph(i->second, MAKE_CSTRING("G[" << FILL(4,'0') << i->first << "].stretching"));
119  }
120 
121  out.Write();
122  out.Close();
123 }
Utility class to parse command line options.
Definition: JParser.hh:1500
int npar
number of fit parameters
Data structure for graph data.
Definition: JGraph.hh:21
#define STATUS(A)
Definition: JMessage.hh:63
#define MAKE_CSTRING(A)
Make C-string.
Definition: JPrint.hh:151
then for HISTOGRAM in h0 h1
Definition: JMatrixNZ.sh:71
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
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:1961
int debug
debug level
Definition: JSirene.cc:66
Auxiliary data structure for sequence of same character.
Definition: JManip.hh:328
#define FATAL(A)
Definition: JMessage.hh:67
int nhit
number of hits
double UNIXTimeStart
start time
General purpose class for object reading from a list of file names.
const JLimit & getLimit() const
Get limit.
Definition: JLimit.hh:73
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:62
Double_t g1(const Double_t x)
Function.
Definition: JQuantiles.cc:25