Jpp 20.0.0-27-g39925593c-D
the software that should make you happy
Loading...
Searching...
No Matches
Functions
JKaboom.cc File Reference
#include <string>
#include <iostream>
#include <fstream>
#include <iomanip>
#include <vector>
#include "TROOT.h"
#include "TFile.h"
#include "TTree.h"
#include "JDAQ/JDAQEvaluator.hh"
#include "km3net-dataformat/offline/Evt.hh"
#include "JDetector/JDetector.hh"
#include "JDetector/JDetectorToolkit.hh"
#include "JROOT/JROOTClassSelector.hh"
#include "JSupport/JMultipleFileScanner.hh"
#include "JSupport/JTreeScanner.hh"
#include "JSupport/JAutoTreeScanner.hh"
#include "JSupport/JSupport.hh"
#include "Jeep/JPrint.hh"
#include "Jeep/JParser.hh"
#include "Jeep/JMessage.hh"
#include "JTrigger/JDAQHitToTSelector.hh"
#include "JSupernova.hh"

Go to the source code of this file.

Functions

int initialize_root_tree (JCoincidenceSN &obs, int *frame_number, TTree *tr)
 Initialize the output ROOT file with CCSN candidate information.
 
int main (int argc, char **argv)
 

Detailed Description

Authors
selhedri,iagoos

Application to analyse buffered L1 timeslices for the online CCSN analysis

Definition in file JKaboom.cc.

Function Documentation

◆ initialize_root_tree()

int initialize_root_tree ( JCoincidenceSN & obs,
int * frame_number,
TTree * tr )
inline

Initialize the output ROOT file with CCSN candidate information.

Parameters
obsSingle-DOM observables
frame_numberTimeslice frame number
trOutput ROOT tree

Definition at line 50 of file JKaboom.cc.

50 {
51 tr->Branch("dom_id", &obs.moduleID, "dom_id/I");
52 tr->Branch("vetoIndex", &obs.vetoIndex, "vetoIndex/I");
53 tr->Branch("event_time", &obs.time, "event_time/D");
54 tr->Branch("timeslice_time_s", &obs.timeslice_time, "timeslice_time_s/D");
55 tr->Branch("frame_number", frame_number, "frame_number/I");
56 tr->Branch("multiplicity", &(obs.multiplicity), "multiplicity/I");
57 return 0;
58}
Template definition of a multi-dimensional oscillation probability interpolation table.

◆ main()

int main ( int argc,
char ** argv )

Definition at line 60 of file JKaboom.cc.

61{
62
63 typedef JRange<int> JRange_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;
74 JROOTClassSelector selector;
75 int debug;
77
78 // Check input parameters
79
80 try {
81
82 JParser<> zap("Example application to study supernova detection background");
83
84 zap['f'] = make_field(inputFile);
85 zap['o'] = make_field(outputFile) = "kaboom.root";
86 zap['n'] = make_field(numberOfEvents) = JLimit::max();
88 zap['T'] = make_field(TMax_ns) = 10.0;
89 zap['V'] = make_field(TVeto_ns) = 1000.0;
91 zap['d'] = make_field(debug) = 1;
94
95 zap(argc, argv);
96 }
97 catch(const exception &error) {
98 FATAL(error.what() << endl);
99 }
100
101 cout.tie(&cerr);
102
104
105 try {
107 }
108 catch(const JException& error) {
109 FATAL(error);
110 }
111
112 TTree *tr = new TTree("singledom", "SN coincidence info");
113 JCoincidenceSN SNdata(0,0,0);
114 int frame_number;
116
117 // -----------------------------------
118 // STEP 0: prepare
119 // -----------------------------------
120
121 // configure input streams
123
125
127
128 pts->configure(inputFile);
129
130 // detect input stream size
131
132
133 int fEnd = pts->rbegin()->getFrameIndex();
134 int fStart = pts->begin( )->getFrameIndex();
135
136
137 // crop histogram if -n is specified
138 if (fEnd > inputFile.getUpperLimit()) {
139 fEnd = fStart + inputFile.getUpperLimit();
140 }
141
142 int fLength = 1 + fEnd - fStart;
143
144 NOTICE("begin | end | length = " << fStart << " | " << fEnd << " | " << fLength << endl);
145
146 // detector routers
147
148 const JModuleRouter moduleRouter(detector);
149
151
152
153 // -----------------------------------
154 // STEP 1: building vetoes from events
155 // -----------------------------------
156
157 int runNumber = 0;
158
159 TH1D* h_vtr = new TH1D("VetoTimeRange","VetoTimeRange", 10000, 0, 10000);
160
162
164
165 for (; evIn.hasNext(); ) {
166
167 STATUS("event: " << setw(10) << evIn.getCounter() << '\r'); DEBUG(endl);
168
169 JDAQEvent* event = evIn.next();
170
171 if (!runNumber) { runNumber = event->getRunNumber(); }
172
173 JVeto vp(*event, hitRouter);
174
175 triggeredEvents[event->getFrameIndex()].push_back(vp);
176
177 h_vtr->Fill(vp.getLength());
178
179 }
180
181 STATUS(triggeredEvents.size() << " JDAQEvents loaded in veto buffer." << endl);
182
183 //-----------------------------
184 // STEP 2: timeslice processing
185 // ----------------------------
186
187 counter_type counter = 0;
188
189 TFile *out = new TFile(outputFile.c_str(), "RECREATE");
190 tr->SetDirectory(out);
191
192 for ( ; pts->hasNext() && counter != inputFile.getLimit(); ++counter) {
193
194 STATUS("timeslice: " << setw(10) << counter << '\r'); DEBUG(endl);
195
196 const JDAQTimeslice* timeslice = pts->next();
197
198 frame_number = timeslice->getFrameIndex();
199
200 if (counter == 0) {
201 NOTICE("Start frame index = " << frame_number << " @ T " << timeslice->getTimesliceStart() << endl);
202 }
203
205
206 if (triggeredEvents.count(frame_number)) {
208 }
209
210 // STANDARD PROCESSING
211
212 JDataSN preTrigger(TMax_ns, preTriggerThreshold);// totSelector_ns);
213
214 preTrigger(timeslice, moduleRouter, totSelector_ns);
215
217
219
220 // generate and configure event-based veto
221
222 // count trigger
223
224 JRange_t A = JRange_t(preTriggerThreshold, NUMBER_OF_PMTS);
225
226 // select single-module clusters with max multiplicity in range A
227 JSNFilterM trgA1(A, 1);
228
229 // select non-vetoed clusters with max multiplicity in range A
231
233 double event_time = tstart.getTimeNanoSecond() * 1e-9;
234
235 for(auto &p: trigger){
236 SNdata = p.getPeak();
237 SNdata.vetoIndex = trgA1(p) ? (trgAV(p) ? 0 : 2) : (trgAV(p) ? 1 : 3); // 0 if no veto, 1 if only DOM correlation veto, 2 if only muon veto, 3 if both
238 SNdata.timeslice_time = event_time;
239 tr->Fill();
240 }
241 }
242
243 //-------------------------------------
244 // STEP 3: Writing data
245 // ------------------------------------
246
247 tr->Write();
248 out->Close();
249}
string outputFile
int initialize_root_tree(JCoincidenceSN &obs, int *frame_number, TTree *tr)
Initialize the output ROOT file with CCSN candidate information.
Definition JKaboom.cc:50
#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
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
Auxiliary class to store reduced information of a coincidence on an optical module This class allows ...
Definition JSupernova.hh:52
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.
JDAQUTCExtended getTimesliceStart() const
Get start of timeslice.
int getFrameIndex() const
Get frame index.
Data structure for UTC time.
void load(const std::string &file_name, JDetector &detector)
Load detector from input file.
Type definition of range.
Definition JHead.hh:43
Detector file.
Definition JHead.hh:227
Auxiliary class to select ROOT class based on 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.