Jpp test-rotations-old
the software that should make you happy
Loading...
Searching...
No Matches
JRipple.cc File Reference
#include <string>
#include <iostream>
#include <fstream>
#include <iomanip>
#include <vector>
#include "TROOT.h"
#include "TFile.h"
#include "TH1D.h"
#include "TH2F.h"
#include "TF1.h"
#include "km3net-dataformat/online/JDAQClock.hh"
#include "JDAQ/JDAQEvaluator.hh"
#include "JDAQ/JDAQSummarysliceIO.hh"
#include "JDetector/JDetector.hh"
#include "JDetector/JDetectorToolkit.hh"
#include "Jeep/JPrint.hh"
#include "Jeep/JParser.hh"
#include "Jeep/JMessage.hh"
#include "JTrigger/JDAQHitToTSelector.hh"
#include "JTrigger/JHitR0.hh"
#include "JTrigger/JMatchL0.hh"
#include "JTrigger/JBuildL0.hh"
#include "JTrigger/JSummaryRouter.hh"
#include "JTrigger/JSuperFrame1D.hh"
#include "JTrigger/JSuperFrame2D.hh"
#include "JTrigger/JTimesliceRouter.hh"
#include "JTools/JCombinatorics.hh"
#include "JROOT/JManager.hh"
#include "JLang/JObjectMultiplexer.hh"
#include "JMath/JMathToolkit.hh"
#include "JROOT/JROOTClassSelector.hh"
#include "JSupernova/JSupernova.hh"
#include "JSupport/JMultipleFileScanner.hh"
#include "JSupport/JTreeScanner.hh"
#include "JSupport/JAutoTreeScanner.hh"
#include "JSupport/JSupport.hh"
#include "JSupport/JSupportToolkit.hh"
#include "JRipple.hh"

Go to the source code of this file.

Functions

int main (int argc, char **argv)
 

Detailed Description

Author
mlincett

Example application to examine hit rates (and active channels) as a function of time.

Definition in file JRipple.cc.

Function Documentation

◆ main()

int main ( int argc,
char ** argv )

Definition at line 66 of file JRipple.cc.

67{
68
69 JMultipleFileScanner<> inputFile;
70 JLimit_t& numberOfEvents = inputFile.getLimit();
71 string outputFile;
72 string detectorFile;
73 double TMax_ns;
74 JROOTClassSelector selector;
75 int debug;
76 int binWidth_ms;
77 bool backVeto;
78 bool globalOutputOnly;
79 bool mode2D;
80 JRange<int> multiplicityRange;
81 JDAQHitToTSelector totSelector_ns;
82 // Check input parameters
83
84 try {
85
86 JParser<> zap("Example program to examine rates as a function of time on ms-level timescales.");
87
88 zap['f'] = make_field(inputFile);
89 zap['o'] = make_field(outputFile) = "ripple.root";
90 zap['n'] = make_field(numberOfEvents) = JLimit::max();
91 zap['a'] = make_field(detectorFile);
93 zap['T'] = make_field(TMax_ns, "Time window for local coincidences (if 0 run in L0 mode)") = 0.0;
94 zap['B'] = make_field(binWidth_ms, "Bin width (experimental)") = 1;
95 zap['V'] = make_field(backVeto, "Apply retroactive veto.");
96 zap['D'] = make_field(mode2D, "L1 mode: create 2D histogram with time differences of coincidences (heavy memory usage, ignored if TMax_ns = 0)");
97 zap['O'] = make_field(globalOutputOnly, "Write only aggregate histograms");
98 zap['M'] = make_field(multiplicityRange, "L1 mode: multiplicity range (ignored if TMax_ns = 0)") = JRange<int>(2,31);
99 zap['t'] = make_field(totSelector_ns, "ToT acceptance range") = JDAQHitToTSelector(3, 255);
100 zap['d'] = make_field(debug) = 1;
101
102
103 zap(argc, argv);
104 }
105 catch(const exception &error) {
106 FATAL(error.what() << endl);
107 }
108
109
110 if (((int)(getFrameTime() / 1e6)) % binWidth_ms) {
111 FATAL("Frame time must be an integer multiple of bin width");
112 }
113
115
116 try {
117 load(detectorFile, detector);
118 }
119 catch(const JException& error) {
120 FATAL(error);
121 }
122
123 // Configure input streams
124
126
128
129 JTreeScannerInterface<JDAQTimeslice>* pts = zmap[selector];
130
131 pts->configure(inputFile);
132
133 int fEnd = pts->rbegin()->getFrameIndex();
134 int fStart = pts->begin( )->getFrameIndex();
135
136 if (fEnd > inputFile.getUpperLimit()) {
137 fEnd = fStart + inputFile.getUpperLimit();
138 }
139
140 double tEnd_ms = fEnd * getFrameTime() / 1.0e6;
141 double tStart_ms = (fStart - 1) * getFrameTime() / 1.0e6;
142
143 int runNumber = pts->begin()->getRunNumber();
144
145 TString runTag = Form("%d" , runNumber);
146
147 double tRun_ms = tEnd_ms - tStart_ms;
148
149 NOTICE("START/END/DURATION [s]: " << tStart_ms / 1000 << " " << tEnd_ms / 1000 << " " << tRun_ms / 1000 << endl);
150
151 // hit containers and data routers
152
153 typedef JHitR0 hit_type;
154 typedef JSuperFrame2D<hit_type> JSuperFrame2D_t;
155 typedef JSuperFrame1D<hit_type> JSuperFrame1D_t;
156
157 const JMatchL0<hit_type> match(TMax_ns); // time window self-coincidences
158
159 const JModuleRouter moduleRouter(detector);
160 const JDAQHitRouter hitRouter(detector);
161
162 JSummaryRouter summaryRouter;
163
164 // output histograms
165
166 typedef JManager<string, TH1F> JManager_t;
167
168 int nx = (tEnd_ms - tStart_ms) / binWidth_ms;
169
170 JManager_t RT_DOM(new TH1F("RT_%", NULL, nx , tStart_ms, tEnd_ms));
171 JManager_t NC_DOM(new TH1F("NC_%", NULL, nx / 100, tStart_ms, tEnd_ms));
172
173 TString rt_tag = runTag + TString(".RT_DET_%");
174 TString nc_tag = runTag + TString(".NC_DET_%");
175
176 JManager_t RT_DET(new TH1F(rt_tag, NULL, nx, tStart_ms, tEnd_ms));
177 JManager_t NC_DET(new TH1F(nc_tag, NULL, nx / 100, tStart_ms, tEnd_ms));
178
179 h2d_t* hT = NULL;
180
181
182 if (mode2D) {
183
184 const int ny = 50;
185
186 const int max_size = floor(h2d_limit / (sizeof(h2d_bintype) * ny));
187
188 if (nx >= max_size) {
189 FATAL("2D histogram size not supported by ROOT file output; limit input size (-n) below " << floor(max_size / 100.0) << endl);
190 }
191
192 TString rt2d_tag = runTag + TString(".RT2D_DET");
193
194 hT = new h2d_t(rt2d_tag, NULL, nx, tStart_ms, tEnd_ms, ny, -TMax_ns, +TMax_ns);
195 hT->SetDirectory(0);
196
197 }
198
199 //---------------------
200 // Timeslice processing
201 //---------------------
202
203 counter_type counter = 0;
204
205 // preload first timeslice
206
207 JDAQTimeslice curr;
208
209 if (pts->hasNext()) {
210 curr = *(pts->next());
211 ++counter;
212 } else {
213 FATAL("Input file is too short.");
214 }
215
216 // loop over timeslices
217
218 for ( ; pts->hasNext() && counter != inputFile.getLimit(); ++counter) {
219
220 STATUS("timeslice: " << setw(10) << counter << '\r'); DEBUG(endl);
221
222 JDAQTimeslice next = *(pts->next());
223
224 const int ic = curr.getFrameIndex();
225
226 const int in = next.getFrameIndex();
227
228 const double tTimeslice_ms = getTimeOfRTS(ic) / 1.0e6;
229
230 // look up next summaryslice
231
232 JDAQSummaryslice* nextSummary = NULL;
233
234 if (backVeto && ((in - ic) == 1)) {
235 nextSummary = new JDAQSummaryslice(next);
236 summaryRouter.update(nextSummary);
237 }
238
239 // loop over frames
240
241 for (JDAQTimeslice::const_iterator frame = curr.begin(); frame != curr.end(); ++frame) {
242
243 const int moduleID = frame->getModuleID();
244 const JModule& module = moduleRouter.getModule(moduleID);
245 const string moduleLabel = getLabel(module);
246
247 TH1F* RD = RT_DOM[moduleLabel];
248
249 vector<bool> veto(NUMBER_OF_PMTS, false);
250
251 // current veto
252
253 if (frame->testHighRateVeto() || frame->testFIFOStatus()) {
254 for (int pmt = 0; pmt < NUMBER_OF_PMTS; pmt++) {
255 veto[pmt] = ( frame->testHighRateVeto(pmt) || frame->testFIFOStatus(pmt) );
256 }
257 }
258
259 // retroactive veto
260
261 if (nextSummary != NULL) {
262
263 JDAQSummaryFrame nextFrame = summaryRouter.getSummaryFrame(moduleID);
264
265 if (nextFrame.testHighRateVeto() || nextFrame.testFIFOStatus()) {
266 for (int pmt = 0; pmt < NUMBER_OF_PMTS; pmt++) {
267 veto[pmt] = veto[pmt] || ( nextFrame.testHighRateVeto(pmt) || nextFrame.testFIFOStatus(pmt) );
268 }
269 }
270
271 }
272
273 // track number of active channel
274
275 NC_DOM[moduleLabel]->Fill(tTimeslice_ms, count(veto.begin(), veto.end(), false));
276
277 // veto
278
279 JSuperFrame2D_t& buffer2D = JSuperFrame2D_t::demultiplex(*frame, module, totSelector_ns);
280
281 for (JSuperFrame2D_t::iterator i = buffer2D.begin(); i != buffer2D.end(); ++i) {
282 if (veto[i->getPMTAddress()]) {
283 i->reset();
284 }
285 }
286
287 buffer2D.preprocess(JPreprocessor::join_t, match);
288
289 JSuperFrame1D_t& data = JSuperFrame1D_t::multiplex(buffer2D);
290
291
292 // data process
293
294 if (data.size() > 1) {
295
296 // sort for L1 processing
297
298 if (TMax_ns > 0) {
299 sort(data.begin(), data.end());
300 }
301
302 // loop over hits in frame
303
304 for (vector<hit_type>::const_iterator p = data.begin(); p != data.end(); p++) {
305
306 const double tHit_ms = tTimeslice_ms + (p->getT() / 1.0e6);
307
308 if (TMax_ns == 0.0) { // L0 mode
309
310 RD->Fill(tHit_ms);
311
312 } else { // L1 mode
313
314 vector<hit_type>::const_iterator q = p;
315
316 while (++q != data.end() && q->getT() - p->getT() <= TMax_ns) {}
317
318 int M = distance(p, q);
319
320 if (multiplicityRange(M)) {
321
322 RD->Fill(tHit_ms);
323
324 if (mode2D) { // L1 2D mode
325
326 const double W = factorial(M, 2);
327
328 for (vector<hit_type>::const_iterator __p = p; __p != q; ++__p) {
329
330 for (vector<hit_type>::const_iterator __q = __p; ++__q != q; ) {
331
332 double dt = JCombinatorics::getSign(__p->getPMT(), __q->getPMT()) * (__q->getT() - __p->getT());
333
334 hT->Fill(tHit_ms, dt, 1.0 / W);
335
336 }
337 }
338 }
339
340 }
341
342 p = (q - 1);
343
344 }
345
346 }
347 }
348 }
349
350 // cleanup
351
352 if (nextSummary != NULL) { delete nextSummary; }
353
354 // update pointers
355
356 curr = next;
357
358 }
359
360
361 //--------------------------------------
362 // Histogram processing
363 //--------------------------------------
364
365 NOTICE("Processing histograms." << endl);
366
367 for (JDetector::const_iterator module = detector.begin(); module != detector.end(); ++module) {
368
369 string moduleLabel = getLabel(*module);
370
371 if (RT_DOM.count(moduleLabel) && NC_DOM.count(moduleLabel)) {
372
373 TH1F* rt = RT_DOM.at(moduleLabel);
374 TH1F* nc = NC_DOM.at(moduleLabel);
375
376 for (int b = 1; b <= rt->GetXaxis()->GetNbins(); b++) {
377
378 double r = rt->GetBinContent(b);
379 double t = rt->GetBinCenter( b);
380
381 RT_DET["SUM"]->Fill(t, r);
382
383 }
384
385
386 for (int b = 1; b <= nc->GetXaxis()->GetNbins(); b++) {
387
388 double n = nc->GetBinContent(b);
389 double t = nc->GetBinCenter( b);
390
391 NC_DET["SUM"]->Fill(t, n);
392
393 }
394
395 } else {
396
397 DEBUG(moduleLabel << " not active." << endl);
398
399 }
400
401 }
402
403 NOTICE("Writing output file" << endl);
404
405 if (outputFile != "") {
406
407 TFile out(outputFile.c_str(), "RECREATE");
408
409 NOTICE("Writing 1D histograms" << endl);
410
411 RT_DET.Write(out);
412 NC_DET.Write(out);
413
414 NOTICE("Writing 2D histogram" << endl);
415
416 if (hT != NULL) {
417 hT->Write();
418 }
419
420 if (!globalOutputOnly) {
421
422 TString dir_tag = runTag + TString(".Modules");
423
424 NOTICE("Writing individual modules histograms" << endl);
425
426 TDirectory* dir = out.mkdir(dir_tag);
427 RT_DOM.Write(*dir);
428 NC_DOM.Write(*dir);
429 }
430
431 NOTICE("Closing file" << endl);
432
433 out.Close();
434 }
435
436}
string outputFile
#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.
Data structure for a composite optical module.
Definition JModule.hh:75
General exception.
Definition JException.hh:24
Utility class to parse command line options.
Definition JParser.hh:1698
Auxiliary class to manage set of compatible ROOT objects (e.g. histograms) using unique keys.
Definition JManager.hh:47
Router for fast addressing of summary data in KM3NETDAQ::JDAQSummaryslice data structure as a functio...
const JDAQSummaryFrame & getSummaryFrame(const JDAQModuleIdentifier &module) const
Get summary frame.
void update(const JDAQSummaryslice *ps)
Update router.
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.
int getFrameIndex() const
Get frame index.
bool testFIFOStatus() const
Test FIFO status.
bool testHighRateVeto() const
Test high-rate veto status.
Data storage class for rate measurements of all PMTs in one module.
std::string getLabel(const JLocation &location)
Get module label for monitoring and other applications.
Definition JLocation.hh:247
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.
std::set< JROOTClassSelector > getROOTClassSelection(const bool option=false)
Get ROOT class selection.
Long64_t counter_type
Type definition for counter.
Float_t h2d_bintype
Definition JRipple.hh:6
const unsigned int h2d_limit
Definition JRipple.hh:11
TH2F h2d_t
Definition JRipple.hh:7
const int n
Definition JPolint.hh:791
double getFrameTime()
Get frame time duration.
Definition JDAQClock.hh:162
double getTimeOfRTS(const JDAQChronometer &chronometer)
Get time of last RTS in ns since start of run for a given chronometer.
static const int NUMBER_OF_PMTS
Total number of PMTs in module.
Definition JDAQ.hh:26
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 to select DAQ hits based on time-over-treshold value.
@ join_t
join consecutive hits according match criterion