Jpp 19.3.0-rc.2
the software that should make you happy
Loading...
Searching...
No Matches
JTurbot2D.cc File Reference

Example program to search for out of sync shifts around integral timeslices evolving over time. More...

#include <string>
#include <iostream>
#include <iomanip>
#include <limits>
#include <map>
#include <vector>
#include "TROOT.h"
#include "TFile.h"
#include "TH1D.h"
#include "TH2D.h"
#include "TF1.h"
#include "JDAQ/JDAQEventIO.hh"
#include "JDAQ/JDAQTimesliceIO.hh"
#include "JDAQ/JDAQEvaluator.hh"
#include "JDetector/JDetector.hh"
#include "JDetector/JDetectorToolkit.hh"
#include "JDetector/JDAQHitRouter.hh"
#include "JDetector/JPMTParametersMap.hh"
#include "JTrigger/JBuildL0.hh"
#include "JTrigger/JBuildL1.hh"
#include "JLang/JSinglePointer.hh"
#include "JROOT/JManager.hh"
#include "JROOT/JROOTClassSelector.hh"
#include "JSupport/JSingleFileScanner.hh"
#include "JSupport/JTreeScannerInterface.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"

Go to the source code of this file.

Functions

int main (int argc, char **argv)
 

Detailed Description

Example program to search for out of sync shifts around integral timeslices evolving over time.

Author
shallmann

Definition in file JTurbot2D.cc.

Function Documentation

◆ main()

int main ( int argc,
char ** argv )

Definition at line 44 of file JTurbot2D.cc.

45{
46 using namespace std;
47 using namespace JPP;
48 using namespace KM3NETDAQ;
49
50 JSingleFileScanner<> inputFile;
51 JLimit_t& numberOfEvents = inputFile.getLimit();
52 string outputFile;
53 string detectorFile;
54 double TMax_ns;
55 int numberOfTimeslices;
56 double binWidth_ns;
57 double binWidth_min_timeEvolution;
58 JROOTClassSelector selector;
59 int debug;
60
61 try {
62
63 JParser<> zap("Example program to search for out of sync shifts around integral timeslices evolving over time.");
64
65 zap['f'] = make_field(inputFile);
66 zap['o'] = make_field(outputFile) = "";
67 zap['a'] = make_field(detectorFile);
68 zap['n'] = make_field(numberOfEvents) = JLimit::max();
69 zap['T'] = make_field(TMax_ns) = 20.0;
70 zap['N'] = make_field(numberOfTimeslices) = 40;
71 zap['W'] = make_field(binWidth_ns) = 10e3;
72 zap['M'] = make_field(binWidth_min_timeEvolution) = 2;
74 zap['d'] = make_field(debug) = 2;
75
76 zap(argc, argv);
77 }
78 catch(const exception& error) {
79 FATAL(error.what() << endl);
80 }
81
82
83
84 map<int, int> MASK; // mask PMTs with permament high-rate veto.
85
86 const int WR = 0x80000000; // White Rabbit
87
88 MASK[808969848] = 0x00000020; // floor 03, pmt 05
89 MASK[809544061] = 0x00000080; // floor 05, pmt 07
90 MASK[808432835] = 0x00004000; // floor 18, pmt 14
91
92
94
95 try {
96 load(detectorFile, detector);
97 }
98 catch(const JException& error) {
99 FATAL(error);
100 }
101
102 const JDAQHitRouter router(detector);
103
104
105 typedef double hit_type;
106
107 JBuildL0<hit_type> buildL0;
108 JBuildL1<hit_type> buildL1(JBuildL1Parameters(TMax_ns, true));
109 vector <hit_type> data;
110
112
114
116
117 if (selector == "") {
118
119 if ((ps = new JTreeScanner< JAssertConversion<JDAQTimesliceSN, JDAQTimeslice> >(inputFile))->getEntries() != 0 ||
120 (ps = new JTreeScanner< JAssertConversion<JDAQTimesliceL2, JDAQTimeslice> >(inputFile))->getEntries() != 0 ||
121 (ps = new JTreeScanner< JAssertConversion<JDAQTimesliceL1, JDAQTimeslice> >(inputFile))->getEntries() != 0 ||
122 (ps = new JTreeScanner< JAssertConversion<JDAQTimesliceL0, JDAQTimeslice> >(inputFile))->getEntries() != 0 ||
123 (ps = new JTreeScanner< JAssertConversion<JDAQTimeslice, JDAQTimeslice> >(inputFile))->getEntries() != 0) {
124 } else {
125 FATAL("No timeslice data." << endl);
126 }
127
128 NOTICE("Selected class " << ps->getClassName() << endl);
129
130 } else {
131
132 ps = zmap[selector];
133
134 ps->configure(inputFile);
135 }
136
137 typedef JManager<int, TH2D> JManager_t;
138
139 // shift in integer time slices
140 const Double_t ymin = -(numberOfTimeslices + 0.5);
141 const Double_t ymax = +(numberOfTimeslices + 0.5);
142 const Int_t ny = (Int_t) (ymax - ymin);
143 // time in M minute binning
144 const Double_t xmin = inputFile.getLimit().min()*getFrameTime()*1e-9;// time of first event - offset
145 const Double_t xmax = min(inputFile.getLimit().max(), ps->getEntries())*getFrameTime()*1e-9; // time of last event
146 const Int_t nx = (Int_t) ((xmax - xmin) / (60*binWidth_min_timeEvolution)) + 1;
147
148 JManager_t manager(new TH2D("M2D_%", ";time [s];shift [timeslices]", nx, xmin, xmax, ny, ymin, ymax));
149
150 typedef map<int, vector<double> > map_type; // frame index -> event times
151
152 map_type buffer;
153
154 for (JTreeScanner<JDAQEvent> in(inputFile.getFilename(), inputFile.getLimit()); in.hasNext(); ) {
155
156 STATUS("event: " << setw(10) << in.getCounter() << '\r'); DEBUG(endl);
157
158 JDAQEvent* event = in.next();
159
160 // Average time of triggered hits
161
162 double t0 = 0.0;
163
164 for (JDAQEvent::const_iterator<JDAQTriggeredHit> hit = event->begin<JDAQTriggeredHit>(); hit != event->end<JDAQTriggeredHit>(); ++hit) {
165 t0 += getTime(*hit, router.getPMT(*hit));
166 }
167
168 t0 /= event->size<JDAQTriggeredHit>();
169 t0 += event->getFrameIndex() * getFrameTime();
170
171 buffer[event->getFrameIndex()].push_back(t0);
172 }
173 STATUS(endl);
174
175
176 while (ps->hasNext()) {
177
178 STATUS("event: " << setw(10) << ps->getCounter() << '\r'); DEBUG(endl);
179
180 JDAQTimeslice* timeslice = ps->next();
181
182 map_type::const_iterator p = buffer.lower_bound(timeslice->getFrameIndex() - numberOfTimeslices);
183 map_type::const_iterator q = buffer.upper_bound(timeslice->getFrameIndex() + numberOfTimeslices);
184
185 int number_of_events = 0;
186
187 for (map_type::const_iterator i = p; i != q; ++i) {
188 number_of_events += i->second.size();
189 }
190
191 if (number_of_events == 0) {
192 continue;
193 }
194
195 for (JDAQTimeslice::const_iterator frame = timeslice->begin(); frame != timeslice->end(); ++frame) {
196
197 // Test that none of the PMTs has high-rate veto.
198
199 if ((frame->getStatus() & ~MASK[frame->getModuleID()] & ~WR) == 0) {
200
201 data.clear();
202
203 buildL1(*frame, router.getModule(frame->getModuleID()), back_inserter(data));
204
205 TH2D* h2 = manager[frame->getModuleID()];
206
207 for (vector<hit_type>::const_iterator hit = data.begin(); hit != data.end(); ++hit) {
208
209 const double t1 = *hit + frame->getFrameIndex() * getFrameTime();
210
211 for (map_type::const_iterator i = p; i != q; ++i) {
212 for (map_type::mapped_type::const_iterator j = i->second.begin(); j != i->second.end(); ++j) {
213
214 const double t0 = *j;
215
216 // fill only shift is near integer time slice
217 const double timeFromInt = (t1 - t0)/getFrameTime() - round((t1 - t0)/getFrameTime());
218 if ( abs(timeFromInt)*getFrameTime() < binWidth_ns/2. ){
219 h2->Fill(timeslice->getFrameIndex() * getFrameTime() * 1e-9, round(t1 - t0)/getFrameTime());
220 }
221 }
222 }
223 }
224 }
225 }
226 }
227 STATUS(endl);
228
229 if (outputFile != "") {
230 manager.Write(outputFile.c_str());
231 }
232}
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
Simple wrapper around JModuleRouter class for direct addressing of PMT data in detector data structur...
Detector data structure.
Definition JDetector.hh:96
Auxialiary class to assert type conversion.
General exception.
Definition JException.hh:24
The template JSinglePointer class can be used to hold a pointer to an object.
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
Template definition for direct access of elements in ROOT TChain.
Template L0 hit builder.
Definition JBuildL0.hh:38
Template L1 hit builder.
Definition JBuildL1.hh:90
int getFrameIndex() const
Get frame index.
Template const_iterator.
Definition JDAQEvent.hh:68
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.
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
std::set< JROOTClassSelector > getROOTClassSelection(const bool option=false)
Get ROOT class selection.
const char * getTime()
Get current local time conform ISO-8601 standard.
int j
Definition JPolint.hh:801
KM3NeT DAQ data structures and auxiliaries.
Definition DataQueue.cc:39
double getFrameTime()
Get frame time duration.
Definition JDAQClock.hh:162
std::map< int, range_type > map_type
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 data structure for L1 build parameters.
Definition JBuildL1.hh:38