Jpp  master_rocky-40-g5f0272dcd
the software that should make you happy
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 "JROOT/JMinimizer.hh"
#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

◆ main()

int main ( int  argc,
char **  argv 
)

Definition at line 55 of file JHitL1.cc.

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