Jpp  18.0.0-rc.4
the software that should make you happy
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Functions
JFilter.cc File Reference

Example program to test performance of various hit filters based on JTRIGGER::clusterize method. More...

#include <string>
#include <iostream>
#include <iomanip>
#include <vector>
#include <set>
#include <map>
#include <algorithm>
#include <iterator>
#include "TROOT.h"
#include "TFile.h"
#include "TH1D.h"
#include "TProfile.h"
#include "km3net-dataformat/offline/Head.hh"
#include "km3net-dataformat/offline/Evt.hh"
#include "km3net-dataformat/offline/Hit.hh"
#include "km3net-dataformat/tools/time_converter.hh"
#include "JAAnet/JHead.hh"
#include "JAAnet/JHeadToolkit.hh"
#include "JAAnet/JAAnetToolkit.hh"
#include "JDAQ/JDAQEventIO.hh"
#include "JDAQ/JDAQTimesliceIO.hh"
#include "km3net-dataformat/online/JDAQClock.hh"
#include "JSirene/JPulse.hh"
#include "JPhysics/JConstants.hh"
#include "JDetector/JDetector.hh"
#include "JDetector/JDetectorToolkit.hh"
#include "JDetector/JModuleRouter.hh"
#include "JDetector/JPMTRouter.hh"
#include "JDetector/JTimeRange.hh"
#include "JTrigger/JHit.hh"
#include "JTrigger/JFrame.hh"
#include "JTrigger/JTimeslice.hh"
#include "JTrigger/JHitL0.hh"
#include "JTrigger/JHitL1.hh"
#include "JTrigger/JBuildL1.hh"
#include "JTrigger/JBuildL2.hh"
#include "JTrigger/JAlgorithm.hh"
#include "JTrigger/JMatch3B.hh"
#include "JTrigger/JMatch3D.hh"
#include "JTrigger/JMatch1D.hh"
#include "JTrigger/JBind2nd.hh"
#include "JSupport/JTriggeredFileScanner.hh"
#include "JSupport/JMonteCarloFileSupportkit.hh"
#include "JSupport/JSupport.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

Example program to test performance of various hit filters based on JTRIGGER::clusterize method.

Author
mdejong

Definition in file JFilter.cc.

Function Documentation

int main ( int  argc,
char **  argv 
)

Definition at line 212 of file JFilter.cc.

213 {
214  using namespace std;
215  using namespace JPP;
216  using namespace KM3NETDAQ;
217 
218  JTriggeredFileScanner<> inputFile;
219  JLimit_t& numberOfEvents = inputFile.getLimit();
220  string outputFile;
221  string detectorFile;
222  double Tmax_ns;
223  double roadWidth_m;
224  double ctMin;
225  double thetaRad;
226  char cluster;
227  int histogram;
228  int debug;
229 
230  try {
231 
232  JParser<> zap("Example program to test performance of various hit filters.");
233 
234  zap['f'] = make_field(inputFile);
235  zap['o'] = make_field(outputFile) = "filter.root";
236  zap['a'] = make_field(detectorFile);
237  zap['n'] = make_field(numberOfEvents) = JLimit::max();
238  zap['T'] = make_field(Tmax_ns) = 20.0; // [ns]
239  zap['R'] = make_field(roadWidth_m) = 150.0; // [m]
240  zap['C'] = make_field(ctMin) = 0.0; //
241  zap['Q'] = make_field(thetaRad) = 0.2; //
242  zap['c'] = make_field(cluster) = 'A', 'B', 'C', 'D', 'E', 'F';
243  zap['H'] = make_field(histogram) = 1, 2;
244  zap['d'] = make_field(debug) = 1;
245 
246  zap(argc, argv);
247  }
248  catch(const exception &error) {
249  FATAL(error.what() << endl);
250  }
251 
252 
253  using namespace KM3NETDAQ;
254 
255 
257 
258  try {
259  load(detectorFile, detector);
260  }
261  catch(const JException& error) {
262  FATAL(error);
263  }
264 
265  const JModuleRouter moduleRouter(detector);
266  const JPMTRouter pmtRouter (detector);
267 
268 
269  detector -= get<JPosition3D>(getHeader(inputFile));
270 
271 
272  TFile out(outputFile.c_str(), "recreate");
273 
274  TProfile* he;
275  TProfile* hp;
276 
277  if (histogram == 1) {
278 
280 
281  double x = -0.5;
282 
283  for ( ; x < 30.0; x += 1.0)
284  X.push_back(x);
285 
286  for ( ; x < 50.0; x += 2.0)
287  X.push_back(x);
288 
289  for ( ; x < 100.0; x += 5.0)
290  X.push_back(x);
291 
292  for ( ; x < 200.0; x += 10.0)
293  X.push_back(x);
294 
295  he = new TProfile("he", NULL, X.size() - 1, X.data());
296  hp = new TProfile("hp", NULL, X.size() - 1, X.data());
297 
298  } else {
299 
300  he = new TProfile("he", NULL, 28, 0.0, 7.0);
301  hp = new TProfile("hp", NULL, 28, 0.0, 7.0);
302  }
303 
304  TH1D ht1("ht1", NULL, 550, -50.0, +500.0);
305  TH1D ht2("ht2", NULL, 550, -50.0, +500.0);
306  TH1D hd1("hd1", NULL, 100, -1.0, +1.0);
307  TH1D hd2("hd2", NULL, 100, -1.0, +1.0);
308  TH1D hx1("hx1", NULL, 100, 0.0, +250.0);
309  TH1D hx2("hx2", NULL, 100, 0.0, +250.0);
310  TH1D hw1("hw1", NULL, 100, -0.5, +99.5);
311  TH1D hw2("hw2", NULL, 100, -0.5, +99.5);
312 
313 
314  typedef map<int, JBuffer> JMap_t;
315  typedef vector<JHitL1> JDataL1_t;
316 
317  JMap_t zmap;
318  const JBuildL2<JHitL1> buildL2(JL2Parameters(2, Tmax_ns, ctMin));
319  const JMatch3B<JHitL1> match3B(roadWidth_m, Tmax_ns);
320  const JMatch3D<JHitL1> match3D(Tmax_ns);
321  const JMatch1D<JHitL1> match1D(roadWidth_m, Tmax_ns);
322  const JHorizon horizon(0.2, 2);
323 
324  const JVersor3D gui(getSinThetaC(), 0.0, getCosThetaC()); // photon emission direction
325 
326 
327  while (inputFile.hasNext()) {
328 
329  STATUS("event: " << setw(10) << inputFile.getCounter() << '\r'); DEBUG(endl);
330 
331  JTriggeredFileScanner<>::multi_pointer_type ps = inputFile.next();
332 
333  const JDAQEvent* tev = ps;
334  const Evt* event = ps;
335 
336  if (has_neutrino(*event)) {
337 
338  const Trk& neutrino = get_neutrino(*event);
339  const double E = neutrino.E; // GeV
340 
341 
342  vector<Trk>::const_iterator muon = event->mc_trks.end();
343 
344  {
345  double Emax = 0.0;
346 
347  for (vector<Trk>::const_iterator track = event->mc_trks.begin(); track != event->mc_trks.end(); ++track) {
348 
349  // (anti) muon
350 
351  if (is_muon(*track) && track->E > Emax) {
352  muon = track;
353  Emax = track->E;
354  }
355  }
356  }
357 
358  if (muon != event->mc_trks.end()) {
359 
360  JPosition3D pos = getPosition (*muon);
361  double t0 = muon->t;
362 
363  const JRotation3D R(getDirection(*muon));
364 
365  pos.rotate(R);
366 
367  const time_converter converter(*event, *tev);
368 
369  // Monte Carlo hits
370 
371  zmap.clear();
372 
373  const vector<Hit>& hitlist = event->mc_hits;
374 
375  for (vector<Hit>::const_iterator i = hitlist.begin(); i != hitlist.end(); ++i) {
376  if (!is_noise(*i)) {
377  zmap[pmtRouter.getParentModuleID(i->pmt_id)].push_back(*i);
378  }
379  }
380 
381 
382  // L1 Monte Carlo true data
383 
384  for (JMap_t::iterator i = zmap.begin(); i != zmap.end(); ) {
385 
386  i->second.compile(Tmax_ns, 2);
387 
388  if (i->second.empty())
389  zmap.erase(i++);
390  else
391  ++i;
392  }
393 
394  const int L1mc = zmap.size();
395 
396 
397  // L1 triggered data
398 
399  JDataL1_t dataL1;
400  JDataL1_t zbuf;
401 
402  JDAQTimeslice timeslice(*tev, true);
403 
404  for (JDAQTimeslice::const_iterator i = timeslice.begin(); i != timeslice.end(); ++i) {
405 
406  zbuf.clear();
407 
408  buildL2(*i, moduleRouter.getModule(i->getModuleID()), back_inserter(zbuf));
409 
410  if (!zbuf.empty()) {
411 
412  sort(zbuf.begin(), zbuf.end(), timeSorter<JHitL1>);
413 
414  dataL1.push_back(*zbuf.begin());
415  }
416  }
417 
418 
419  // cluster
420 
421  JDataL1_t::iterator __end = dataL1.begin();
422 
423  switch (cluster) {
424 
425  case 'A':
426  //__end = dataL1.end();
427  __end = clusterize(dataL1.begin(), dataL1.end(), make_match(matchAll));
428  break;
429 
430  case 'B':
431  __end = clusterize(dataL1.begin(), dataL1.end(), weightSorter<JHitL1>, match3D);
432  break;
433 
434  case 'C':
435  __end = clusterize(dataL1.begin(), dataL1.end(), match3B);
436  break;
437 
438  case 'D':
439  __end = clusterizeWeight(dataL1.begin(), dataL1.end(), match3B);
440  break;
441 
442  case 'E':
443 
444  for (JDataL1_t::iterator i = dataL1.begin(); i != dataL1.end(); ++i)
445  i->rotate(R);
446 
447  __end = clusterizeWeight(dataL1.begin(), dataL1.end(), match1D);
448 
449  for (JDataL1_t::iterator i = dataL1.begin(); i != dataL1.end(); ++i)
450  i->rotate_back(R);
451 
452  break;
453 
454  case 'F':
455 
456  __end = clusterizeWeight(dataL1.begin(), dataL1.end(), match3B);
457 
458  for (JDataL1_t::iterator i = dataL1.begin(); i != dataL1.end(); ++i)
459  i->rotate(R);
460 
461  __end = partition(dataL1.begin(), __end, horizon);
462  __end = clusterizeWeight(dataL1.begin(), __end, match1D);
463 
464  for (JDataL1_t::iterator i = dataL1.begin(); i != dataL1.end(); ++i)
465  i->rotate_back(R);
466 
467  break;
468 
469  default:
470  __end = dataL1.end();
471  break;
472  }
473 
474 
475  // L1 monitor
476 
477  JDataL1_t::iterator __q = __end;
478 
479  set<int> L1;
480  set<int> L1ok;
481 
482  for (vector<JHitL1>::iterator hit = dataL1.begin(); hit != __q; ) {
483 
484  L1.insert(hit->getModuleID());
485 
486  const double t1 = converter.getTime(hit->getT());
487 
488  if (zmap[hit->getModuleID()].has(JTimeRange(t1 - 0.5*Tmax_ns, t1 + 0.5*Tmax_ns))) {
489 
490  L1ok.insert(hit->getModuleID());
491 
492  ++hit;
493 
494  } else {
495 
496  iter_swap(hit, --__q);
497  }
498  }
499 
500 
501  if (L1mc != 0 && !L1.empty()) {
502 
503  Double_t x = numeric_limits<Double_t>::max();
504 
505  if (histogram == 1)
506  x = L1mc;
507  else
508  x = log10(E);
509 
510  he->Fill(x, (Double_t) L1ok.size() / (Double_t) L1mc);
511  hp->Fill(x, (Double_t) L1ok.size() / (Double_t) L1.size());
512 
513 
514  // transform hits to muon system
515 
516  for (JDataL1_t::iterator hit = dataL1.begin(); hit != __end; ++hit) {
517  hit->transform(R, pos);
518  }
519 
520 
521  for (JDataL1_t::iterator hit = dataL1.begin(); hit != __end; ++hit) {
522 
523  const double t1 = t0 + getInverseSpeedOfLight() * (hit->getZ() + hit->getX() * getTanThetaC());
524 
525  double dot = +1.0;
526 
527  for (JHitL1::const_iterator i = hit->begin(); i != hit->end(); ++i) {
528  if (i->getDot(gui) < dot)
529  dot = i->getDot(gui);
530  }
531 
532  const double w = (hit->rbegin()->getT() - hit->begin()->getT()); // / hit->size();
533 
534  if (distance(hit, __q) > 0) {
535  ht1.Fill(converter.getTime(hit->getT()) - t1);
536  hd1.Fill(dot);
537  hx1.Fill(hit->getX());
538  hw1.Fill(w);
539  } else {
540  ht2.Fill(converter.getTime(hit->getT()) - t1);
541  hd2.Fill(dot);
542  hx2.Fill(hit->getX());
543  hw2.Fill(w);
544  }
545  }
546  }
547  }
548  }
549  }
550  STATUS(endl);
551 
552  out.Write();
553  out.Close();
554 }
Router for direct addressing of PMT data in detector data structure.
Definition: JPMTRouter.hh:35
Utility class to parse command line options.
Definition: JParser.hh:1514
General exception.
Definition: JException.hh:23
data_type w[N+1][M+1]
Definition: JPolint.hh:778
1D match criterion.
Definition: JMatch1D.hh:31
static struct JTRIGGER::clusterizeWeight clusterizeWeight
then usage $script< input file >[option[primary[working directory]]] nWhere option can be E
Definition: JMuonPostfit.sh:40
std::vector< T >::difference_type distance(typename std::vector< T >::const_iterator first, typename PhysicsEvent::const_iterator< T > second)
Specialisation of STL distance.
#define STATUS(A)
Definition: JMessage.hh:63
bool has_neutrino(const Evt &evt)
Test whether given event has an incoming neutrino.
Detector data structure.
Definition: JDetector.hh:89
bool is_muon(const Trk &track)
Test whether given track is a (anti-)muon.
Router for direct addressing of module data in detector data structure.
Rotation matrix.
Definition: JRotation3D.hh:111
Auxiliary class to synchronously read DAQ events and Monte Carlo events (and optionally other events)...
JMatchHelper< JHit_t > make_match(bool(*match)(const JHit_t &, const JHit_t &))
Auxiliary method to make JMatch object based on pointer to match function.
Definition: JMatch.hh:115
Auxiliary class to convert DAQ hit time to/from Monte Carlo hit time.
string outputFile
double E
Energy [GeV] (either MC truth or reconstructed)
Definition: Trk.hh:20
bool is_noise(const Hit &hit)
Verify hit origin.
Head getHeader(const JMultipleFileScanner_t &file_list)
Get Monte Carlo header.
double getCosThetaC()
Get average cosine of Cherenkov angle of water corresponding to group velocity.
Auxiliary class for defining the range of iterations of objects.
Definition: JLimit.hh:41
Template L2 builder.
Definition: JBuildL2.hh:45
Detector file.
Definition: JHead.hh:226
JDirection3D getDirection(const Vec &dir)
Get direction.
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1989
set_variable E_E log10(E_{fit}/E_{#mu})"
JPosition3D getPosition(const Vec &pos)
Get position.
Data time slice.
#define FATAL(A)
Definition: JMessage.hh:67
Data structure for L2 parameters.
then JCookie sh JDataQuality D $DETECTOR_ID R
Definition: JDataQuality.sh:41
void load(const std::string &file_name, JDetector &detector)
Load detector from input file.
const double getInverseSpeedOfLight()
Get inverse speed of light.
no fit printf nominal n $STRING awk v X
double getTanThetaC()
Get average tangent of Cherenkov angle of water corresponding to group velocity.
Data structure for position in three dimensions.
Definition: JPosition3D.hh:36
const JLimit & getLimit() const
Get limit.
Definition: JLimit.hh:84
const Trk & get_neutrino(const Evt &evt)
Get incoming neutrino.
do set_variable DETECTOR_TXT $WORKDIR detector
Data structure for normalised vector in three dimensions.
Definition: JVersor3D.hh:26
General purpose class for multiple pointers.
The Trk class represents a Monte Carlo (MC) particle as well as a reconstructed track/shower.
Definition: Trk.hh:14
JPosition3D & rotate(const JRotation3D &R)
Rotate.
Definition: JPosition3D.hh:186
static struct JTRIGGER::clusterize clusterize
3D match criterion.
Definition: JMatch3D.hh:29
int debug
debug level
The Evt class respresent a Monte Carlo (MC) event as well as an offline event.
Definition: Evt.hh:20
JTOOLS::JRange< double > JTimeRange
Type definition for time range (unit [s]).
double getSinThetaC()
Get average sine of Cherenkov angle of water corresponding to group velocity.
3D match criterion with road width.
Definition: JMatch3B.hh:34
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:62