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

Example program to determine time slewing of L1 hits. More...

#include <string>
#include <iostream>
#include <iomanip>
#include <vector>
#include <algorithm>
#include "TROOT.h"
#include "TFile.h"
#include "TH1D.h"
#include "TF1.h"
#include "km3net-dataformat/offline/Head.hh"
#include "km3net-dataformat/offline/Evt.hh"
#include "km3net-dataformat/offline/Hit.hh"
#include "JDAQ/JDAQEventIO.hh"
#include "JDAQ/JDAQTimesliceIO.hh"
#include "JAAnet/JHead.hh"
#include "JAAnet/JHeadToolkit.hh"
#include "JDetector/JDetector.hh"
#include "JDetector/JDetectorToolkit.hh"
#include "JDetector/JModuleRouter.hh"
#include "JTrigger/JHit.hh"
#include "JTrigger/JFrame.hh"
#include "JTrigger/JTimeslice.hh"
#include "JTrigger/JSuperFrame2D.hh"
#include "JTrigger/JHitL0.hh"
#include "JTrigger/JHitL1.hh"
#include "JTrigger/JBuildL0.hh"
#include "JTrigger/JBuildL1.hh"
#include "JTrigger/JBuildL2.hh"
#include "JTrigger/JAlgorithm.hh"
#include "JTrigger/JTimeConverter.hh"
#include "JSupport/JTriggeredFileScanner.hh"
#include "JSupport/JMonteCarloFileSupportkit.hh"
#include "JSupport/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 determine time slewing of L1 hits.

Author
mdejong

Definition in file JHitL1.cc.

Function Documentation

int main ( int  argc,
char **  argv 
)

Definition at line 52 of file JHitL1.cc.

53 {
54  using namespace std;
55  using namespace JPP;
56  using namespace KM3NETDAQ;
57 
58  JTriggeredFileScanner<> inputFile;
59  JLimit_t& numberOfEvents = inputFile.getLimit();
60  string outputFile;
61  string detectorFile;
62  double Tmax_ns;
63  double ctMin;
64  int debug;
65 
66  try {
67 
68  JParser<> zap("Example program to determine time slewing of L1 hits.");
69 
70  zap['f'] = make_field(inputFile);
71  zap['o'] = make_field(outputFile) = "hitL1.root";
72  zap['a'] = make_field(detectorFile);
73  zap['n'] = make_field(numberOfEvents) = JLimit::max();
74  zap['T'] = make_field(Tmax_ns) = 15.0; // [ns]
75  zap['C'] = make_field(ctMin) = 0.0; //
76  zap['d'] = make_field(debug) = 1;
77 
78  zap(argc, argv);
79  }
80  catch(const exception& error) {
81  FATAL(error.what() << endl);
82  }
83 
84 
85  cout.tie(&cerr);
86 
87 
89 
90  try {
91  load(detectorFile, detector);
92  }
93  catch(const JException& error) {
94  FATAL(error);
95  }
96 
97  const JModuleRouter router(detector);
98 
99  detector -= get<JPosition3D>(getHeader(inputFile));
100 
101 
102  TFile out(outputFile.c_str(), "recreate");
103 
104  vector<TH1D*> H1D;
105 
106  const int nx = 100;
107  const double xmin = -25.0; // [ns]
108  const double xmax = +25.0; // [ns]
109 
110  for (int i = 0; i <= NUMBER_OF_PMTS; ++i) {
111 
112  ostringstream os;
113 
114  os << "M[" << setfill('0') << setw(2) << i << "]";
115 
116  H1D.push_back(new TH1D(os.str().c_str(), NULL, nx, xmin, xmax));
117  }
118 
119 
120  typedef vector<JHitL0> JDataL0_t;
121  typedef vector<JHitL2> JDataL2_t;
122 
123  JSuperFrame2D <JHit> buffer; // I/O buffer for time-sorted calibrated hits
124  const JBuildL0<JHitL0> buildL0;
125  const JBuildL2<JHitL2> buildL2(JL2Parameters(2, Tmax_ns, ctMin));
126 
127 
128  while (inputFile.hasNext()) {
129 
130  STATUS("event: " << setw(10) << inputFile.getCounter() << '\r'); DEBUG(endl);
131 
132  JTriggeredFileScanner<>::multi_pointer_type ps = inputFile.next();
133 
134  const JDAQEvent* tev = ps;
135  const Evt* event = ps;
136 
137  const JTimeConverter converter(*event, *tev);
138 
139  vector<Trk>::const_iterator muon = find_if(event->mc_trks.begin(), event->mc_trks.end(), is_muon);
140 
141  if (muon != event->mc_trks.end()) {
142 
143  // muon
144 
145  const JTrack3D trk = getTrack(*muon);
146 
147  JDataL0_t dataL0;
148  JDataL2_t dataL2;
149 
150  buildL0(*tev, router, true, back_inserter(dataL0));
151  buildL2(*tev, router, true, back_inserter(dataL2));
152 
153  // select first hit
154 
155  sort(dataL2.begin(), dataL2.end(), timeSorter<JHitL2>);
156 
157  dataL2.erase(unique(dataL2.begin(), dataL2.end(), equal_to<JDAQModuleIdentifier>()), dataL2.end());
158 
159  for (JDataL0_t::const_iterator hit = dataL0.begin(); hit != dataL0.end(); ++hit) {
160 
161  const double t0 = trk.getT(*hit);
162  const double t1 = converter.getTime(hit->getT());
163 
164  H1D[1]->Fill(t1 - t0);
165  }
166 
167 
168  for (JDataL2_t::const_iterator hit = dataL2.begin(); hit != dataL2.end(); ++hit) {
169 
170  const double t0 = trk.getT(*hit);
171  const double t1 = converter.getTime(hit->begin()->getT());
172  const size_t index = min(hit->size(), H1D.size() - 1);
173 
174  H1D[index]->Fill(t1 - t0);
175  }
176  }
177  }
178  STATUS(endl);
179 
180 
181  TF1 f1("f1", "[0]*exp(-0.5*(x-[1])*(x-[1])/([2]*[2])) + [3]");
182 
183  string option = "LL";
184 
185  if (debug < JEEP::debug_t && option.find('Q') == string::npos) {
186  option += "Q";
187  }
188 
189 
190  for (vector<TH1D*>::iterator h1 = H1D.begin(); h1 != H1D.end(); ++h1) {
191 
192  // start values
193 
194  Double_t ymax = 0.0;
195  Double_t x0 = 0.0; // [ns]
196  Double_t sigma = 2.0; // [ns]
197  Double_t ymin = 0.0;
198 
199  for (int i = 1; i != (*h1)->GetNbinsX(); ++i) {
200 
201  const Double_t x = (*h1)->GetBinCenter (i);
202  const Double_t y = (*h1)->GetBinContent(i);
203 
204  if (y > ymax) {
205  ymax = y;
206  x0 = x;
207  }
208  }
209 
210  f1.SetParameter(0, ymax);
211  f1.SetParameter(1, x0);
212  f1.SetParameter(2, sigma);
213  f1.FixParameter(3, ymin);
214 
215 
216  // fit
217 
218  (*h1)->Fit(&f1, option.c_str(), "same", x0 - 3 * sigma, x0 + 3 * sigma);
219 
220 
221  cout << "\tt0.push_back(" << showpos << FIXED(4,2) << f1.GetParameter(1) << ");" << endl;
222  }
223 
224  out.Write();
225  out.Close();
226 }
Utility class to parse command line options.
Definition: JParser.hh:1493
General exception.
Definition: JException.hh:23
debug
Definition: JMessage.hh:29
JTrack3E getTrack(const Trk &track)
Get track.
Template specialisation of L0 builder for JHitL0 data type.
Definition: JBuildL0.hh:102
#define STATUS(A)
Definition: JMessage.hh:63
Detector data structure.
Definition: JDetector.hh:80
bool is_muon(const Trk &track)
Test whether given track is a (anti-)muon.
Router for direct addressing of module data in detector data structure.
Auxiliary class to synchronously read DAQ events and Monte Carlo events (and optionally other events)...
then for HISTOGRAM in h0 h1
Definition: JMatrixNZ.sh:69
Auxiliary data structure for floating point format specification.
Definition: JPrint.hh:481
string outputFile
Head getHeader(const JMultipleFileScanner_t &file_list)
Get Monte Carlo header.
Auxiliary class for defining the range of iterations of objects.
Definition: JLimit.hh:41
Detector file.
Definition: JHead.hh:130
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1954
Template specialisation of L2 builder for JHitL2 data type.
Definition: JBuildL2.hh:186
void load(const JString &file_name, JDetector &detector)
Load detector from input file.
int debug
debug level
Definition: JSirene.cc:61
double getT(const JVector3D &pos) const
Get arrival time of Cherenkov light at given position.
Definition: JTrack3D.hh:126
#define FATAL(A)
Definition: JMessage.hh:67
Data structure for L2 parameters.
2-dimensional frame with time calibrated data from one optical module.
Auxiliary class to convert DAQ/trigger hit time to/from Monte Carlo hit time.
const JLimit & getLimit() const
Get limit.
Definition: JLimit.hh:73
static const int NUMBER_OF_PMTS
Total number of PMTs in module.
Definition: JDAQ.hh:26
General purpose class for multiple pointers.
The Evt class respresent a Monte Carlo (MC) event as well as an offline event.
Definition: Evt.hh:19
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:62