Jpp  master_rocky-40-g5f0272dcd
the software that should make you happy
Functions
JCalibrateMuon.cc File Reference

Program to perform detector calibration using reconstructed muon trajectories. More...

#include <string>
#include <iostream>
#include <iomanip>
#include <vector>
#include <set>
#include <algorithm>
#include <memory>
#include "TROOT.h"
#include "TFile.h"
#include "TProfile.h"
#include "TH2D.h"
#include "km3net-dataformat/offline/Head.hh"
#include "km3net-dataformat/offline/Evt.hh"
#include "JDAQ/JDAQEventIO.hh"
#include "JDAQ/JDAQTimesliceIO.hh"
#include "JDAQ/JDAQSummarysliceIO.hh"
#include "JDetector/JDetector.hh"
#include "JDetector/JDetectorToolkit.hh"
#include "JDetector/JModuleRouter.hh"
#include "JDynamics/JDynamics.hh"
#include "JLang/JPredicate.hh"
#include "JTrigger/JHit.hh"
#include "JTrigger/JFrame.hh"
#include "JTrigger/JTimeslice.hh"
#include "JTrigger/JHitL0.hh"
#include "JTrigger/JBuildL0.hh"
#include "JTrigger/JTriggerParameters.hh"
#include "JSupport/JSupport.hh"
#include "JSupport/JParallelFileScanner.hh"
#include "JSupport/JSummaryFileRouter.hh"
#include "JSupport/JMeta.hh"
#include "JReconstruction/JEvt.hh"
#include "JReconstruction/JMuonGandalf.hh"
#include "JReconstruction/JMuonGandalfParameters_t.hh"
#include "JTools/JAbstractHistogram.hh"
#include "JROOT/JRootToolkit.hh"
#include "JROOT/JManager.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

Program to perform detector calibration using reconstructed muon trajectories.

It is recommended to perform the detector calibration after JMuonGandalf.cc and to accordingly set option JMuonGandalfParameters_t::numberOfPrefits = 1, so that the best track from the previous output is used. Several histograms will be produced which can subsequently be processed with JHobbit.cc for the determination of the time calibration of optical modules and JFrodo.cc for the time-slewing correction.

Author
mdejong

Definition in file JCalibrateMuon.cc.

Function Documentation

◆ main()

int main ( int  argc,
char **  argv 
)

Definition at line 66 of file JCalibrateMuon.cc.

67 {
68  using namespace std;
69  using namespace JPP;
70  using namespace KM3NETDAQ;
71 
72  typedef JParallelFileScanner< JTypeList<JDAQEvent, JFIT::JEvt> > JParallelFileScanner_t;
73  typedef JParallelFileScanner_t::multi_pointer_type multi_pointer_type;
75  JACOUSTICS::JEvt>::typelist calibration_types;
76  typedef JMultipleFileScanner<calibration_types> JCalibration_t;
77  typedef JAbstractHistogram<double> histogram_type;
78 
79  JParallelFileScanner_t inputFile;
80  string outputFile;
81  JLimit_t& numberOfEvents = inputFile.getLimit();
82  string detectorFile;
83  JCalibration_t calibrationFile;
84  double Tmax_s;
85  string pdfFile;
86  JMuonGandalfParameters_t parameters;
87  histogram_type calibrate;
88  int debug;
89 
90  try {
91 
92  parameters.numberOfPrefits = 1;
93 
94  JParser<> zap("Program to perform detector calibration using reconstructed muon trajectories.");
95 
96  zap['f'] = make_field(inputFile);
97  zap['o'] = make_field(outputFile) = "gandalf.root";
98  zap['a'] = make_field(detectorFile);
99  zap['+'] = make_field(calibrationFile) = JPARSER::initialised();
100  zap['T'] = make_field(Tmax_s) = 100.0;
101  zap['n'] = make_field(numberOfEvents) = JLimit::max();
102  zap['P'] = make_field(pdfFile);
103  zap['c'] = make_field(calibrate, "Histogram for time calibration per optical module.");
104  zap['@'] = make_field(parameters) = JPARSER::initialised();
105  zap['d'] = make_field(debug) = 1;
106 
107  zap(argc, argv);
108  }
109  catch(const exception& error) {
110  FATAL(error.what() << endl);
111  }
112 
113 
114  if (!calibrate.is_valid()) {
115  FATAL("Invalid calibration data " << calibrate << endl);
116  }
117 
118  if (parameters.numberOfPrefits != 1) {
119  WARNING("Number of prefits " << parameters.numberOfPrefits << " != " << 1 << endl);
120  }
121 
123 
124  try {
125  load(detectorFile, detector);
126  }
127  catch(const JException& error) {
128  FATAL(error);
129  }
130 
131  unique_ptr<JDynamics> dynamics;
132 
133  try {
134 
135  dynamics.reset(new JDynamics(detector, Tmax_s));
136 
137  dynamics->load(calibrationFile);
138  }
139  catch(const exception& error) {
140  if (!calibrationFile.empty()) {
141  FATAL(error.what());
142  }
143  }
144 
145  const JModuleRouter router(dynamics ? dynamics->getDetector() : detector);
146 
147  JSummaryFileRouter summary(inputFile, parameters.R_Hz);
148 
149  JMuonGandalf fit(parameters, router, summary, pdfFile, debug);
150 
151  typedef vector<JHitL0> JDataL0_t;
152  typedef vector<JHitW0> JDataW0_t;
153 
154  const JBuildL0<JHitL0> buildL0;
155 
156 
157  typedef JManager<string, TProfile> JManager_t;
158 
159  TH2D ha("ha", NULL, 256, -0.5, 255.5,
160  calibrate.getNumberOfBins(), calibrate.getLowerLimit(), calibrate.getUpperLimit());
161 
162  TProfile hb("hb", NULL, 256, -0.5, 255.5);
163 
164  TProfile hr("hr", NULL, 60, 0.0, 150.0);
165 
166  TH2D h2("h2", NULL, detector.size(), -0.5, detector.size() - 0.5,
167  calibrate.getNumberOfBins(), calibrate.getLowerLimit(), calibrate.getUpperLimit());
168 
169  JManager_t g1(new TProfile("%", NULL, detector.size(), -0.5, detector.size() - 0.5, -1.0, +1.0));
170 
171 
172  while (inputFile.hasNext()) {
173 
174  STATUS("event: " << setw(10) << inputFile.getCounter() << '\r'); DEBUG(endl);
175 
176  multi_pointer_type ps = inputFile.next();
177 
178  JDAQEvent* tev = ps;
179  JFIT::JEvt* in = ps;
180 
181  summary.update(*tev);
182 
183  if (dynamics) {
184  dynamics->update(*tev);
185  }
186 
187  in->select(parameters.numberOfPrefits, qualitySorter);
188 
189  JDataL0_t dataL0;
190 
191  buildL0(*tev, router, true, back_inserter(dataL0));
192 
193  for (JFIT::JEvt::const_iterator track = in->begin(); track != in->end(); ++track) {
194 
195  const JRotation3D R (getDirection(*track));
196  const JLine1Z tz(getPosition (*track).rotate(R), track->getT());
197  JRange<double> Z_m;
198 
199  if (track->hasW(JSTART_LENGTH_METRES)) {
200  Z_m = JZRange(0.0, track->getW(JSTART_LENGTH_METRES)) + JRange<double>(parameters.ZMin_m, parameters.ZMax_m);
201  }
202 
203  const JFIT::JModel<JLine1Z> match(tz, parameters.roadWidth_m, JRange<double>(parameters.TMin_ns, parameters.TMax_ns), Z_m);
204 
205  // hit selection based on start value
206 
207  JDataW0_t data;
208 
209  for (JDataL0_t::const_iterator i = dataL0.begin(); i != dataL0.end(); ++i) {
210 
211  double rate_Hz = summary.getRate(*i);
212 
213  if (rate_Hz <= 0.0) {
214  rate_Hz = summary.getRate();
215  }
216 
217  JHitW0 hit(*i, rate_Hz);
218 
219  hit.rotate(R);
220 
221  if (match(hit)) {
222  data.push_back(hit);
223  }
224  }
225 
226  // select first hit
227 
228  sort(data.begin(), data.end(), JHitL0::compare);
229 
230  JDataW0_t::iterator __end = unique(data.begin(), data.end(), equal_to<JDAQPMTIdentifier>());
231 
232 
233  // set fit parameters
234 
235  if (track->getE() > 0.1)
236  fit.JRegressor_t::E_GeV = track->getE();
237  else
238  fit.JRegressor_t::E_GeV = parameters.E_GeV;
239 
240  set<int> buffer;
241 
242  for (JDataW0_t::const_iterator hit = data.begin(); hit != __end; ++hit) {
243  buffer.insert(hit->getModuleID());
244  }
245 
246  const JLine3Z ta(tz);
247 
248  for (set<int>::const_iterator id = buffer.begin(); id != buffer.end(); ++id) {
249 
250  JDataW0_t::iterator q = partition(data.begin(), __end, make_predicate(&JHitW0::getModuleID, *id, JComparison::ne()));
251 
252  if (distance(data.begin(), q) - fit.parameters.size() > 0) {
253 
254  fit(ta, data.begin(), q); // re-fit
255 
256  for (JDataW0_t::const_iterator hit = q; hit != __end; ++hit) {
257 
258  const int index = router.getIndex(*id);
259  const double t1 = fit.value.getT(hit->getPosition());
260  JTrack3D gradient = fit(fit.value, *hit).gradient;
261 
262  gradient.rotate_back(R);
263 
264  ha.Fill(hit->getToT(), hit->getT1() - t1);
265  hb.Fill(hit->getToT(), gradient.getT());
266 
267  hr.Fill(fit.value.getDistance(*hit), gradient.getT());
268 
269  h2.Fill(index, hit->getT() - t1);
270 
271  g1["T"]->Fill(index, gradient.getT());
272  g1["X"]->Fill(index, gradient.getX());
273  g1["Y"]->Fill(index, gradient.getY());
274  g1["Z"]->Fill(index, gradient.getZ());
275  }
276  }
277  }
278  }
279  }
280  STATUS(endl);
281 
282  TFile out(outputFile.c_str(), "recreate");
283 
284  putObject(out, JMeta(argc, argv));
285 
286  out << ha;
287  out << hb;
288  out << hr;
289  out << h2;
290  out << g1;
291 
292  out.Write();
293  out.Close();
294 }
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 WARNING(A)
Definition: JMessage.hh:65
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:2142
Double_t g1(const Double_t x)
Function.
Definition: JQuantiles.cc:25
std::vector< T >::difference_type distance(typename std::vector< T >::const_iterator first, typename PhysicsEvent::const_iterator< T > second)
Specialisation of STL distance.
Detector data structure.
Definition: JDetector.hh:96
Router for direct addressing of module data in detector data structure.
Data structure for set of track fit results.
void select(const JSelector_t &selector)
Select fits.
Data structure for fit of straight line paralel to z-axis.
Definition: JLine1Z.hh:29
Data structure for fit of straight line in positive z-direction.
Definition: JLine3Z.hh:40
JAxis3D & rotate_back(const JRotation3D &R)
Rotate back axis.
Definition: JAxis3D.hh:240
JPosition3D & rotate(const JRotation3D &R)
Rotate.
Definition: JPosition3D.hh:186
Rotation matrix.
Definition: JRotation3D.hh:114
double getT(const JVector3D &pos) const
Get arrival time of Cherenkov light at given position.
Definition: JTrack3D.hh:126
double getY() const
Get y position.
Definition: JVector3D.hh:104
double getZ() const
Get z position.
Definition: JVector3D.hh:115
double getX() const
Get x position.
Definition: JVector3D.hh:94
General exception.
Definition: JException.hh:24
Utility class to parse command line options.
Definition: JParser.hh:1698
Auxiliary class for a hit with background rate value.
Definition: JHitW0.hh:23
Auxiliary class to manage set of compatible ROOT objects (e.g. histograms) using unique keys.
Definition: JManager.hh:47
General purpose class for object reading from a list of file names.
General purpose class for parallel reading of objects from a single file or multiple files.
File router for fast addressing of summary data.
Template specialisation of L0 builder for JHitL0 data type.
Definition: JBuildL0.hh:105
static const int JSTART_LENGTH_METRES
distance between first and last hits in metres from JStart.cc
JDirection3D getDirection(const Vec &dir)
Get direction.
JPosition3D getPosition(const Vec &pos)
Get position.
void load(const std::string &file_name, JDetector &detector)
Load detector from input file.
JTOOLS::JRange< double > JZRange
Definition: JFit/JModel.hh:21
JPredicate< JResult_t T::*, JComparison::eq > make_predicate(JResult_t T::*member, const JResult_t value)
Helper method to create predicate for data member.
Definition: JPredicate.hh:128
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
bool qualitySorter(const JFit &first, const JFit &second)
Comparison of fit results.
bool putObject(TDirectory &dir, const TObject &object)
Write object to ROOT directory.
KM3NeT DAQ data structures and auxiliaries.
Definition: DataQueue.cc:39
Definition: JSTDTypes.hh:14
Detector file.
Definition: JHead.hh:227
Acoustic event fit.
Orientation of module.
Dynamic detector calibration.
Definition: JDynamics.hh:84
Template specialisation of class JModel to match hit with muon trajectory along z-axis.
Definition: JFit/JModel.hh:36
Auxiliary class for recursive type list generation.
Definition: JTypeList.hh:351
Empty structure for specification of parser element that is initialised (i.e. does not require input)...
Definition: JParser.hh:68
double TMin_ns
minimal time w.r.t. Cherenkov hypothesis [ns]
double TMax_ns
maximal time w.r.t. Cherenkov hypothesis [ns]
Wrapper class to make final fit of muon trajectory.
Definition: JMuonGandalf.hh:75
Auxiliary class for defining the range of iterations of objects.
Definition: JLimit.hh:45
const JLimit & getLimit() const
Get limit.
Definition: JLimit.hh:84
Auxiliary class for ROOT I/O of application specific meta data.
Definition: JMeta.hh:72
Simple data structure for histogram binning.
Auxiliary data structure for sorting of hits.
Definition: JHitL0.hh:85