Jpp 20.0.0-72-g597b30bc9
the software that should make you happy
Loading...
Searching...
No Matches
JPulsar.cc File Reference

Application for dark room time calibration using external laser. More...

#include <string>
#include <iostream>
#include <iomanip>
#include <map>
#include "TROOT.h"
#include "TFile.h"
#include "TH1D.h"
#include "TF1.h"
#include "TFitResult.h"
#include "JROOT/JMinimizer.hh"
#include "JDAQ/JDAQTimesliceIO.hh"
#include "JDetector/JDetector.hh"
#include "JDetector/JDetectorToolkit.hh"
#include "JDetector/JModuleRouter.hh"
#include "JDetector/JPMTIdentifier.hh"
#include "JTrigger/JHit.hh"
#include "JTrigger/JFrame.hh"
#include "JTrigger/JTimeslice.hh"
#include "JTrigger/JHitR0.hh"
#include "JTrigger/JBuildL0.hh"
#include "JCalibrate/JTDC_t.hh"
#include "JSupport/JMultipleFileScanner.hh"
#include "JSupport/JSupport.hh"
#include "JSupport/JMeta.hh"
#include "JTools/JRange.hh"
#include "JLang/JObjectMultiplexer.hh"
#include "JROOT/JROOTClassSelector.hh"
#include "Jeep/JProperties.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

Application for dark room time calibration using external laser.

The time calibration is made for the specified reference PMTs (option -! "\<module identifier\> \<PMT channel\>".
In this, the wild card value -1 can be used. A ROOT histogram corresponding to each reference PMT is created, filled, fitted and stored to the specified output file.
The detector can be updated using option -A
The laser should be synchronised with the same clock system as the detector and produce single pulses at a fixed frequency (option -l <laserFrequency_Hz>.

In addition, the transit-time distribution of the reference PMTs can be monitored (see JLegolas.cc).
For this, the expectation values should be less than one and the time-over-threshold distribution nominal.

Author
mbouwhuis

Definition in file JPulsar.cc.

Function Documentation

◆ main()

int main ( int argc,
char ** argv )

Definition at line 86 of file JPulsar.cc.

87{
88 using namespace std;
89 using namespace JPP;
90 using namespace KM3NETDAQ;
91
93
95 JLimit_t& numberOfEvents = inputFile.getLimit();
96 string outputFile;
97 string detectorFile;
98 JTDC_t TDC;
99 double laserFrequency_Hz;
100 bool overwriteDetector;
101 JROOTClassSelector selector;
102 string option;
103 double Wmin = 1e3; // Minimal signal intensity
104 JRange_t T_ns;
105 JRange_t X;
106 int debug;
107
108 try {
109
110 JProperties properties;
111
112 properties.insert(gmake_property(Wmin));
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).");
117 zap['o'] = make_field(outputFile, "output file.") = "pulsar.root";
118 zap['a'] = make_field(detectorFile, "detector file.");
119 zap['!'] = make_field(TDC,
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.");
122 zap['n'] = make_field(numberOfEvents) = JLimit::max();
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";
126 zap['@'] = make_field(properties) = JPARSER::initialised();
127 zap['C'] = make_field(selector, "timeslice selector, e.g. JDAQTimesliceL1.") = getROOTClassSelection<JDAQTimesliceTypes_t>();
128 zap['T'] = make_field(T_ns, "time window for time-over-threshold monitor") = JRange_t(-10.0, +10.0);
129 zap['x'] = make_field(X, "time window for histogram") = JRange_t();
130 zap['d'] = make_field(debug, "debug.") = 1;
131
132 zap(argc, argv);
133 }
134 catch(const exception& error) {
135 FATAL(error.what() << endl);
136 }
137
138
139 TDC.is_valid(true);
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'; }
149 if (option.find('Q') == string::npos && debug < JEEP::debug_t) { option += 'Q'; }
150
151
153
154 try {
155 load(detectorFile, detector);
156 }
157 catch(const JException& error) {
158 FATAL(error);
159 }
160
161 const JModuleRouter moduleRouter(detector);
162
163
164 // ROOT I/O
165
166 TFile out(outputFile.c_str(), "recreate");
167
168 putObject(out, JMeta(argc, argv));
169
171
172 map_type zmap;
173
174 const double xmin = (X == JRange_t() ? -0.5 * laserPeriod_ns : X.getLowerLimit());
175 const double xmax = (X == JRange_t() ? +0.5 * laserPeriod_ns : X.getUpperLimit());
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
187 const JTDC_t::range_type range = TDC.equal_range(module->getID());
188
189 for (JTDC_t::const_iterator i = range.first; i != range.second; ++i) {
190
191 NOTICE("Reference PMT at"
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
197 const JPMTIdentifier id(module->getID(),i->second);
198
199 ostringstream os;
200
201 os << getLabel(module->getLocation()) << ' ' << getLabel(id);
202
203 zmap.insert(make_pair(id, new TH1D(os.str().c_str(), NULL, nx, xmin, xmax)));
204 }
205 }
206 }
207
208
209 typedef vector<JHitR0> JDataL0_t;
210
211 const JBuildL0<JHitR0> buildL0;
212
213
215
216 counter_type counter = 0;
217
218 for ( ; in.hasNext() && counter != inputFile.getLimit(); ++counter) {
219
220 STATUS("event: " << setw(10) << counter << '\r'); DEBUG(endl);
221
222 const JDAQTimeslice* timeslice = in.next();
223
224 for (JDAQTimeslice::const_iterator frame = timeslice->begin(); frame != timeslice->end(); ++frame) {
225
226 const JTDC_t::range_type range = TDC.equal_range(frame->getModuleID());
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
238 const JPMTIdentifier id(frame->getModuleID(), hit->getPMT());
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 }
261 STATUS(endl);
262
263
264 map<JPMTIdentifier, double> t0; // fitted time offsets
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
271 JPMTIdentifier pmt = it->first;
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 // start values
282
283 Double_t ymax = 0.0;
284 Double_t x0 = 0.0; // [ns]
285 Double_t sigma = 2.0; // [ns]
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) {
294 ymax = y;
295 x0 = x;
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 // fit
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
317 if (debug >= notice_t || !result->IsValid() || f1.GetParameter(0) < Wmin) {
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 // check for multiple peaks
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; // minimal width
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; // reset
348 il = ir;
349 ir = il;
350
351 } else {
352
353 W -= h1->GetBinContent(il) - Y; // shift
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
377 for (map<JPMTIdentifier, int>::const_iterator i = counts.begin(); i != counts.end(); ++i) {
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
391 const JTDC_t::range_type range = TDC.equal_range(module->getID());
392
393 if (range.first != range.second) {
394
395 const JPMTIdentifier id(module->getID(), range.first->second); // select t0 of the first reference PMT in this optical module
396
398
399 if (p0 != t0.end()) {
400
401 const double t1 = p0->second;
402
403 for (int pmt = 0; pmt != NUMBER_OF_PMTS; ++pmt) {
404
405 map<JPMTIdentifier, double>::const_iterator p1 = t0.find(JPMTIdentifier(module->getID(),pmt));
406
407 if (p1 != t0.end())
408 module->getPMT(pmt).subT0(p1->second); // offset with t0 of this reference PMT
409 else
410 module->getPMT(pmt).subT0(t1); // offset with t0 of first reference PMT in same optical module
411 }
412
413 } else {
414
415 if (!module->getPMT(id.getPMTAddress()).has(PMT_DISABLE)) {
416 ++errors;
417 }
418
419 ERROR("Module/PMT "
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
431 detector.comment.add(JMeta(argc, argv));
432
433 store(detectorFile, detector);
434
435 } else {
436
437 FATAL("Number of errors " << errors << " != 0" << endl);
438 }
439 }
440}
string outputFile
TPaveText * p1
#define DEBUG(A)
Message macros.
Definition JMessage.hh:62
#define STATUS(A)
Definition JMessage.hh:63
#define NOTICE(A)
Definition JMessage.hh:64
#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
#define gmake_property(A)
macros to convert (template) parameter to JPropertiesElement object
Detector data structure.
Definition JDetector.hh:96
Router for direct addressing of module data in detector data structure.
Utility class to parse parameter values.
General exception.
Definition JException.hh:24
Auxiliary class for multiplexing object iterators.
Utility class to parse command line options.
Definition JParser.hh:1698
General purpose class for object reading from a list of file names.
Range of values.
Definition JRange.hh:42
T getLowerLimit() const
Get lower limit.
Definition JRange.hh:202
T getUpperLimit() const
Get upper limit.
Definition JRange.hh:213
Template L0 hit builder.
Definition JBuildL0.hh:38
std::string getLabel(const JLocation &location)
Get module label for monitoring and other applications.
Definition JLocation.hh:247
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.
@ debug_t
debug
Definition JMessage.hh:29
@ notice_t
notice
Definition JMessage.hh:32
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.
return result
Definition JPolint.hh:862
KM3NeT DAQ data structures and auxiliaries.
Definition DataQueue.cc:39
double getFrameTime()
Get frame time duration.
Definition JDAQClock.hh:162
double getTimeOfFrame(const int frame_index)
Get start time of frame in ns since start of run for a given frame index.
Definition JDAQClock.hh:185
static const int NUMBER_OF_PMTS
Total number of PMTs in module.
Definition JDAQ.hh:26
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
Definition pmt_status.hh:13
Auxiliary data structure for sequence of same character.
Definition JManip.hh:330
Auxiliary data structure for floating point format specification.
Definition JManip.hh:448
Type definition of range.
Definition JHead.hh:43
Detector file.
Definition JHead.hh:227
Auxiliary class for TDC constraints.
Definition JTDC_t.hh:39
range_type equal_range(const int id)
Get range of constraints for given module.
Definition JTDC_t.hh:101
bool is_valid(const bool option=false) const
Check validity of TDC constrants.
Definition JTDC_t.hh:171
Empty structure for specification of parser element that is initialised (i.e. does not require input)...
Definition JParser.hh:68
Auxiliary class to select ROOT class based on class name.
Auxiliary class for defining the range of iterations of objects.
Definition JLimit.hh:45
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