Jpp
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
JSignalL1.cc
Go to the documentation of this file.
1 
2 #include <string>
3 #include <iostream>
4 #include <iomanip>
5 #include <vector>
6 
7 #include "TROOT.h"
8 #include "TFile.h"
9 #include "TProfile.h"
10 
11 #include "evt/Head.hh"
12 #include "evt/Evt.hh"
13 #include "evt/Hit.hh"
14 #include "JDAQ/JDAQTimeslice.hh"
15 
16 #include "JAAnet/JAAnetToolkit.hh"
18 
19 #include "JDetector/JDetector.hh"
26 
27 #include "JTrigger/JHit.hh"
28 #include "JTrigger/JFrame.hh"
29 #include "JTrigger/JTimeslice.hh"
31 #include "JTrigger/JHitL0.hh"
32 #include "JTrigger/JHitL1.hh"
33 #include "JTrigger/JBuildL1.hh"
34 #include "JTrigger/JBuildL2.hh"
35 
38 #include "JSupport/JSupport.hh"
39 
40 #include "Jeep/JParser.hh"
41 #include "Jeep/JMessage.hh"
42 
43 
44 /**
45  * \file
46  *
47  * Example program to test JTRIGGER::JBuildL1 and JTRIGGER::JBuildL2 hit coincidence building with Monte Carlo events.
48  * \author mdejong
49  */
50 int main(int argc, char **argv)
51 {
52  using namespace std;
53  using namespace JPP;
54 
55  JMultipleFileScanner<Evt> inputFile;
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  cout.tie(&cerr);
85 
86 
87  JDetector detector;
88 
89  try {
90  load(detectorFile, detector);
91  }
92  catch(const JException& error) {
93  FATAL(error);
94  }
95 
96  JDetectorSimulator simbad(detector);
97 
98  simbad.reset(new JPMTDefaultSimulator(pmtParameters, detector));
99  simbad.reset(new JCLBDefaultSimulator());
100 
101 
102  TFile out(outputFile.c_str(), "recreate");
103 
104  TProfile hn("hn", NULL, 31, 0.5, +31.5);
105  TProfile hc("hc", NULL, 21, -1.05, +1.05);
106  TProfile ht("ht", NULL, 20, 0.5, +20.5);
107 
108 
109  const JModuleRouter moduleRouter(detector);
110 
111  typedef double hit_type;
113  typedef JSuperFrame2D<hit_type> JSuperFrame2D_t;
114  typedef JBuildL1 <hit_type> JBuildL1_t;
115 
116 
117  const JBuildL1_t buildL1((hit_type) Tmax_ns, true);
118  JSuperFrame2D_t buffer;
119 
120 
121  for (JMultipleFileScanner<Evt>& in = inputFile; in.hasNext(); ) {
122 
123  STATUS("event: " << setw(10) << in.getCounter() << '\r'); DEBUG(endl);
124 
125  Evt* event = in.next();
126 
127  const int frame_index = 1;
128 
129  event->mc_t = 0.5 * getFrameTime();
130 
131  const JDAQChronometer chronometer(detector.getID(), 1, frame_index, JDAQUTCExtended(getTimeOfFrame(frame_index)));
132 
133  JEventTimeslice timeslice(chronometer, simbad, *event);
134 
135  for (JDAQTimeslice::const_iterator super_frame = timeslice.begin(); super_frame != timeslice.end(); ++super_frame) {
136 
137  if (moduleRouter.hasModule(super_frame->getModuleID())) {
138 
139  // calibration
140 
141  const JModule& module = detector.getModule(moduleRouter.getAddress(super_frame->getModuleID()));
142 
143  buffer(*super_frame, module);
144 
145  JFrameL1_t dataL1;
146 
147  buildL1(buffer, back_inserter(dataL1));
148 
149  if (!dataL1.empty()) {
150 
151  JBuildL2<hit_type> buildL2(1, Tmax_ns, -1.0);
152 
153  JFrameL1_t d1(dataL1);
154  JFrameL1_t d2;
155 
156  for (int i = 1; i <= hn.GetNbinsX(); ++i) {
157 
158  buildL2.numberOfHits = (int) hn.GetBinCenter(i);
159 
160  d2.clear();
161 
162  buildL2(buffer, d1, back_inserter(d2));
163 
164  hn.Fill((double) buildL2.numberOfHits, (double) d2.size() / (double) dataL1.size());
165 
166  d1.swap(d2);
167  }
168  }
169 
170 
171  if (!dataL1.empty()) {
172 
173  JBuildL2<hit_type> buildL2(2, Tmax_ns, -1.0);
174 
175  JFrameL1_t d1(dataL1);
176  JFrameL1_t d2;
177 
178  for (int i = 1; i <= hc.GetNbinsX(); ++i) {
179 
180  buildL2.ctMin = hc.GetBinCenter(i);
181 
182  d2.clear();
183 
184  buildL2(buffer, d1, back_inserter(d2));
185 
186  hc.Fill(buildL2.ctMin, (double) d2.size() / (double) dataL1.size());
187 
188  d1.swap(d2);
189  }
190  }
191 
192 
193  if (!dataL1.empty()) {
194 
195  JBuildL2<hit_type> buildL2(2, 0.0, -1.0);
196 
197  JFrameL1_t d1(dataL1);
198  JFrameL1_t d2;
199 
200  for (int i = ht.GetNbinsX(); i != 0; --i) {
201 
202  buildL2.TMaxLocal_ns = ht.GetBinCenter(i);
203 
204  d2.clear();
205 
206  buildL2(buffer, d1, back_inserter(d2));
207 
208  ht.Fill(buildL2.TMaxLocal_ns, (double) d2.size() / (double) dataL1.size());
209 
210  d1.swap(d2);
211  }
212  }
213  }
214  }
215  }
216  NOTICE(endl);
217 
218  out.Write();
219  out.Close();
220 }
Utility class to parse command line options.
Definition: JParser.hh:1410
Timeslice with Monte Carlo event.
#define STATUS(A)
Definition: JMessage.hh:61
Structure to store the ToT mean and standard deviation of the hits produced by a nanobeacon in a sour...
Empty structure for specification of parser element that is initialised (i.e.
Definition: JParser.hh:64
JSuperFrame2D< hit_type > JSuperFrame2D_t
Definition: JDataFilter.cc:93
string outputFile
Data structure for UTC time.
Data structure for detector geometry and calibration.
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
Basic data structure for L0 hit.
JLimit JLimit_t
Type definition of limit.
Definition: JLimit.hh:214
JBuildL1< hit_type > JBuildL1_t
Definition: JDataFilter.cc:95
Basic data structure for time and time over threshold information of hit.
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1836
double getFrameTime()
Get frame time duration.
Definition: JDAQClock.hh:162
Definition of hit and track types and auxiliary methods for handling Monte Carlo data.
#define NOTICE(A)
Definition: JMessage.hh:62
void load(const JString &file_name, JDetector &detector)
Load detector from input file.
int debug
debug level
Definition: JSirene.cc:59
General purpose messaging.
#define FATAL(A)
Definition: JMessage.hh:65
Scanning of objects from multiple files according a format that follows from the extension of each fi...
Direct access to module in detector data structure.
vector< hit_type > JFrameL1_t
Definition: JDataFilter.cc:91
Utility class to parse command line options.
ROOT TTree parameter settings.
2-dimensional frame with time calibrated data from one optical module.
const JLimit & getLimit() const
Get limit.
Definition: JLimit.hh:72
Basic data structure for L1 hit.
double hit_type
Definition: JDataFilter.cc:89
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:60
int main(int argc, char *argv[])
Definition: Main.cpp:15