Jpp 20.0.0-27-g39925593c-D
the software that should make you happy
Loading...
Searching...
No Matches
JTurbot2D.cc
Go to the documentation of this file.
1
2#include <string>
3#include <iostream>
4#include <iomanip>
5#include <limits>
6#include <map>
7#include <vector>
8
9#include "TROOT.h"
10#include "TFile.h"
11#include "TH1D.h"
12#include "TH2D.h"
13#include "TF1.h"
14
15#include "JDAQ/JDAQEventIO.hh"
17#include "JDAQ/JDAQEvaluator.hh"
22#include "JTrigger/JBuildL0.hh"
23#include "JTrigger/JBuildL1.hh"
25#include "JROOT/JManager.hh"
31#include "JSupport/JSupport.hh"
32
33#include "Jeep/JPrint.hh"
34#include "Jeep/JParser.hh"
35#include "Jeep/JMessage.hh"
36
37
38/**
39 * \file
40 *
41 * Example program to search for out of sync shifts around integral timeslices evolving over time.
42 * \author shallmann
43 */
44int main(int argc, char **argv)
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;
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) = "";
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;
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 {
97 }
98 catch(const JException& error) {
99 FATAL(error);
100 }
101
102 const JDAQHitRouter router(detector);
103
104
105 typedef double hit_type;
106
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
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());
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}
Direct access to PMT data in detector data structure for DAQ hits.
string outputFile
Data structure for detector geometry and calibration.
Dynamic ROOT object management.
General purpose messaging.
#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
Utility class to parse command line options.
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition JParser.hh:2142
I/O formatting auxiliaries.
Scanning of objects from a single file according a format that follows from the extension of each fil...
ROOT TTree parameter settings of various packages.
int main(int argc, char **argv)
Definition JTurbot2D.cc:44
Simple wrapper around JModuleRouter class for direct addressing of PMT data in detector data structur...
const JModule & getModule(const JDAQKeyHit &hit) const
Get module parameters.
const JPMT & getPMT(const JDAQKeyHit &hit) const
Get PMT parameters.
Detector data structure.
Definition JDetector.hh:96
General exception.
Definition JException.hh:24
Template definition of a multi-dimensional oscillation probability interpolation table.
Template definition for direct access of elements in ROOT TChain.
int getFrameIndex() const
Get frame index.
Template const_iterator.
Definition JDAQEvent.hh:68
JTriggerCounter_t next()
Increment trigger counter.
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).
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
Transmission with position.
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 data structure for L1 build parameters.
Definition JBuildL1.hh:38