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

Example program to test JTRIGGER::JBuildL1 and JTRIGGER::JBuildL2 hit coincidence building with Monte Carlo events. More...

#include <string>
#include <iostream>
#include <iomanip>
#include <vector>
#include "TROOT.h"
#include "TFile.h"
#include "TProfile.h"
#include "km3net-dataformat/offline/Head.hh"
#include "km3net-dataformat/offline/Evt.hh"
#include "km3net-dataformat/offline/Hit.hh"
#include "JDAQ/JDAQTimesliceIO.hh"
#include "JAAnet/JAAnetToolkit.hh"
#include "JTimeslice/JEventTimeslice.hh"
#include "JDetector/JDetector.hh"
#include "JDetector/JDetectorToolkit.hh"
#include "JDetector/JModuleRouter.hh"
#include "JDetector/JDetectorSimulator.hh"
#include "JDetector/JPMTParametersMap.hh"
#include "JDetector/JPMTDefaultSimulator.hh"
#include "JDetector/JCLBDefaultSimulator.hh"
#include "JTrigger/JHit.hh"
#include "JTrigger/JFrame.hh"
#include "JTrigger/JTimeslice.hh"
#include "JTrigger/JSuperFrame2D.hh"
#include "JTrigger/JHitL0.hh"
#include "JTrigger/JHitL1.hh"
#include "JTrigger/JBuildL1.hh"
#include "JTrigger/JBuildL2.hh"
#include "JSupport/JMultipleFileScanner.hh"
#include "JSupport/JMonteCarloFileSupportkit.hh"
#include "JSupport/JSupport.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 test JTRIGGER::JBuildL1 and JTRIGGER::JBuildL2 hit coincidence building with Monte Carlo events.

Author
mdejong

Definition in file JSignalL1.cc.

Function Documentation

◆ main()

int main ( int argc,
char ** argv )

Definition at line 50 of file JSignalL1.cc.

51{
52 using namespace std;
53 using namespace JPP;
54
56 JLimit_t& numberOfEvents = inputFile.getLimit();
57 string outputFile;
58 string detectorFile;
59 JPMTParametersMap pmtParameters;
60 double Tmax_ns;
61 int debug;
62
63 try {
64
65 JParser<> zap("Example program to test hit coincidence building with Monte Carlo events.");
66
67 zap['f'] = make_field(inputFile);
68 zap['o'] = make_field(outputFile) = "buildL2.root";
69 zap['n'] = make_field(numberOfEvents) = JLimit::max();
70 zap['a'] = make_field(detectorFile);
71 zap['P'] = make_field(pmtParameters) = JPARSER::initialised();
72 zap['T'] = make_field(Tmax_ns) = 20.0; // [ns]
73 zap['d'] = make_field(debug) = 0;
74
75 zap(argc, argv);
76 }
77 catch(const exception &error) {
78 FATAL(error.what() << endl);
79 }
80
81
82 using namespace KM3NETDAQ;
83
84
86
87 try {
88 load(detectorFile, detector);
89 }
90 catch(const JException& error) {
91 FATAL(error);
92 }
93
95
96 simbad.reset(new JPMTDefaultSimulator(pmtParameters, detector));
97 simbad.reset(new JCLBDefaultSimulator());
98
99
100 TFile out(outputFile.c_str(), "recreate");
101
102 TProfile hn("hn", NULL, 31, 0.5, +31.5);
103 TProfile hc("hc", NULL, 21, -1.05, +1.05);
104 TProfile ht("ht", NULL, 20, 0.5, +20.5);
105
106
107 const JModuleRouter moduleRouter(detector);
108
109 typedef double hit_type;
110 typedef vector <hit_type> JFrameL1_t;
111 typedef JSuperFrame2D<hit_type> JSuperFrame2D_t;
112 typedef JBuildL1 <hit_type> JBuildL1_t;
113
114
115 const JBuildL1_t buildL1(JBuildL1Parameters((hit_type) Tmax_ns, true));
116 JSuperFrame2D_t buffer;
117
118
119 for (JMultipleFileScanner<Evt>& in = inputFile; in.hasNext(); ) {
120
121 STATUS("event: " << setw(10) << in.getCounter() << '\r'); DEBUG(endl);
122
123 Evt* event = in.next();
124
125 const int frame_index = 1;
126
127 event->mc_t = 0.5 * getFrameTime();
128
129 const JDAQChronometer chronometer(detector.getID(), 1, frame_index, JDAQUTCExtended(getTimeOfFrame(frame_index)));
130
131 JEventTimeslice timeslice(chronometer, simbad, *event);
132
133 for (JDAQTimeslice::const_iterator super_frame = timeslice.begin(); super_frame != timeslice.end(); ++super_frame) {
134
135 if (moduleRouter.hasModule(super_frame->getModuleID())) {
136
137 // calibration
138
139 const JModule& module = detector.getModule(moduleRouter.getAddress(super_frame->getModuleID()));
140
141 buffer(*super_frame, module);
142
143 JFrameL1_t dataL1;
144
145 buildL1(buffer, back_inserter(dataL1));
146
147 if (!dataL1.empty()) {
148
149 JBuildL2<hit_type> buildL2(JL2Parameters(1, Tmax_ns, -1.0));
150
151 JFrameL1_t d1(dataL1);
152 JFrameL1_t d2;
153
154 for (int i = 1; i <= hn.GetNbinsX(); ++i) {
155
156 buildL2.numberOfHits = (int) hn.GetBinCenter(i);
157
158 d2.clear();
159
160 buildL2(buffer, d1, back_inserter(d2));
161
162 hn.Fill((double) buildL2.numberOfHits, (double) d2.size() / (double) dataL1.size());
163
164 d1.swap(d2);
165 }
166 }
167
168
169 if (!dataL1.empty()) {
170
171 JBuildL2<hit_type> buildL2(JL2Parameters(2, Tmax_ns, -1.0));
172
173 JFrameL1_t d1(dataL1);
174 JFrameL1_t d2;
175
176 for (int i = 1; i <= hc.GetNbinsX(); ++i) {
177
178 buildL2.ctMin = hc.GetBinCenter(i);
179
180 d2.clear();
181
182 buildL2(buffer, d1, back_inserter(d2));
183
184 hc.Fill(buildL2.ctMin, (double) d2.size() / (double) dataL1.size());
185
186 d1.swap(d2);
187 }
188 }
189
190
191 if (!dataL1.empty()) {
192
193 JBuildL2<hit_type> buildL2(JL2Parameters(2, 0.0, -1.0));
194
195 JFrameL1_t d1(dataL1);
196 JFrameL1_t d2;
197
198 for (int i = ht.GetNbinsX(); i != 0; --i) {
199
200 buildL2.TMaxLocal_ns = ht.GetBinCenter(i);
201
202 d2.clear();
203
204 buildL2(buffer, d1, back_inserter(d2));
205
206 ht.Fill(buildL2.TMaxLocal_ns, (double) d2.size() / (double) dataL1.size());
207
208 d1.swap(d2);
209 }
210 }
211 }
212 }
213 }
214 NOTICE(endl);
215
216 out.Write();
217 out.Close();
218}
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
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
Auxiliary class for map of PMT parameters.
General exception.
Definition JException.hh:24
Utility class to parse command line options.
Definition JParser.hh:1698
General purpose class for object reading from a list of file names.
virtual bool hasNext() override
Check availability of next element.
Template L2 builder.
Definition JBuildL2.hh:49
2-dimensional frame with time calibrated data from one optical module.
Data structure for UTC time.
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).
KM3NeT DAQ data structures and auxiliaries.
Definition DataQueue.cc:39
double getFrameTime()
Get frame time duration.
Definition JDAQClock.hh:162
double getTimeOfFrame(const int frame_index)
Get start time of frame in ns since start of run for a given frame index.
Definition JDAQClock.hh:185
The Evt class respresent a Monte Carlo (MC) event as well as an offline event.
Definition Evt.hh:21
double mc_t
MC: time where the mc-event was put in the timeslice, since start of run (offset+frameidx*timeslice_d...
Definition Evt.hh:47
Detector file.
Definition JHead.hh:227
Acoustic hit.
Definition JBillabong.cc:70
Empty structure for specification of parser element that is initialised (i.e. does not require input)...
Definition JParser.hh:68
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
Data structure for L2 parameters.
Timeslice with Monte Carlo event.