Jpp
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Functions
JPostfit2F.cc File Reference

Example program to compare fit results from two files. More...

#include <string>
#include <iostream>
#include <iomanip>
#include "TROOT.h"
#include "TFile.h"
#include "TH1D.h"
#include "TH2D.h"
#include "evt/Head.hh"
#include "evt/Evt.hh"
#include "JAAnet/JAAnetToolkit.hh"
#include "JDAQ/JDAQEvent.hh"
#include "JSupport/JTriggeredFileScanner.hh"
#include "JSupport/JMonteCarloFileSupportkit.hh"
#include "JSupport/JSupport.hh"
#include "JFit/JEvt.hh"
#include "JFit/JEvtToolkit.hh"
#include "JFit/JFitParameters.hh"
#include "JMath/JMathToolkit.hh"
#include "Jeep/JParser.hh"
#include "Jeep/JMessage.hh"
#include "JHistogramToolkit.hh"

Go to the source code of this file.

Functions

int main (int argc, char **argv)
 

Detailed Description

Example program to compare fit results from two files.

Author
mdejong

Definition in file JPostfit2F.cc.

Function Documentation

int main ( int  argc,
char **  argv 
)

Definition at line 34 of file JPostfit2F.cc.

35 {
36  using namespace std;
37  using namespace JPP;
38 
39  typedef JTriggeredFileScanner<JEvt> JTriggeredFileScanner_t;
40  typedef JTriggeredFileScanner_t::multi_pointer_type multi_pointer_type;
41 
42  JTriggeredFileScanner_t inputFileA;
43  JTriggeredFileScanner_t inputFileB;
44  JLimit_t numberOfEvents;
45  string outputFile;
46  double angle_Deg;
47  int debug;
48 
49  try {
50 
51  JParser<> zap("Example program to compare fit results from two files.");
52 
53  zap['a'] = make_field(inputFileA);
54  zap['b'] = make_field(inputFileB);
55  zap['o'] = make_field(outputFile) = "postfit2f.root";
56  zap['n'] = make_field(numberOfEvents) = JLimit::max();
57  zap['A'] = make_field(angle_Deg) = 1.0;
58  zap['d'] = make_field(debug) = 2;
59 
60  zap(argc, argv);
61  }
62  catch(const exception& error) {
63  FATAL(error.what() << endl);
64  }
65 
66 
67 
68  inputFileA.setLimit(numberOfEvents);
69  inputFileB.setLimit(numberOfEvents);
70 
71 
72  TFile out(outputFile.c_str(), "recreate");
73 
74  JManager manager;
75 
76  manager.insert(make_key(JQUALITY), 101, 0.0, 250.0);
77  manager.insert(make_key(JGANDALF_BETA0_RAD), 51, 0.0, 2.0e-2);
78  manager.insert(make_key(JGANDALF_BETA1_RAD), 51, 0.0, 2.0e-2);
79  manager.insert(make_key(JGANDALF_NUMBER_OF_HITS), 51, 0.0, 100);
80  manager.insert(make_key(JSTART_NPE_MIP), 51, 0.0, 10.0);
81  manager.insert(make_key(JSTART_NPE_MIP_TOTAL), 51, 0.0, 10.0);
82  manager.insert(make_key(JSTART_LENGTH_METRES), 51, 0.0, 100.0);
83  manager.insert(make_key(JENERGY_ENERGY), 51, 1.0, 5.0, true);
84  manager.insert(make_key(JENERGY_CHI2), 51, 0.0, 200.0);
85 
86  TH1D hA("h[A]", NULL, 100, -3.0, +2.3); // [log(deg)]
87  TH1D hB("h[B]", NULL, 100, -3.0, +2.3); // [log(deg)]
88 
89  TH2D h2("h2", NULL,
90  100, -100.0, +100.0, // quality
91  360, -180.0, +180.0); // deg
92 
93 
94  while (inputFileA.hasNext() && inputFileB.hasNext()) {
95 
96  STATUS("event: " << setw(10) << inputFileA.getCounter() << '\r'); DEBUG(endl);
97 
98  multi_pointer_type psA = inputFileA.next();
99  multi_pointer_type psB = inputFileB.next();
100 
101  // find the same corresponding Monte Carlo true event based on the event identifier.
102 
103  while (psA.get<Evt>()->mc_id < psB.get<Evt>()->mc_id && inputFileA.hasNext()) { psA = inputFileA.next(); }
104  while (psB.get<Evt>()->mc_id < psA.get<Evt>()->mc_id && inputFileB.hasNext()) { psB = inputFileB.next(); }
105 
106  if (!psA.is_valid()) { continue; }
107  if (!psB.is_valid()) { continue; }
108 
109  if (psA.get<Evt>()->mc_id == psB.get<Evt>()->mc_id) {
110 
111  Evt* event = psA;
112  JEvt* evtA = psA;
113  JEvt* evtB = psB;
114 
115  vector<Trk>::const_iterator muon = find_if(event->mc_trks.begin(), event->mc_trks.end(), is_muon);
116 
117  if (muon == event->mc_trks.end()) { continue; }
118 
119  if (evtA->empty()) { continue; }
120  if (evtB->empty()) { continue; }
121 
122  JEvt::const_iterator fitA = evtA->begin();
123  JEvt::const_iterator fitB = evtB->begin();
124 
125  const Double_t angleA = getAngle(getDirection(*muon), getDirection(*fitA));
126  const Double_t angleB = getAngle(getDirection(*muon), getDirection(*fitB));
127 
128  hA.Fill(log10(angleA));
129  hB.Fill(log10(angleB));
130 
131  h2.Fill(fitB->getQ() - fitA->getQ(), angleB - angleA);
132 
133  if (angleA > angle_Deg) {
134  manager.Fill(*fitA, *fitB, angleB < angle_Deg);
135  }
136  }
137  }
138  STATUS(endl);
139 
140  out.Write();
141  out.Close();
142 }
void insert(const JKey &key, const Int_t nx, const Double_t xmin, const Double_t xmax, const double logx=false)
Insert set of histograms.
#define make_key(PARAMETER)
Make key.
Utility class to parse command line options.
Definition: JParser.hh:1410
chi2 from JEnergy.cc
double getAngle(const JFirst_t &first, const JSecond_t &second)
Get space angle between objects.
Auxiliary class to manage set of histograms.
#define STATUS(A)
Definition: JMessage.hh:61
bool is_muon(const Trk &track)
Test whether given track is a (anti-)muon.
Structure to store the ToT mean and standard deviation of the hits produced by a nanobeacon in a sour...
string outputFile
angular resolution [rad] from JGandalf.cc
JLimit JLimit_t
Type definition of limit.
Definition: JLimit.hh:214
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1836
static const int JQUALITY
fit quality identifier
angular resolution [rad] from JGandalf.cc
int debug
debug level
Definition: JSirene.cc:59
uncorrected energy [GeV] from JEnergy.cc
#define FATAL(A)
Definition: JMessage.hh:65
number of hits from JGandalf.cc
number of photo-electrons up to the barycentre from JStart.cc
Data structure for set of track fit results.
Definition: JEvt.hh:312
distance between first and last hits in metres from JStart.cc
JDirection3D getDirection(const Vec &v)
Get direction.
number of photo-electrons along the whole track from JStart.cc
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:60