Jpp  15.0.1-rc.1-highQE
the software that should make you happy
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 
17 #include "JDAQ/JDAQEventIO.hh"
18 #include "JDAQ/JDAQTimesliceIO.hh"
19 
20 #include "JAAnet/JHead.hh"
21 #include "JAAnet/JHeadToolkit.hh"
22 #include "JAAnet/JAAnetToolkit.hh"
23 
24 #include "JDetector/JDetector.hh"
27 
28 #include "JTrigger/JHit.hh"
29 #include "JTrigger/JFrame.hh"
30 #include "JTrigger/JTimeslice.hh"
32 #include "JTrigger/JHitL0.hh"
33 #include "JTrigger/JHitL1.hh"
34 #include "JTrigger/JBuildL0.hh"
35 #include "JTrigger/JBuildL1.hh"
36 #include "JTrigger/JBuildL2.hh"
37 #include "JTrigger/JAlgorithm.hh"
38 
41 #include "JSupport/JSupport.hh"
42 
43 #include "Jeep/JPrint.hh"
44 #include "Jeep/JParser.hh"
45 #include "Jeep/JMessage.hh"
46 
47 
48 /**
49  * \file
50  * Example program to determine time slewing of L1 hits.
51  * \author mdejong
52  */
53 int main(int argc, char **argv)
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 
86  cout.tie(&cerr);
87 
88 
90 
91  try {
92  load(detectorFile, detector);
93  }
94  catch(const JException& error) {
95  FATAL(error);
96  }
97 
98  const JModuleRouter router(detector);
99 
100  detector -= get<JPosition3D>(getHeader(inputFile));
101 
102 
103  TFile out(outputFile.c_str(), "recreate");
104 
105  vector<TH1D*> H1D;
106 
107  const int nx = 100;
108  const double xmin = -25.0; // [ns]
109  const double xmax = +25.0; // [ns]
110 
111  for (int i = 0; i <= NUMBER_OF_PMTS; ++i) {
112 
113  ostringstream os;
114 
115  os << "M[" << setfill('0') << setw(2) << i << "]";
116 
117  H1D.push_back(new TH1D(os.str().c_str(), NULL, nx, xmin, xmax));
118  }
119 
120 
121  typedef vector<JHitL0> JDataL0_t;
122  typedef vector<JHitL2> JDataL2_t;
123 
124  JSuperFrame2D <JHit> buffer; // I/O buffer for time-sorted calibrated hits
125  const JBuildL0<JHitL0> buildL0;
126  const JBuildL2<JHitL2> buildL2(JL2Parameters(2, Tmax_ns, ctMin));
127 
128 
129  while (inputFile.hasNext()) {
130 
131  STATUS("event: " << setw(10) << inputFile.getCounter() << '\r'); DEBUG(endl);
132 
133  JTriggeredFileScanner<>::multi_pointer_type ps = inputFile.next();
134 
135  const JDAQEvent* tev = ps;
136  const Evt* event = ps;
137 
138  const time_converter converter(*event, *tev);
139 
140  vector<Trk>::const_iterator muon = find_if(event->mc_trks.begin(), event->mc_trks.end(), is_muon);
141 
142  if (muon != event->mc_trks.end()) {
143 
144  // muon
145 
146  const JTrack3D trk = getTrack(*muon);
147 
148  JDataL0_t dataL0;
149  JDataL2_t dataL2;
150 
151  buildL0(*tev, router, true, back_inserter(dataL0));
152  buildL2(*tev, router, true, back_inserter(dataL2));
153 
154  // select first hit
155 
156  sort(dataL2.begin(), dataL2.end(), timeSorter<JHitL2>);
157 
158  dataL2.erase(unique(dataL2.begin(), dataL2.end(), equal_to<JDAQModuleIdentifier>()), dataL2.end());
159 
160  for (JDataL0_t::const_iterator hit = dataL0.begin(); hit != dataL0.end(); ++hit) {
161 
162  const double t0 = trk.getT(*hit);
163  const double t1 = converter.getTime(hit->getT());
164 
165  H1D[1]->Fill(t1 - t0);
166  }
167 
168 
169  for (JDataL2_t::const_iterator hit = dataL2.begin(); hit != dataL2.end(); ++hit) {
170 
171  const double t0 = trk.getT(*hit);
172  const double t1 = converter.getTime(hit->begin()->getT());
173  const size_t index = min(hit->size(), H1D.size() - 1);
174 
175  H1D[index]->Fill(t1 - t0);
176  }
177  }
178  }
179  STATUS(endl);
180 
181 
182  TF1 f1("f1", "[0]*exp(-0.5*(x-[1])*(x-[1])/([2]*[2])) + [3]");
183 
184  string option = "LL";
185 
186  if (debug < JEEP::debug_t && option.find('Q') == string::npos) {
187  option += "Q";
188  }
189 
190 
191  for (vector<TH1D*>::iterator h1 = H1D.begin(); h1 != H1D.end(); ++h1) {
192 
193  // start values
194 
195  Double_t ymax = 0.0;
196  Double_t x0 = 0.0; // [ns]
197  Double_t sigma = 2.0; // [ns]
198  Double_t ymin = 0.0;
199 
200  for (int i = 1; i != (*h1)->GetNbinsX(); ++i) {
201 
202  const Double_t x = (*h1)->GetBinCenter (i);
203  const Double_t y = (*h1)->GetBinContent(i);
204 
205  if (y > ymax) {
206  ymax = y;
207  x0 = x;
208  }
209  }
210 
211  f1.SetParameter(0, ymax);
212  f1.SetParameter(1, x0);
213  f1.SetParameter(2, sigma);
214  f1.FixParameter(3, ymin);
215 
216 
217  // fit
218 
219  (*h1)->Fit(&f1, option.c_str(), "same", x0 - 3 * sigma, x0 + 3 * sigma);
220 
221 
222  cout << "\tt0.push_back(" << showpos << FIXED(4,2) << f1.GetParameter(1) << ");" << endl;
223  }
224 
225  out.Write();
226  out.Close();
227 }
Utility class to parse command line options.
Definition: JParser.hh:1500
General exception.
Definition: JException.hh:23
debug
Definition: JMessage.hh:29
int main(int argc, char *argv[])
Definition: Main.cc:15
ROOT TTree parameter settings of various packages.
JTrack3E getTrack(const Trk &track)
Get track.
Algorithms for hit clustering and sorting.
Synchronously read DAQ events and Monte Carlo events (and optionally other events).
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)...
then for HISTOGRAM in h0 h1
Definition: JMatrixNZ.sh:71
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
Basic data structure for time and time over threshold information of hit.
string outputFile
Data structure for detector geometry and calibration.
Head getHeader(const JMultipleFileScanner_t &file_list)
Get Monte Carlo header.
Basic data structure for L0 hit.
Auxiliary class for defining the range of iterations of objects.
Definition: JLimit.hh:41
I/O formatting auxiliaries.
Detector file.
Definition: JHead.hh:196
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1961
Template specialisation of L2 builder for JHitL2 data type.
Definition: JBuildL2.hh:219
Definition of hit and track types and auxiliary methods for handling Monte Carlo data.
int debug
debug level
Definition: JSirene.cc:63
double getT(const JVector3D &pos) const
Get arrival time of Cherenkov light at given position.
Definition: JTrack3D.hh:126
General purpose messaging.
Auxiliary include file for time conversion between DAQ/trigger hit and Monte Carlo hit...
#define FATAL(A)
Definition: JMessage.hh:67
Direct access to module in detector data structure.
Data structure for L2 parameters.
void load(const std::string &file_name, JDetector &detector)
Load detector from input file.
Utility class to parse command line options.
then usage $script< input_file >< detector_file > fi set_variable OUTPUT_DIR set_variable SELECTOR JDAQTimesliceL1 set_variable DEBUG case set_variable DEBUG
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.
double getTime() const
Get DAQ/trigger time minus Monte Carlo time.
Basic data structure for L1 hit.
The Evt class respresent a Monte Carlo (MC) event as well as an offline event.
Definition: Evt.hh:19