Jpp  18.3.0
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 
105 
106  try {
107  load(detectorFile, detector);
108  }
109  catch(const JException& error) {
110  FATAL(error);
111  }
112 
113  // -----------------------------------
114  // STEP 0: prepare
115  // -----------------------------------
116 
117  // configure input streams
118 
120 
122 
123  JTreeScannerInterface<JDAQTimeslice>* pts = zmap[selector];
124 
125  pts->configure(inputFile);
126 
127  // detect input stream size
128 
129 
130  int fEnd = pts->rbegin()->getFrameIndex();
131  int fStart = pts->begin( )->getFrameIndex();
132 
133 
134  // crop histogram if -n is specified
135  if (fEnd > inputFile.getUpperLimit()) {
136  fEnd = fStart + inputFile.getUpperLimit();
137  }
138 
139  int fLength = 1 + fEnd - fStart;
140 
141  NOTICE("begin | end | length = " << fStart << " | " << fEnd << " | " << fLength << endl);
142 
143 
144  if (fStart > fEnd) {
145  FATAL("FATAL: inconsistent TTree indexing" << endl);
146  }
147 
148  // detector routers
149 
150  const JModuleRouter moduleRouter(detector);
151 
152  const JDAQHitRouter hitRouter(detector);
153 
154  // data structures
155 
156  const int nx = fLength;
157  const double xmin = fStart;
158  const double xmax = fEnd + 1;
159 
160  const int ny = NUMBER_OF_PMTS + 1;
161  const double ymin = -0.5;
162  const double ymax = ymin + ny;
163 
164  typedef JManager <string, TH2D> JManager2D_t;
165  typedef JManager <string, TH1D> JManager1D_t;
166 
167  JManager2D_t MT(new TH2D(mul_p , NULL, nx, xmin, xmax, ny, ymin, ymax));
168  JManager1D_t ST(new TH1D(status_p, NULL, nx, xmin, xmax));
169 
170  // -----------------------------------
171  // STEP 1: building vetoes from events
172  // -----------------------------------
173 
174  int runNumber = 0;
175 
176  TH1D* h_vtr = new TH1D("VetoTimeRange","VetoTimeRange", 10000, 0, 10000);
177 
179 
180  map<int, JVetoSet> triggeredEvents;
181 
182  for (; evIn.hasNext(); ) {
183 
184  STATUS("event: " << setw(10) << evIn.getCounter() << '\r'); DEBUG(endl);
185 
186  JDAQEvent* event = evIn.next();
187 
188  if (!runNumber) { runNumber = event->getRunNumber(); }
189 
190  JVeto vp(*event, hitRouter);
191 
192  triggeredEvents[event->getFrameIndex()].push_back(vp);
193 
194  h_vtr->Fill(vp.getLength());
195 
196  }
197 
198  STATUS(triggeredEvents.size() << " JDAQEvents loaded in veto buffer." << endl);
199 
200  TParameter<int> RUNNR("RUNNR", runNumber);
201  // TParameter<int> TSTART("TSTART", tStart)
202 
203  //-----------------------------
204  // STEP 2: timeslice processing
205  // ----------------------------
206 
207  counter_type counter = 0;
208 
209  for ( ; pts->hasNext() && counter != inputFile.getLimit(); ++counter) {
210 
211  STATUS("timeslice: " << setw(10) << counter << '\r'); DEBUG(endl);
212 
213  const JDAQTimeslice* timeslice = pts->next();
214 
215  int fIndex = timeslice->getFrameIndex();
216 
217  if (counter == 0) {
218  NOTICE("Start frame index = " << fIndex << " @ T " << timeslice->getTimesliceStart() << endl);
219  }
220 
221  JVetoSet veto;
222 
223  if (triggeredEvents.count(fIndex)) {
224  veto = triggeredEvents.at(fIndex);
225  }
226 
227  // L0-L1 processing
228 
229  typedef JHitR0 hit_type;
230 
231  typedef JSuperFrame2D<hit_type> JSuperFrame2D_t;
232 
233  typedef JSuperFrame1D<hit_type> JSuperFrame1D_t;
234 
235  const JMatchL0<hit_type> match(TMax_ns);
236 
237  double fractionActivePMTs = 0;
238 
239  int nDOMs = timeslice->size();
240 
241  for (JDAQTimeslice::const_iterator frame = timeslice->begin(); frame != timeslice->end(); ++frame) {
242 
243  int moduleID = frame->getModuleID();
244 
245  fractionActivePMTs += ((double) frame->countActiveChannels());
246 
247  JSuperFrame2D_t& buffer = JSuperFrame2D_t::demultiplex(*frame,
248  moduleRouter.getModule(moduleID));
249  // totSelector_ns);
250 
251  buffer.preprocess(JPreprocessor::join_t, match);
252 
253  JSuperFrame1D_t& data = JSuperFrame1D_t::multiplex(buffer);
254 
255  // sort(data.begin(), data.end()); // JSuperFrame1D is sorted
256 
257  // histogram for time differences, may not be used
258 
259  TH1D* h2dt = new TH1D("H2DT", "H2DT", 100, -TMax_ns, +TMax_ns);
260 
261  for (vector<JHitR0>::const_iterator p = data.begin(); p != data.end(); ) {
262 
264 
265  while (++q != data.end() && q->getT() - p->getT() <= TMax_ns) {}
266 
267  int m = distance(p, q);
268 
269  // raw multiplicity count //
270 
271  MT["RAW"]->Fill(fIndex, m);
272 
273  // time difference distribution
274 
275  if (selector != "JDAQTimesliceSN" && timeDifferences && m > 1) {
276  const double W = factorial(m, 2);
277  for (vector<JHitR0>::const_iterator __p = p; __p != q; ++__p) {
278  for (vector<JHitR0>::const_iterator __q = __p; ++__q != q; ) {
279  double dt = JCombinatorics::getSign(__p->getPMT(), __q->getPMT()) * (__q->getT() - __p->getT());
280  h2dt->Fill(dt, 1.0/W);
281  }
282  }
283  }
284 
285  // slide iterator
286 
287  p = q;
288 
289  } // end hit processing
290 
291  if (h2dt->GetEntries() > 0) {
292 
293  TF1 f("f", "[0]*exp(-0.5*(x-[1])*(x-[1])/([2]*[2]))/(TMath::Sqrt(2*TMath::Pi())*[2]) + [3]");
294 
295  f.SetParameter(0, h2dt->GetMaximum());
296  f.SetParameter(1, h2dt->GetMean());
297  f.SetParameter(2, h2dt->GetRMS() * 0.25);
298  f.SetParameter(3, h2dt->GetMinimum());
299 
300  h2dt->Fit(&f, "Q", "same");
301 
302  double nb = h2dt->GetNbinsX();
303  double bg_v = f.GetParameter(3) * nb;
304  double sg = h2dt->GetSumOfWeights() - bg_v;
305 
306  // inclusive 2-fold after subtraction of fitted background //
307 
308  MT["FIT"]->Fill(fIndex, 2, sg);
309 
310  }
311 
312  delete h2dt;
313 
314  }
315 
316  // end timeslice processing
317 
318  fractionActivePMTs /= (NUMBER_OF_PMTS * nDOMs);
319 
320  ST["PMT"]->Fill(fIndex, fractionActivePMTs);
321 
322  // STANDARD PROCESSING
323 
324  JDataSN preTrigger(TMax_ns, preTriggerThreshold);// totSelector_ns);
325 
326  preTrigger(timeslice, moduleRouter, totSelector_ns);
327 
328  JTriggerSN trigger(TVeto_ns);
329 
330  trigger(preTrigger);
331 
332  // generate and configure event-based veto
333 
334  // count trigger
335 
336  JRange_t A = JRange_t(4,31);
337 
338  // select single-module clusters with max multiplicity in range A
339  JSNFilterM trgA1(A, 1);
340 
341  // select non-vetoed clusters with max multiplicity in range A
342  JSNFilterMV trgAV(A, veto);
343 
344  vector<double> m_a1 = trigger.getMultiplicities(trgA1);
345 
346  for (vector<double>::const_iterator p = m_a1.begin(); p != m_a1.end(); p++) {
347  MT["TA1"]->Fill(fIndex, *p);
348  }
349 
350  vector<double> m_av = trigger.getMultiplicities(trgAV);
351 
352  for (vector<double>::const_iterator p = m_av.begin(); p != m_av.end(); p++) {
353  MT["TAV"]->Fill(fIndex, *p);
354  }
355 
356  }
357 
358  //-------------------------------------
359  // STEP 3: generating trigger summaries
360  // ------------------------------------
361 
362  if (outputFile != "") {
363 
364  TFile out(outputFile.c_str(), "RECREATE");
365 
366  putObject(out, JMeta(argc,argv));
367 
368  RUNNR.Write();
369  // TSTART.Write();
370 
371  h_vtr->Write();
372 
373  MT.Write(out);
374  ST.Write(out);
375 
376  out.Close();
377 
378  }
379 
380 }
Auxiliary class for ROOT I/O of application specific meta data.
Definition: JMeta.hh:70
const double xmax
Definition: JQuadrature.cc:24
Utility class to parse command line options.
Definition: JParser.hh:1514
General exception.
Definition: JException.hh:24
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:89
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.
Auxiliary interface for direct access of elements in ROOT TChain.
o $QUALITY_ROOT d $DEBUG!CHECK_EXIT_CODE JPlot1D f
Definition: JDataQuality.sh:76
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.
long long int factorial(const long long int n)
Determine factorial.
Definition: JMathToolkit.hh:42
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
Type definition of range.
Definition: JHead.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:1989
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
bool putObject(TDirectory &dir, const TObject &object)
Write object to ROOT directory.
Auxiliary class to select JTreeScanner based on ROOT class name.
#define FATAL(A)
Definition: JMessage.hh:67
const double xmin
Definition: JQuadrature.cc:23
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 &dev null eval JShellParser o a A
int debug
debug level
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