Jpp
JHitL1.cc
Go to the documentation of this file.
1 
2 #include <string>
3 #include <iostream>
4 #include <iomanip>
5 #include <vector>
6 #include <algorithm>
7 
8 #include "TROOT.h"
9 #include "TFile.h"
10 #include "TH1D.h"
11 #include "TF1.h"
12 
13 #include "evt/Head.hh"
14 #include "evt/Evt.hh"
15 #include "evt/Hit.hh"
16 #include "JDAQ/JDAQEvent.hh"
17 #include "JDAQ/JDAQTimeslice.hh"
18 
19 #include "JAAnet/JHead.hh"
20 #include "JAAnet/JHeadToolkit.hh"
21 
22 #include "JDetector/JDetector.hh"
25 
26 #include "JTrigger/JHit.hh"
27 #include "JTrigger/JFrame.hh"
28 #include "JTrigger/JTimeslice.hh"
30 #include "JTrigger/JHitL0.hh"
31 #include "JTrigger/JHitL1.hh"
32 #include "JTrigger/JBuildL0.hh"
33 #include "JTrigger/JBuildL1.hh"
34 #include "JTrigger/JBuildL2.hh"
35 #include "JTrigger/JAlgorithm.hh"
37 
40 #include "JSupport/JSupport.hh"
41 
42 #include "Jeep/JPrint.hh"
43 #include "Jeep/JParser.hh"
44 #include "Jeep/JMessage.hh"
45 
46 
47 /**
48  * \file
49  * Example program to determine time slewing of L1 hits.
50  * \author mdejong
51  */
52 int main(int argc, char **argv)
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 
88  JDetector detector;
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 }
KM3NETDAQ::JDAQEvent
DAQ Event.
Definition: JDAQEvent.hh:34
JSuperFrame2D.hh
FIXED
Auxiliary data structure for floating point format specification.
Definition: JPrint.hh:481
JMessage.hh
JHead.hh
JPrint.hh
JTimeslice.hh
JDETECTOR::load
void load(const JString &file_name, JDetector &detector)
Load detector from input file.
Definition: JDetectorToolkit.hh:456
std::vector< TH1D * >
JSUPPORT::JLimit_t
JLimit JLimit_t
Type definition of limit.
Definition: JLimit.hh:215
JAANET::getTrack
JTrack3E getTrack(const Trk &track)
Get track.
Definition: JAAnetToolkit.hh:256
JTimeConverter.hh
KM3NETDAQ::NUMBER_OF_PMTS
static const int NUMBER_OF_PMTS
Total number of PMTs in module.
Definition: JDAQ.hh:26
JPARSER::JParser
Utility class to parse command line options.
Definition: JParser.hh:1493
JTriggeredFileScanner.hh
JDAQTimeslice.hh
JAANET::is_muon
bool is_muon(const Trk &track)
Test whether given track is a (anti-)muon.
Definition: JAAnetToolkit.hh:367
JHeadToolkit.hh
JSupport.hh
JPP
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
Definition: JAAnetToolkit.hh:37
JBuildL1.hh
JSUPPORT::getHeader
Head getHeader(const JMultipleFileScanner_t &file_list)
Get Monte Carlo header.
Definition: JMonteCarloFileSupportkit.hh:425
JBuildL0.hh
debug
int debug
debug level
Definition: JSirene.cc:59
JSUPPORT::JLimit::getLimit
const JLimit & getLimit() const
Get limit.
Definition: JLimit.hh:73
main
int main(int argc, char **argv)
Definition: JHitL1.cc:52
JHit.hh
JModuleRouter.hh
STATUS
#define STATUS(A)
Definition: JMessage.hh:63
JBuildL2.hh
JParser.hh
JDetectorToolkit.hh
JHitL1.hh
make_field
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1954
DEBUG
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:62
std
Definition: jaanetDictionary.h:36
KM3NETDAQ
KM3NeT DAQ data structures and auxiliaries.
Definition: DataQueue.cc:39
JDAQEvent.hh
JMonteCarloFileSupportkit.hh
JDetector.hh
JHitL0.hh
FATAL
#define FATAL(A)
Definition: JMessage.hh:67
outputFile
string outputFile
Definition: JDAQTimesliceSelector.cc:37
JEEP::debug_t
debug
Definition: JMessage.hh:29
JFrame.hh
JAlgorithm.hh