87{
91
93
95 JLimit_t& numberOfEvents = inputFile.getLimit();
97 string detectorFile;
99 double laserFrequency_Hz;
100 bool overwriteDetector;
102 string option;
103 double Wmin = 1e3;
107
108 try {
109
111
113
114 JParser<> zap(
"Application for dark room time calibration.");
115
116 zap[
'f'] =
make_field(inputFile,
"input file (time slice data from laser calibration).");
118 zap[
'a'] =
make_field(detectorFile,
"detector file.");
120 "Set reference PMTs, e.g."
121 "\n-! \"808969848 0 808982077 23\" sets PMT 0 of module 808969848 and PMT 23 of module 808982077 as references.");
123 zap[
'l'] =
make_field(laserFrequency_Hz,
"laser frequency [Hz]") = 10000;
124 zap[
'A'] =
make_field(overwriteDetector,
"overwrite detector file provided through '-a' with correct time offsets.");
125 zap[
'O'] =
make_field(option,
"ROOT fit option, see TH1::Fit.") =
"LS";
128 zap[
'T'] =
make_field(T_ns,
"time window for time-over-threshold monitor") =
JRange_t(-10.0, +10.0);
131
132 zap(argc, argv);
133 }
134 catch(const exception& error) {
135 FATAL(error.what() << endl);
136 }
137
138
140
141 if (laserFrequency_Hz <= 0.0) {
142 FATAL(
"Invalid laser frequency " << laserFrequency_Hz << endl);
143 }
144
145 const double laserPeriod_ns = 1.0e9 / laserFrequency_Hz;
146
147 if (option.find('R') == string::npos) { option += 'R'; }
148 if (option.find('S') == string::npos) { option += 'S'; }
150
151
153
154 try {
156 }
159 }
160
162
163
164
165
167
169
171
172 map_type zmap;
173
176 const int nx = 2 * (int) (xmax - xmin);
177
178 TH1D h0("h0", NULL, nx, xmin, xmax);
179 TH1D h1("h1", NULL, 256, -0.5, +255.5);
180
182
183 for (JDetector::iterator module =
detector.begin(); module !=
detector.end(); ++module) {
184
185 if (!module->empty()) {
186
188
189 for (JTDC_t::const_iterator i = range.first; i != range.second; ++i) {
190
192 << " string " << setw(3) << module->getString()
193 << " floor " << setw(2) << module->getFloor()
194 << " module " << setw(8) << module->getID()
195 << " channel " << setw(2) << i->second << endl);
196
198
199 ostringstream os;
200
202
203 zmap.insert(make_pair(id, new TH1D(os.str().c_str(), NULL, nx, xmin, xmax)));
204 }
205 }
206 }
207
208
210
212
213
215
217
218 for ( ; in.hasNext() && counter != inputFile.getLimit(); ++counter) {
219
220 STATUS(
"event: " << setw(10) << counter <<
'\r');
DEBUG(endl);
221
223
224 for (JDAQTimeslice::const_iterator frame = timeslice->begin(); frame != timeslice->end(); ++frame) {
225
227
228 if (range.first != range.second) {
229
230 const double t0 = get_time(
getTimeOfFrame(frame->getFrameIndex()), laserPeriod_ns);
231
232 JDataL0_t dataL0;
233
234 buildL0(*frame, moduleRouter.getModule(frame->getModuleID()), back_inserter(dataL0));
235
236 for (JDataL0_t::const_iterator hit = dataL0.begin(); hit != dataL0.end(); ++hit) {
237
239
240 map_type::const_iterator p = zmap.find(id);
241
242 if (p != zmap.end()) {
243
244 const double t1 = get_time(t0 + hit->getT(), laserPeriod_ns);
245
246 p->second->Fill(t1);
247
248 h0.Fill(t1);
249
250 if (T_ns(t1)) {
251
252 h1.Fill(hit->getToT());
253
254 counts[id] += 1;
255 }
256 }
257 }
258 }
259 }
260 }
262
263
265
266
267 TF1 f1("f1", "[0]*exp(-0.5*(x-[1])*(x-[1])/([2]*[2]))/(TMath::Sqrt(2.0*TMath::Pi())*[2]) + [3]");
268
269 for (map_type::iterator it = zmap.begin(); it != zmap.end(); ++it) {
270
272 TH1D* h1 = it->second;
273
274 if (h1->GetEntries() == 0) {
275 WARNING(
"Histogram " << h1->GetName() <<
" empty" << endl);
276 continue;
277 }
278
279 STATUS(
"--- PMT = " << pmt <<
"; histogram " << h1->GetName() << endl);
280
281
282
283 Double_t ymax = 0.0;
284 Double_t x0 = 0.0;
285 Double_t sigma = 2.0;
286 Double_t ymin = 0.0;
287
288 for (int i = 1; i != h1->GetNbinsX(); ++i) {
289
290 const Double_t
x = h1->GetBinCenter (i);
291 const Double_t
y = h1->GetBinContent(i);
292
293 if (y > ymax) {
296 }
297 }
298
299 f1.SetParameter(0, ymax/sigma);
300 f1.SetParameter(1, x0);
301 f1.SetParameter(2, sigma);
302 f1.SetParameter(3, ymin);
303
304 for (Int_t i = 0; i != f1.GetNpar(); ++i) {
305 f1.SetParError(i, 0.0);
306 }
307
308
309
310
311 TFitResultPtr
result = h1->Fit(&f1, option.c_str(),
"same", x0 - 3 * sigma, x0 + 3 * sigma);
312
313 if (
result.Get() == NULL) {
314 FATAL(
"Invalid TFitResultPtr " << h1->GetName() << endl);
315 }
316
318 cout <<
"Histogram " << h1->GetName() <<
" fit " << (
result->IsValid() ?
"ok" :
"failed") << endl;
319 cout <<
"\tw = " <<
FIXED(12,3) << f1.GetParameter(0) << endl;
320 cout <<
"\tx0 = " <<
FIXED(12,3) << x0 << endl;
321 cout <<
"\tt0 = " <<
FIXED(12,3) << f1.GetParameter(1) << endl;
322 }
323
324
325
326
327 int number_of_peaks = 0;
328
329 Double_t dx = 2.0 * f1.GetParameter(2);
330 Double_t W = 0.0;
331 Double_t Y = f1.GetParameter(3);
332
333 if (dx < (xmax - xmin) / nx) {
334 dx = (xmax - xmin) / nx;
335 }
336
337 for (Int_t il = 1, ir = il; ir <= nx; ) {
338
339 for ( ; ir <= nx && h1->GetBinCenter(ir) <= h1->GetBinCenter(il) + dx; ++ir) {
340 W += h1->GetBinContent(ir) - Y;
341 }
342
343 if (W >= Wmin) {
344
345 number_of_peaks += 1;
346
347 W = 0.0;
348 il = ir;
349 ir = il;
350
351 } else {
352
353 W -= h1->GetBinContent(il) - Y;
354 il += 1;
355 }
356 }
357
358 if (number_of_peaks != 1) {
359 WARNING(
"Number of peaks " << h1->GetName() <<
' ' << number_of_peaks <<
" != 1" << endl);
360 }
361
362 if (
result->IsValid() && f1.GetParameter(0) >= Wmin) {
363 t0[pmt] = f1.GetParameter(1);
364 }
365 }
366
367 out.Write();
368 out.Close();
369
370
371 if (counter != 0) {
372
373 const double W = laserFrequency_Hz * counter *
getFrameTime() * 1.0e-9;
374
375 NOTICE(
"Expection values [npe]" << endl);
376
378 NOTICE(i->first <<
' ' <<
FIXED(7,3) << i->second / W << endl);
379 }
380 }
381
382
383 if (overwriteDetector) {
384
385 int errors = 0;
386
387 for (JDetector::iterator module =
detector.begin(); module !=
detector.end(); ++module) {
388
389 if (!module->empty()) {
390
392
393 if (range.first != range.second) {
394
396
398
399 if (p0 != t0.end()) {
400
401 const double t1 = p0->second;
402
404
406
408 module->getPMT(pmt).subT0(p1->second);
409 else
410 module->getPMT(pmt).subT0(t1);
411 }
412
413 } else {
414
415 if (!module->getPMT(
id.getPMTAddress()).has(
PMT_DISABLE)) {
416 ++errors;
417 }
418
420 << setw(10) << module->getID() <<
"/" <<
FILL(2,
'0') <<
id.getPMTAddress() <<
FILL() <<
' '
421 << "missing or insufficient signal." << endl);
422 }
423 }
424 }
425 }
426
427 if (errors == 0) {
428
429 NOTICE(
"Store calibration data on file " << detectorFile << endl);
430
432
434
435 } else {
436
437 FATAL(
"Number of errors " << errors <<
" != 0" << endl);
438 }
439 }
440}
#define DEBUG(A)
Message macros.
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
#define gmake_property(A)
macros to convert (template) parameter to JPropertiesElement object
Router for direct addressing of module data in detector data structure.
Utility class to parse parameter values.
Auxiliary class for multiplexing object iterators.
Utility class to parse command line options.
General purpose class for object reading from a list of file names.
std::string getLabel(const JLocation &location)
Get module label for monitoring and other applications.
void load(const std::string &file_name, JDetector &detector)
Load detector from input file.
void store(const std::string &file_name, const JDetector &detector)
Store detector to output file.
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
bool putObject(TDirectory &dir, const TObject &object)
Write object to ROOT directory.
std::set< JROOTClassSelector > getROOTClassSelection(const bool option=false)
Get ROOT class selection.
Long64_t counter_type
Type definition for counter.
KM3NeT DAQ data structures and auxiliaries.
double getFrameTime()
Get frame time duration.
double getTimeOfFrame(const int frame_index)
Get start time of frame in ns since start of run for a given frame index.
static const int NUMBER_OF_PMTS
Total number of PMTs in module.
std::map< int, range_type > map_type
static const int PMT_DISABLE
ile KM3NeT Data Definitions v3.6.2 https://git.km3net.de/common/km3net-dataformat
Auxiliary data structure for sequence of same character.
Auxiliary data structure for floating point format specification.
Type definition of range.
Auxiliary class for TDC constraints.
range_type equal_range(const int id)
Get range of constraints for given module.
bool is_valid(const bool option=false) const
Check validity of TDC constrants.
Empty structure for specification of parser element that is initialised (i.e. does not require input)...
Auxiliary class to select ROOT class based on class name.
Auxiliary class for defining the range of iterations of objects.
static counter_type max()
Get maximum counter value.