Jpp
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Functions
JPulsar.cc File Reference

Application for dark room time calibration. 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 "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 "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/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.

Author
mbouwhuis

Definition in file JPulsar.cc.

Function Documentation

int main ( int  argc,
char **  argv 
)

Definition at line 97 of file JPulsar.cc.

98 {
99  using namespace std;
100  using namespace JPP;
101  using namespace KM3NETDAQ;
102 
103  typedef JRange<double> JRange_t;
104 
106  JLimit_t& numberOfEvents = inputFile.getLimit();
107  string outputFile;
108  string detectorFile;
109  JTDC_t TDC;
110  double laserFrequency_Hz;
111  bool overwriteDetector;
112  JROOTClassSelector selector;
113  string option;
114  JRange_t T_ns;
115  int debug;
116 
117  try {
118 
119  JParser<> zap("Application for dark room time calibration.");
120 
121  zap['f'] = make_field(inputFile, "input file (time slice data from laser calibration).");
122  zap['o'] = make_field(outputFile, "output file.") = "pulsar.root";
123  zap['a'] = make_field(detectorFile, "detector file.");
124  zap['!'] = make_field(TDC,
125  "Set reference PMTs, e.g."
126  "\n-! \"808969848 0 808982077 23\" sets PMT 0 of module 808969848 and PMT 23 of module 808982077 as references.");
127  zap['n'] = make_field(numberOfEvents) = JLimit::max();
128  zap['l'] = make_field(laserFrequency_Hz, "laser frequency [Hz]") = 10000;
129  zap['A'] = make_field(overwriteDetector, "overwrite detector file provided through '-a' with correct time offsets.");
130  zap['O'] = make_field(option, "ROOT fit option, see TH1::Fit.") = "LS";
131  zap['C'] = make_field(selector, "timeslice selector, e.g. JDAQTimesliceL1.") = getROOTClassSelection<JDAQTimesliceTypes_t>();
132  zap['T'] = make_field(T_ns) = JRange_t(-10.0, +10.0);
133  zap['d'] = make_field(debug, "debug.") = 1;
134 
135  zap(argc, argv);
136  }
137  catch(const exception& error) {
138  FATAL(error.what() << endl);
139  }
140 
141 
142  for (JTDC_t::iterator tdc = TDC.begin(); tdc != TDC.end(); ) {
143 
144  if (tdc->second == WILD_CARD) {
145 
146  const int id = tdc->first;
147 
148  TDC.erase(tdc);
149 
150  for (int pmt = 0; pmt != NUMBER_OF_PMTS; ++pmt) {
151  tdc = TDC.insert(make_pair(id,pmt));
152  }
153 
154  } else {
155 
156  ++tdc;
157  }
158  }
159 
160  for (JTDC_t::const_iterator tdc = TDC.begin(); tdc != TDC.end(); ++tdc) {
161  if (tdc->second < 0 || tdc->second >= NUMBER_OF_PMTS) {
162  FATAL("Illegal TDC (= readout channel) identifier: " << tdc->second << " {0, .., " << NUMBER_OF_PMTS - 1 << "}" << endl);
163  }
164  }
165 
166  if (laserFrequency_Hz <= 0.0) {
167  FATAL("Invalid laser frequency " << laserFrequency_Hz << endl);
168  }
169 
170 
171  const double laserPeriod_ns = 1.0e9 / laserFrequency_Hz;
172 
173  cout.tie(&cerr);
174 
176 
177  try {
178  load(detectorFile, detector);
179  }
180  catch(const JException& error) {
181  FATAL(error);
182  }
183 
184  const JModuleRouter moduleRouter(detector);
185 
186 
187  // ROOT I/O
188 
189  TFile out(outputFile.c_str(), "recreate");
190 
191  putObject(&out, JMeta(argc, argv));
192 
193  typedef map<JPMTIdentifier, TH1D*> map_type;
194 
195  map_type zmap;
196 
197  const int nx = 2 * (int) laserPeriod_ns;
198  const double xmin = -0.5 * laserPeriod_ns;
199  const double xmax = 0.5 * laserPeriod_ns;
200 
201  TH1D h0("h0", NULL, nx, xmin, xmax);
202  TH1D h1("h1", NULL, 256, -0.5, +255.5);
203 
205 
206  for (JDetector::iterator module = detector.begin(); module != detector.end(); ++module) {
207 
208  const range_type range = get_range(TDC, module->getID());
209 
210  for (JTDC_t::const_iterator i = range.first; i != range.second; ++i) {
211 
212  NOTICE("Reference PMT at"
213  << " string " << setw(3) << module->getString()
214  << " floor " << setw(2) << module->getFloor()
215  << " module " << setw(8) << module->getID()
216  << " channel " << setw(2) << i->second << endl);
217 
218  const JPMTIdentifier id(module->getID(),i->second);
219 
220  ostringstream os;
221 
222  os << getLabel(module->getLocation()) << ' ' << getLabel(id);
223 
224  zmap.insert(make_pair(id, new TH1D(os.str().c_str(), NULL, nx, xmin, xmax)));
225  }
226  }
227 
228 
229  typedef vector<JHitR0> JDataL0_t;
230 
231  const JBuildL0<JHitR0> buildL0;
232 
233 
235 
236  counter_type counter = 0;
237 
238  for ( ; in.hasNext() && counter != inputFile.getLimit(); ++counter) {
239 
240  STATUS("event: " << setw(10) << counter << '\r'); DEBUG(endl);
241 
242  const JDAQTimeslice* timeslice = in.next();
243 
244  for (JDAQTimeslice::const_iterator frame = timeslice->begin(); frame != timeslice->end(); ++frame) {
245 
246  const range_type range = get_range(TDC, frame->getModuleID());
247 
248  if (range.first != range.second) {
249 
250  const double t0 = get_time(getTimeOfFrame(frame->getFrameIndex()), laserPeriod_ns);
251 
252  JDataL0_t dataL0;
253 
254  buildL0(*frame, moduleRouter.getModule(frame->getModuleID()), back_inserter(dataL0));
255 
256  for (JDataL0_t::const_iterator hit = dataL0.begin(); hit != dataL0.end(); ++hit) {
257 
258  const JPMTIdentifier id(frame->getModuleID(), hit->getPMT());
259 
260  map_type::const_iterator p = zmap.find(id);
261 
262  if (p != zmap.end()) {
263 
264  const double t1 = get_time(t0 + hit->getT(), laserPeriod_ns);
265 
266  p->second->Fill(t1);
267 
268  h0.Fill(t1);
269 
270  if (T_ns(t1)) {
271 
272  h1.Fill(hit->getToT());
273 
274  counts[id] += 1;
275  }
276  }
277  }
278  }
279  }
280  }
281  STATUS(endl);
282 
283 
284  map<JPMTIdentifier, double> t0; // fitted time offsets
285 
286 
287  TF1 f1("f1", "[0]*exp(-0.5*(x-[1])*(x-[1])/([2]*[2])) + [3]");
288 
289  for (map_type::iterator it = zmap.begin(); it != zmap.end(); ++it) {
290 
291  JPMTIdentifier pmt = it->first;
292  TH1D* h1 = it->second;
293 
294  if (h1->GetEntries() == 0) {
295  WARNING("Histogram " << h1->GetName() << " empty" << endl);
296  continue;
297  }
298 
299  STATUS("--- PMT = " << pmt << "; histogram " << h1->GetName() << endl);
300 
301  // start values
302 
303  Double_t ymax = 0.0;
304  Double_t x0 = 0.0; // [ns]
305  Double_t sigma = 2.0; // [ns]
306  Double_t ymin = 0.0;
307 
308  for (int i = 1; i != h1->GetNbinsX(); ++i) {
309 
310  const Double_t x = h1->GetBinCenter (i);
311  const Double_t y = h1->GetBinContent(i);
312 
313  if (y > ymax) {
314  ymax = y;
315  x0 = x;
316  }
317  }
318 
319  f1.SetParameter(0, ymax);
320  f1.SetParameter(1, x0);
321  f1.SetParameter(2, sigma);
322  f1.SetParameter(3, ymin);
323 
324 
325  // fit
326 
327  TFitResultPtr result = h1->Fit(&f1, option.c_str(), "same", x0 - 3 * sigma, x0 + 3 * sigma);
328 
329 
330  if (debug >= notice_t || !result->IsValid()) {
331  cout << "Histogram " << h1->GetName() << " fit " << (result->IsValid() ? "" : "failed") << endl;
332  cout << "\tx0 = " << FIXED(12,3) << x0 << endl;
333  cout << "\tt0 = " << FIXED(12,3) << f1.GetParameter(1) << endl;
334  }
335 
336  t0[pmt] = f1.GetParameter(1);
337  }
338 
339 
340  if (counter != 0) {
341 
342  const double W = laserFrequency_Hz * counter * getFrameTime() * 1.0e-9;
343 
344  NOTICE("Expection values [npe]" << endl);
345 
346  for (map<JPMTIdentifier, int>::const_iterator i = counts.begin(); i != counts.end(); ++i) {
347  NOTICE(i->first << ' ' << FIXED(7,3) << i->second / W << endl);
348  }
349  }
350 
351 
352  if (overwriteDetector) {
353 
354  NOTICE("Store calibration data on file " << detectorFile << endl);
355 
356  detector.comment.add(JMeta(argc, argv));
357 
358  for (JDetector::iterator module = detector.begin(); module != detector.end(); ++module) {
359 
360  const range_type range = get_range(TDC, module->getID());
361 
362  if (range.first != range.second) {
363 
364  const double t1 = t0[JPMTIdentifier(range.first->first,
365  range.first->second)]; // t0 of the first reference PMT in this optical module
366 
367  for (int pmt = 0; pmt != NUMBER_OF_PMTS; ++pmt) {
368 
369  map<JPMTIdentifier, double>::const_iterator p = t0.find(JPMTIdentifier(module->getID(),pmt));
370 
371  if (p != t0.end())
372  module->getPMT(pmt).subT0(p->second); // offset with t0 of this reference PMT
373  else
374  module->getPMT(pmt).subT0(t1); // offset with t0 of first reference PMT in same optical module
375  }
376  }
377  }
378 
379  store(detectorFile, detector);
380  }
381 
382  out.Write();
383  out.Close();
384 }
Auxiliary class for ROOT I/O of application specific meta data.
Definition: JMeta.hh:71
Utility class to parse command line options.
Definition: JParser.hh:1493
General exception.
Definition: JException.hh:23
#define WARNING(A)
Definition: JMessage.hh:65
std::string getLabel(const JLocation &location)
Get module label for monitoring and other applications.
Definition: JLocation.hh:246
#define STATUS(A)
Definition: JMessage.hh:63
Detector data structure.
Definition: JDetector.hh:80
Auxiliary class to select ROOT class based on class name.
Router for direct addressing of module data in detector data structure.
esac print_variable DETECTOR INPUT_FILE OUTPUT_FILE CDF for TYPE in
Definition: JSirene.sh:45
notice
Definition: JMessage.hh:32
then for HISTOGRAM in h0 h1
Definition: JMatrixNZ.sh:69
Long64_t counter_type
Type definition for counter.
Auxiliary data structure for floating point format specification.
Definition: JPrint.hh:481
Auxiliary class for multiplexing object iterators.
string outputFile
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
esac $JPP_DIR examples JDetector JTransitTime o $OUTPUT_FILE n N $NPE T $TTS_NS d $DEBUG for HISTOGRAM in tts tt2 pmt
Definition: JTransitTime.sh:36
Auxiliary class for defining the range of iterations of objects.
Definition: JLimit.hh:41
JKey_t first
Definition: JPair.hh:128
Template specialisation of L0 builder for JHitR0 data type.
Definition: JBuildL0.hh:174
Detector file.
Definition: JHead.hh:130
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1954
double getFrameTime()
Get frame time duration.
Definition: JDAQClock.hh:162
return result
Definition: JPolint.hh:695
#define NOTICE(A)
Definition: JMessage.hh:64
Data time slice.
void load(const JString &file_name, JDetector &detector)
Load detector from input file.
int debug
debug level
Definition: JSirene.cc:61
#define FATAL(A)
Definition: JMessage.hh:67
JRange< Double_t > JRange_t
Definition: JFitToT.hh:34
General purpose class for object reading from a list of file names.
static const char WILD_CARD
Definition: JDAQTags.hh:34
void store(const JString &file_name, const JDetector &detector)
Store detector to output file.
const JLimit & getLimit() const
Get limit.
Definition: JLimit.hh:73
static const int NUMBER_OF_PMTS
Total number of PMTs in module.
Definition: JDAQ.hh:26
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:62
bool putObject(TDirectory *dir, const T &object)
Write object to ROOT directory.