Jpp 19.3.0-rc.3
the software that should make you happy
Loading...
Searching...
No Matches
JCompareKatoomba.cc File Reference

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

#include <iostream>
#include <iomanip>
#include "TROOT.h"
#include "TFile.h"
#include "TH1D.h"
#include "TH2D.h"
#include "TProfile.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 compare acoustic fit results.

Author
mdejong

Definition in file JCompareKatoomba.cc.

Function Documentation

◆ main()

int main ( int argc,
char ** argv )

Definition at line 31 of file JCompareKatoomba.cc.

32{
33 using namespace std;
34 using namespace JPP;
35
37 JLimit_t& numberOfEvents = inputFile.getLimit();
38 string outputFile;
39 int debug;
40
41 try {
42
43 JParser<> zap("Example program to compare acoustic fit results.");
44
45 zap['f'] = make_field(inputFile, "input file (output of JKatoomba[.sh])");
46 zap['n'] = make_field(numberOfEvents) = JLimit::max();
47 zap['o'] = make_field(outputFile) = "comparison.root";
48 zap['d'] = make_field(debug) = 2;
49
50 zap(argc, argv);
51 }
52 catch(const exception &error) {
53 FATAL(error.what() << endl);
54 }
55
56
57 if (inputFile.size() != 2u) {
58 FATAL("Invalid number of input files; need 2 files for comparison." << endl);
59 }
60
61 JManager<int, TH2D> H2(new TH2D("string[%]", NULL, 500, -50.0, +50.0, 500, -50.0, +50.0));
62
65
66 while (inA.hasNext() && inB.hasNext()) {
67
68 STATUS("event: " << setw(10) << inA.getCounter() << '\r'); DEBUG(endl);
69
70 JEvt* pA = inA.next();
71 JEvt* pB = inB.next();
72
73 // find the same event based on the start and stop times of the events
74
75 while (pA->UNIXTimeStop < pB->UNIXTimeStart && inA.hasNext()) { pA = inA.next(); }
76 while (pB->UNIXTimeStop < pA->UNIXTimeStart && inB.hasNext()) { pB = inB.next(); }
77
78 if (pA->UNIXTimeStart < pB->UNIXTimeStop &&
79 pB->UNIXTimeStart < pA->UNIXTimeStop) {
80
81 for (JEvt::const_iterator iA = pA->begin(); iA != pA->end(); ++iA) {
82 for (JEvt::const_iterator iB = pB->begin(); iB != pB->end(); ++iB) {
83
84 if (iA->id == iB->id) {
85
86 const double tx = (iA->tx - iB->tx) * 1.0e3; // [mrad]
87 const double ty = (iA->ty - iB->ty) * 1.0e3; // [mrad]
88
89 H2 ->Fill(tx, ty);
90 H2[iA->id]->Fill(tx, ty);
91
92 break;
93 }
94 }
95 }
96 }
97 }
98
99
100 TFile out(outputFile.c_str(), "recreate");
101
102 out << H2 << *H2;
103
104 out.Write();
105 out.Close();
106}
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:72
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition JParser.hh:2142
Utility class to parse command line options.
Definition JParser.hh:1698
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.
Template definition for direct access of elements in ROOT TChain.
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
Acoustic event fit.
double UNIXTimeStop
stop time
double UNIXTimeStart
start time
Auxiliary class for defining the range of iterations of objects.
Definition JLimit.hh:45
static counter_type max()
Get maximum counter value.
Definition JLimit.hh:128