Jpp - the software that should make you happy
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Functions
JKexing2D.cc File Reference
#include <string>
#include <iostream>
#include <fstream>
#include <iomanip>
#include <vector>
#include "TROOT.h"
#include "TFile.h"
#include "TH2D.h"
#include "TF1.h"
#include "TParameter.h"
#include "JDAQ/JDAQEvaluator.hh"
#include "JDetector/JDetector.hh"
#include "JDetector/JDetectorToolkit.hh"
#include "JROOT/JManager.hh"
#include "JROOT/JROOTClassSelector.hh"
#include "JSupport/JMultipleFileScanner.hh"
#include "JSupport/JTreeScanner.hh"
#include "JSupport/JAutoTreeScanner.hh"
#include "JSupport/JSupport.hh"
#include "JSupport/JMeta.hh"
#include "JLang/JObjectMultiplexer.hh"
#include "JMath/JMathToolkit.hh"
#include "Jeep/JPrint.hh"
#include "Jeep/JParser.hh"
#include "Jeep/JMessage.hh"
#include "JTools/JCombinatorics.hh"
#include "JTools/JQuantile.hh"
#include "km3net-dataformat/online/JDAQHit.hh"
#include "JTrigger/JHit.hh"
#include "JTrigger/JHitR0.hh"
#include "JTrigger/JSuperFrame2D.hh"
#include "JTrigger/JDAQHitToTSelector.hh"
#include "JSupernova/JSupernova.hh"
#include "JSupernova/JKexing2D.hh"

Go to the source code of this file.

Functions

int main (int argc, char **argv)
 

Detailed Description

Author
mlincett

Example application to analyse intra-run supernova background

Definition in file JKexing2D.cc.

Function Documentation

int main ( int  argc,
char **  argv 
)

Definition at line 60 of file JKexing2D.cc.

61 {
62 
63  typedef JRange<int> JRange_t;
64  typedef JRange<JDAQHit::JTOT_t> JToTRange_t;
65 
66  JMultipleFileScanner<> inputFile;
67  JLimit_t& numberOfEvents = inputFile.getLimit();
68  string outputFile;
69  string detectorFile;
70  double TMax_ns;
71  double TVeto_ns;
72  JToTRange_t totRange_ns;
73  JDAQHitToTSelector totSelector_ns;
74  JROOTClassSelector selector;
75  bool timeDifferences;
76  int debug;
77  int preTriggerThreshold;
78 
79  // Check input parameters
80 
81  try {
82 
83  JParser<> zap("Example application to study supernova detection background");
84 
85  zap['f'] = make_field(inputFile);
86  zap['o'] = make_field(outputFile) = "kexing2D.root";
87  zap['n'] = make_field(numberOfEvents) = JLimit::max();
88  zap['a'] = make_field(detectorFile);
89  zap['T'] = make_field(TMax_ns) = 10.0;
90  zap['V'] = make_field(TVeto_ns) = 1000.0;
91  zap['D'] = make_field(timeDifferences);
92  zap['C'] = make_field(selector) = getROOTClassSelection<JDAQTimesliceTypes_t>();
93  zap['d'] = make_field(debug) = 1;
94  zap['P'] = make_field(preTriggerThreshold) = 4;
95  zap['t'] = make_field(totSelector_ns) = JDAQHitToTSelector(4, 255);
96 
97  zap(argc, argv);
98  }
99  catch(const exception &error) {
100  FATAL(error.what() << endl);
101  }
102 
103  cout.tie(&cerr);
104 
106 
107  try {
108  load(detectorFile, detector);
109  }
110  catch(const JException& error) {
111  FATAL(error);
112  }
113 
114  // -----------------------------------
115  // STEP 0: prepare
116  // -----------------------------------
117 
118  // configure input streams
119 
121 
123 
124  JTreeScannerInterface<JDAQTimeslice>* pts = zmap[selector];
125 
126  pts->configure(inputFile);
127 
128  // detect input stream size
129 
130 
131  int fEnd = pts->rbegin()->getFrameIndex();
132  int fStart = pts->begin( )->getFrameIndex();
133 
134 
135  // crop histogram if -n is specified
136  if (fEnd > inputFile.getUpperLimit()) {
137  fEnd = fStart + inputFile.getUpperLimit();
138  }
139 
140  int fLength = 1 + fEnd - fStart;
141 
142  NOTICE("begin | end | length = " << fStart << " | " << fEnd << " | " << fLength << endl);
143 
144 
145  if (fStart > fEnd) {
146  FATAL("FATAL: inconsistent TTree indexing" << endl);
147  }
148 
149  // detector routers
150 
151  const JModuleRouter moduleRouter(detector);
152 
153  const JDAQHitRouter hitRouter(detector);
154 
155  // data structures
156 
157  const int nx = fLength;
158  const double xmin = fStart;
159  const double xmax = fEnd + 1;
160 
161  const int ny = NUMBER_OF_PMTS + 1;
162  const double ymin = -0.5;
163  const double ymax = ymin + ny;
164 
165  typedef JManager <string, TH2D> JManager2D_t;
166  typedef JManager <string, TH1D> JManager1D_t;
167 
168  JManager2D_t MT(new TH2D(mul_p , NULL, nx, xmin, xmax, ny, ymin, ymax));
169  JManager1D_t ST(new TH1D(status_p, NULL, nx, xmin, xmax));
170 
171  // -----------------------------------
172  // STEP 1: building vetoes from events
173  // -----------------------------------
174 
175  int runNumber = 0;
176 
177  TH1D* h_vtr = new TH1D("VetoTimeRange","VetoTimeRange", 10000, 0, 10000);
178 
180 
181  map<int, JVetoSet> triggeredEvents;
182 
183  for (; evIn.hasNext(); ) {
184 
185  STATUS("event: " << setw(10) << evIn.getCounter() << '\r'); DEBUG(endl);
186 
187  JDAQEvent* event = evIn.next();
188 
189  if (!runNumber) { runNumber = event->getRunNumber(); }
190 
191  JVeto vp(*event, hitRouter);
192 
193  triggeredEvents[event->getFrameIndex()].push_back(vp);
194 
195  h_vtr->Fill(vp.getLength());
196 
197  }
198 
199  STATUS(triggeredEvents.size() << " JDAQEvents loaded in veto buffer." << endl);
200 
201  TParameter<int> RUNNR("RUNNR", runNumber);
202  // TParameter<int> TSTART("TSTART", tStart)
203 
204  //-----------------------------
205  // STEP 2: timeslice processing
206  // ----------------------------
207 
208  counter_type counter = 0;
209 
210  for ( ; pts->hasNext() && counter != inputFile.getLimit(); ++counter) {
211 
212  STATUS("timeslice: " << setw(10) << counter << '\r'); DEBUG(endl);
213 
214  const JDAQTimeslice* timeslice = pts->next();
215 
216  int fIndex = timeslice->getFrameIndex();
217 
218  if (counter == 0) {
219  NOTICE("Start frame index = " << fIndex << " @ T " << timeslice->getTimesliceStart() << endl);
220  }
221 
222  JVetoSet veto;
223 
224  if (triggeredEvents.count(fIndex)) {
225  veto = triggeredEvents.at(fIndex);
226  }
227 
228  // L0-L1 processing
229 
230  typedef JHitR0 hit_type;
231 
232  typedef JSuperFrame2D<hit_type> JSuperFrame2D_t;
233 
234  typedef JSuperFrame1D<hit_type> JSuperFrame1D_t;
235 
236  const JMatchL0<hit_type> match(TMax_ns);
237 
238  double fractionActivePMTs = 0;
239 
240  int nDOMs = timeslice->size();
241 
242  for (JDAQTimeslice::const_iterator frame = timeslice->begin(); frame != timeslice->end(); ++frame) {
243 
244  int moduleID = frame->getModuleID();
245 
246  fractionActivePMTs += ((double) frame->countActiveChannels());
247 
248  JSuperFrame2D_t& buffer = JSuperFrame2D_t::demultiplex(*frame,
249  moduleRouter.getModule(moduleID));
250  // totSelector_ns);
251 
252  buffer.preprocess(JPreprocessor::join_t, match);
253 
254  JSuperFrame1D_t& data = JSuperFrame1D_t::multiplex(buffer); data.pop_back();
255 
256  // sort(data.begin(), data.end()); // JSuperFrame1D is sorted
257 
258  // histogram for time differences, may not be used
259 
260  TH1D* h2dt = new TH1D("H2DT", "H2DT", 100, -TMax_ns, +TMax_ns);
261 
262  for (vector<JHitR0>::const_iterator p = data.begin(); p != data.end(); ) {
263 
265 
266  while (++q != data.end() && q->getT() - p->getT() <= TMax_ns) {}
267 
268  int m = distance(p, q);
269 
270  // raw multiplicity count //
271 
272  MT["RAW"]->Fill(fIndex, m);
273 
274  // time difference distribution
275 
276  if (selector != "JDAQTimesliceSN" && timeDifferences && m > 1) {
277  const double W = factorial(m, 2);
278  for (vector<JHitR0>::const_iterator __p = p; __p != q; ++__p) {
279  for (vector<JHitR0>::const_iterator __q = __p; ++__q != q; ) {
280  double dt = JCombinatorics::getSign(__p->getPMT(), __q->getPMT()) * (__q->getT() - __p->getT());
281  h2dt->Fill(dt, 1.0/W);
282  }
283  }
284  }
285 
286  // slide iterator
287 
288  p = q;
289 
290  } // end hit processing
291 
292  if (h2dt->GetEntries() > 0) {
293 
294  TF1 f("f", "[0]*exp(-0.5*(x-[1])*(x-[1])/([2]*[2]))/(TMath::Sqrt(2*TMath::Pi())*[2]) + [3]");
295 
296  f.SetParameter(0, h2dt->GetMaximum());
297  f.SetParameter(1, h2dt->GetMean());
298  f.SetParameter(2, h2dt->GetRMS() * 0.25);
299  f.SetParameter(3, h2dt->GetMinimum());
300 
301  h2dt->Fit(&f, "Q", "same");
302 
303  double nb = h2dt->GetNbinsX();
304  double bg_v = f.GetParameter(3) * nb;
305  double sg = h2dt->GetSumOfWeights() - bg_v;
306 
307  // inclusive 2-fold after subtraction of fitted background //
308 
309  MT["FIT"]->Fill(fIndex, 2, sg);
310 
311  }
312 
313  delete h2dt;
314 
315  }
316 
317  // end timeslice processing
318 
319  fractionActivePMTs /= (NUMBER_OF_PMTS * nDOMs);
320 
321  ST["PMT"]->Fill(fIndex, fractionActivePMTs);
322 
323  // STANDARD PROCESSING
324 
325  JDataSN preTrigger(TMax_ns, preTriggerThreshold);// totSelector_ns);
326 
327  preTrigger(timeslice, moduleRouter, totSelector_ns);
328 
329  JTriggerSN trigger(TVeto_ns);
330 
331  trigger(preTrigger);
332 
333  // generate and configure event-based veto
334 
335  // count trigger
336 
337  JRange_t A = JRange_t(4,31);
338 
339  // select single-module clusters with max multiplicity in range A
340  JSNFilterM trgA1(A, 1);
341 
342  // select non-vetoed clusters with max multiplicity in range A
343  JSNFilterMV trgAV(A, veto);
344 
345  vector<double> m_a1 = trigger.getMultiplicities(trgA1);
346 
347  for (vector<double>::const_iterator p = m_a1.begin(); p != m_a1.end(); p++) {
348  MT["TA1"]->Fill(fIndex, *p);
349  }
350 
351  vector<double> m_av = trigger.getMultiplicities(trgAV);
352 
353  for (vector<double>::const_iterator p = m_av.begin(); p != m_av.end(); p++) {
354  MT["TAV"]->Fill(fIndex, *p);
355  }
356 
357  }
358 
359  //-------------------------------------
360  // STEP 3: generating trigger summaries
361  // ------------------------------------
362 
363  if (outputFile != "") {
364 
365  TFile out(outputFile.c_str(), "RECREATE");
366 
367  putObject(&out, JMeta(argc,argv));
368 
369  RUNNR.Write();
370  // TSTART.Write();
371 
372  h_vtr->Write();
373 
374  MT.Write(out);
375  ST.Write(out);
376 
377  out.Close();
378 
379  }
380 
381 }
Auxiliary class for ROOT I/O of application specific meta data.
Definition: JMeta.hh:70
Utility class to parse command line options.
Definition: JParser.hh:1500
General exception.
Definition: JException.hh:23
JDAQUTCExtended getTimesliceStart() const
Get start of timeslice.
Auxiliary class to define a veto time window on a set of optical modules.
Definition: JSupernova.hh:80
std::vector< T >::difference_type distance(typename std::vector< T >::const_iterator first, typename PhysicsEvent::const_iterator< T > second)
Specialisation of STL distance.
L0 match criterion.
Definition: JMatchL0.hh:27
static const char * mul_p
Definition: JKexing2D.hh:11
#define STATUS(A)
Definition: JMessage.hh:63
Detector data structure.
Definition: JDetector.hh:80
static const char * status_p
Definition: JKexing2D.hh:12
Auxiliary class to select ROOT class based on class name.
SN filter based on veto window.
Definition: JSupernova.hh:364
Router for direct addressing of module data in detector data structure.
do set_array DAQHEADER JPrintDAQHeader f
Definition: JTuneHV.sh:79
Long64_t counter_type
Type definition for counter.
Auxiliary class for a type holder.
Definition: JType.hh:19
string outputFile
Template definition for direct access of elements in ROOT TChain.
Definition: JTreeScanner.hh:91
long long int factorial(const long long int n)
Determine factorial.
Definition: JMathToolkit.hh:42
Auxiliary interface for direct access of elements in ROOT TChain.
Simple wrapper around JModuleRouter class for direct addressing of PMT data in detector data structur...
int getFrameIndex() const
Get frame index.
1-dimensional frame with time calibrated data from one optical module.
Auxiliary class for defining the range of iterations of objects.
Definition: JLimit.hh:41
Auxiliary class to manage set of compatible ROOT objects (e.g. histograms) using unique keys...
Definition: JManager.hh:43
Reduced data structure for L0 hit.
Definition: JHitR0.hh:25
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1961
Auxiliary class to apply the supernova trigger to SN data.
Definition: JSupernova.hh:391
#define NOTICE(A)
Definition: JMessage.hh:64
Data time slice.
SN filter based on multiplicity selection optional suppression of multi-module coincidences WARNING: ...
Definition: JSupernova.hh:329
int debug
debug level
Definition: JSirene.cc:63
Auxiliary class to select JTreeScanner based on ROOT class name.
#define FATAL(A)
Definition: JMessage.hh:67
Auxiliary class to build the supernova trigger dataset.
Definition: JSupernova.hh:134
void load(const std::string &file_name, JDetector &detector)
Load detector from input file.
2-dimensional frame with time calibrated data from one optical module.
Auxiliary class to select DAQ hits based on time-over-treshold value.
do set_variable DETECTOR_TXT $WORKDIR detector
static const int NUMBER_OF_PMTS
Total number of PMTs in module.
Definition: JDAQ.hh:26
source $JPP_DIR setenv csh $JPP_DIR eval JShellParser o a A
JTriggerCounter_t next()
Increment trigger counter.
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:62
Auxiliary class to manage a set of vetoes.
Definition: JSupernova.hh:257
bool putObject(TDirectory *dir, const T &object)
Write object to ROOT directory.