Jpp  17.3.1
the software that should make you happy
 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 "km3net-dataformat/tools/time_converter.hh"
#include "JDAQ/JDAQEventIO.hh"
#include "JDAQ/JDAQTimesliceIO.hh"
#include "JAAnet/JHead.hh"
#include "JAAnet/JHeadToolkit.hh"
#include "JAAnet/JAAnetToolkit.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 "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 53 of file JHitL1.cc.

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