Jpp test-rotations-new
the software that should make you happy
Loading...
Searching...
No Matches
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 "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 "JCalibrate/JCalibrateMuon.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 65 of file JCalibrateMuon.cc.

66{
67 using namespace std;
68 using namespace JPP;
69 using namespace KM3NETDAQ;
70
71 typedef JParallelFileScanner< JTypeList<JDAQEvent, JFIT::JEvt> > JParallelFileScanner_t;
72 typedef JParallelFileScanner_t::multi_pointer_type multi_pointer_type;
74 JACOUSTICS::JEvt>::typelist calibration_types;
75 typedef JMultipleFileScanner<calibration_types> JCalibration_t;
76 typedef JAbstractHistogram<double> histogram_type;
77
78 JParallelFileScanner_t inputFile;
79 string outputFile;
80 JLimit_t& numberOfEvents = inputFile.getLimit();
81 string detectorFile;
82 JCalibration_t calibrationFile;
83 double Tmax_s;
84 string pdfFile;
85 JMuonGandalfParameters_t parameters;
86 histogram_type calibrate;
87 string target;
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.");
104 zap['@'] = make_field(parameters) = JPARSER::initialised();
105 zap['R'] = make_field(target, "option for time calibration.") = module_t, string_t;
106 zap['d'] = make_field(debug) = 1;
107
108 zap(argc, argv);
109 }
110 catch(const exception& error) {
111 FATAL(error.what() << endl);
112 }
113
114
115 if (!calibrate.is_valid()) {
116 FATAL("Invalid calibration data " << calibrate << endl);
117 }
118
119 if (parameters.numberOfPrefits != 1) {
120 WARNING("Number of prefits " << parameters.numberOfPrefits << " != " << 1 << endl);
121 }
122
124
125 try {
126 load(detectorFile, detector);
127 }
128 catch(const JException& error) {
129 FATAL(error);
130 }
131
132 unique_ptr<JDynamics> dynamics;
133
134 if (!calibrationFile.empty()) {
135
136 try {
137
138 dynamics.reset(new JDynamics(detector, Tmax_s));
139
140 dynamics->load(calibrationFile);
141 }
142 catch(const exception& error) {
143 FATAL(error.what());
144 }
145 }
146
147 const JRouter_t* rpm = NULL;
148
149 if (target == module_t)
150 rpm = new JModuleRouter_t(detector);
151 else if (target == string_t)
152 rpm = new JStringRouter_t(detector);
153 else
154 FATAL("No valid target; check option -R" << endl);
155
156 const JModuleRouter router(dynamics ? dynamics->getDetector() : detector);
157
158 JSummaryFileRouter summary(inputFile);
159
160 const JMuonGandalf::storage_type storage(pdfFile, JTimeRange(parameters.TMin_ns, parameters.TMax_ns), parameters.TTS_ns);
161
162 JMuonGandalf fit(parameters, storage, debug);
163
164 typedef vector<JHitL0> JDataL0_t;
165 typedef vector<JHitW0> JDataW0_t;
166
167 const JBuildL0<JHitL0> buildL0;
168
169
170 typedef JManager<string, TProfile> JManager_t;
171
172 TH2D ha("ha", NULL, 256, -0.5, 255.5,
173 calibrate.getNumberOfBins(), calibrate.getLowerLimit(), calibrate.getUpperLimit());
174
175 TProfile hb("hb", NULL, 256, -0.5, 255.5);
176
177 TProfile hr("hr", NULL, 60, 0.0, 150.0);
178
179 TH2D h2("h2", NULL, rpm->getN(), -0.5, rpm->getN() - 0.5,
180 calibrate.getNumberOfBins(), calibrate.getLowerLimit(), calibrate.getUpperLimit());
181
182 JManager_t g1(new TProfile("%", NULL, rpm->getN(), -0.5, rpm->getN() - 0.5, -1.0, +1.0));
183
184
185 while (inputFile.hasNext()) {
186
187 STATUS("event: " << setw(10) << inputFile.getCounter() << '\r'); DEBUG(endl);
188
189 multi_pointer_type ps = inputFile.next();
190
191 JDAQEvent* tev = ps;
192 JFIT::JEvt* in = ps;
193
194 summary.update(*tev);
195
196 if (dynamics) {
197 dynamics->update(*tev);
198 }
199
200 in->select(parameters.numberOfPrefits, qualitySorter);
201
202 JDataL0_t dataL0;
203
204 buildL0(*tev, router, true, back_inserter(dataL0));
205
206 for (JFIT::JEvt::const_iterator track = in->begin(); track != in->end(); ++track) {
207
208 const JRotation3D R (getDirection(*track));
209 const JLine1Z tz(getPosition (*track).rotate(R), track->getT());
210 JRange<double> Z_m;
211
212 if (track->hasW(JSTART_LENGTH_METRES)) {
213 Z_m = JZRange(0.0, track->getW(JSTART_LENGTH_METRES)) + JRange<double>(parameters.ZMin_m, parameters.ZMax_m);
214 }
215
216 const JFIT::JModel<JLine1Z> match(tz, parameters.roadWidth_m, JRange<double>(parameters.TMin_ns, parameters.TMax_ns), Z_m);
217
218 // hit selection based on start value
219
220 JDataW0_t data;
221
222 for (JDataL0_t::const_iterator i = dataL0.begin(); i != dataL0.end(); ++i) {
223
224 JHitW0 hit(*i, summary.getRate(*i, parameters.R_Hz));
225
226 hit.rotate(R);
227
228 if (match(hit)) {
229 data.push_back(hit);
230 }
231 }
232
233 // select first hit
234
235 sort(data.begin(), data.end(), JHitL0::compare);
236
237 JDataW0_t::iterator __end = unique(data.begin(), data.end(), equal_to<JDAQPMTIdentifier>());
238
239
240 // set fit parameters
241
242 if (track->getE() > 0.1)
243 fit.JRegressor_t::E_GeV = track->getE();
244 else
245 fit.JRegressor_t::E_GeV = parameters.E_GeV;
246
247 set<int> buffer;
248
249 for (JDataW0_t::const_iterator hit = data.begin(); hit != __end; ++hit) {
250 buffer.insert(rpm->getIndex(hit->getModuleID()));
251 }
252
253 const JLine3Z ta(tz);
254
255 for (const int index : buffer) {
256
257 JDataW0_t::iterator q = partition(data.begin(), __end, [rpm, index](const JHitW0& hit) { return rpm->getIndex(hit.getModuleID()) != index; });
258
259 if (distance(data.begin(), q) - fit.parameters.size() > 0) {
260
261 fit(ta, data.begin(), q); // re-fit
262
263 for (JDataW0_t::const_iterator hit = q; hit != __end; ++hit) {
264
265 const double t1 = fit.value.getT(hit->getPosition());
266 JTrack3D gradient = fit(fit.value, *hit).gradient;
267
268 gradient.rotate_back(R);
269
270 ha.Fill(hit->getToT(), hit->getT1() - t1);
271 hb.Fill(hit->getToT(), gradient.getT());
272
273 hr.Fill(fit.value.getDistance(*hit), gradient.getT());
274
275 h2.Fill(index, hit->getT() - t1);
276
277 g1["T"]->Fill(index, gradient.getT());
278 g1["X"]->Fill(index, gradient.getX());
279 g1["Y"]->Fill(index, gradient.getY());
280 g1["Z"]->Fill(index, gradient.getZ());
281 }
282 }
283 }
284 }
285 }
286 STATUS(endl);
287
288 TFile out(outputFile.c_str(), "recreate");
289
290 putObject(out, JMeta(argc, argv));
291
292 out << ha;
293 out << hb;
294 out << hr;
295 out << h2;
296 out << g1;
297
298 out.Write();
299 out.Close();
300
301 delete rpm;
302}
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:72
#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.
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.
Range of values.
Definition JRange.hh:42
Template L0 hit builder.
Definition JBuildL0.hh:38
static const int JSTART_LENGTH_METRES
distance between projected positions on the track of optical modules for which the response does not ...
JDirection3D getDirection(const Vec &dir)
Get direction.
JPosition3D getPosition(const Vec &pos)
Get position.
static const std::string string_t
string
Definition JSydney.cc:90
JTOOLS::JRange< double > JTimeRange
Type definition for time range (unit [s]).
const char *const module_t
routing by module
void load(const std::string &file_name, JDetector &detector)
Load detector from input file.
JTOOLS::JRange< double > JZRange
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
Detector file.
Definition JHead.hh:227
Target.
Definition JHead.hh:300
Acoustic event fit.
Interface for routing module identifier to index and vice versa.
virtual int getIndex(const int id) const =0
Get index.
virtual size_t getN() const =0
Get number of indices.
Orientation of module.
Dynamic detector calibration.
Definition JDynamics.hh:86
Auxiliary class to match data points with given model.
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.
Auxiliary class for defining the range of iterations of objects.
Definition JLimit.hh:45
const JLimit & getLimit() const
Get limit.
Definition JLimit.hh:84
static counter_type max()
Get maximum counter value.
Definition JLimit.hh:128
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