Jpp 19.3.0-rc.2
the software that should make you happy
Loading...
Searching...
No Matches
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 "JROOT/JMinimizer.hh"
#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

◆ main()

int main ( int argc,
char ** argv )

Definition at line 62 of file JKexing2D.cc.

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