Jpp
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
JKexing.cc
Go to the documentation of this file.
1 #include <string>
2 #include <iostream>
3 #include <fstream>
4 #include <iomanip>
5 #include <vector>
6 
7 #include "TROOT.h"
8 #include "TFile.h"
9 #include "TH1D.h"
10 #include "TF1.h"
11 
12 #include "JSupernova.hh"
13 
14 #include "JDAQ/JDAQClock.hh"
15 #include "JDAQ/JDAQEvaluator.hh"
16 
17 #include "JDetector/JDetector.hh"
19 
20 #include "JGizmo/JManager.hh"
21 
23 #include "JTools/JQuantile.hh"
24 
26 #include "JSupport/JTreeScanner.hh"
28 #include "JSupport/JSupport.hh"
29 
31 #include "JMath/JMathToolkit.hh"
32 
33 #include "Jeep/JPrint.hh"
34 #include "Jeep/JParser.hh"
35 #include "Jeep/JMessage.hh"
36 
37 using namespace std;
38 using namespace JPP;
39 using namespace KM3NETDAQ;
40 using namespace JSUPERNOVA;
41 
42 /**
43  * \author mlincett
44  * \file
45  * Example application to test supernova trigger on run files.
46  */
47 int main(int argc, char **argv)
48 {
49 
50  typedef JRange<int> JRange_t;
51 
52  JMultipleFileScanner<> inputFile;
53  JLimit_t& numberOfEvents = inputFile.getLimit();
54  string outputFile;
55  string detectorFile;
56  double TMax_ns;
57  double TVeto_ns;
58  JROOTClassSelector selector;
59  int debug;
60  int preTriggerThreshold;
61  JRange_t M;
62 
63 
64  // Check input parameters
65 
66  try {
67 
68  JParser<> zap("Example program test supernova triggers.");
69 
70  zap['f'] = make_field(inputFile);
71  zap['o'] = make_field(outputFile) = "jsn.root";
72  zap['n'] = make_field(numberOfEvents) = JLimit::max();
73  zap['a'] = make_field(detectorFile);
74  zap['T'] = make_field(TMax_ns) = 10.0;
75  zap['V'] = make_field(TVeto_ns) = 1000.0;
76  zap['C'] = make_field(selector) = getROOTClassSelection<JDAQTimesliceTypes_t>();
77  zap['d'] = make_field(debug) = 1;
78  zap['M'] = make_field(M) = JRange<int>(6,10);
79  zap['t'] = make_field(preTriggerThreshold) = 4;
80 
81  zap(argc, argv);
82  }
83  catch(const exception &error) {
84  FATAL(error.what() << endl);
85  }
86 
87  cout.tie(&cerr);
88 
90 
91  try {
92  load(detectorFile, detector);
93  }
94  catch(const JException& error) {
95  FATAL(error);
96  }
97 
98  // Configure input streams
99 
101 
103 
104  JTreeScannerInterface<JDAQTimeslice>* pts = zmap[selector];
105 
106  pts->configure(inputFile);
107 
108  // Configure routers
109 
110  const JModuleRouter moduleRouter(detector);
111 
112  const JDAQHitRouter hitRouter(detector);
113 
114  // -----------------------------------
115  // STEP 1: building vetoes from events
116  // -----------------------------------
117 
118  TH1D* h_vtr = new TH1D("VetoTimeRange","VetoTimeRange", 10000, 0, 10000);
119 
121 
122  map<int, JVetoSet > triggeredEvents;
123 
124  for (; evIn.hasNext(); ) {
125 
126  STATUS("event: " << setw(10) << evIn.getCounter() << '\r'); DEBUG(endl);
127 
128  JDAQEvent* event = evIn.next();
129 
130  JVeto vp(*event, hitRouter);
131 
132  triggeredEvents[event->getFrameIndex()].push_back(vp);
133 
134  h_vtr->Fill(vp.getLength());
135 
136  }
137 
138  STATUS(triggeredEvents.size() << " JDAQEvents loaded in veto buffer." << endl);
139 
140  //--------------------------------
141  // STEP 2: timeslice processing
142  // -------------------------------
143 
144  // 0 -> no filter
145  // 1 -> count once per track
146  // 2 -> suppress track
147  // 3 -> suppress based on trigger veto
148 
149  // Output data structures
150 
151  typedef JManager<int, TH1D> JManager_t;
152 
153  JManager_t SNTH(new TH1D("SNT_F%", NULL, 100, 0.0, 100));
154 
155  JManager_t MUL(new TH1D("MUL_F%", NULL, 1 + 31, -0.5, 31 + 1 - 0.5));
156 
157  const int nStages = 5;
158  vector<vector<int> > trgHistory(nStages);
159 
160  // Loop
161 
162  counter_type counter = 0;
163 
164  for ( ; pts->hasNext() && counter != inputFile.getLimit(); ++counter) {
165 
166  STATUS("timeslice: " << setw(10) << counter << '\r'); DEBUG(endl);
167 
168  const JDAQTimeslice* timeslice = pts->next();
169 
170  int fIndex = timeslice->getFrameIndex();
171 
172  JDataSN preTrigger(TMax_ns, preTriggerThreshold);
173 
174  preTrigger(timeslice, moduleRouter);
175 
176  JTriggerSN trigger(M, TVeto_ns);
177 
178  trigger(preTrigger);
179 
180  // generate and configure event-based veto
181 
182  JVetoSet veto;
183 
184  if (triggeredEvents.count(fIndex)) {
185  veto = triggeredEvents.at(fIndex);
186  }
187 
188  trigger.setVeto(veto);
189 
190  // count trigger
191 
192  JSNFilterM trgOp(M);
193 
194  int rawCount = count_if(preTrigger.begin(), preTrigger.end(), trgOp);
195 
196  vector<int> trgCount = trigger.benchmark();
197 
198  trgHistory[0].push_back(rawCount);
199 
200  for (int i = 0; i < 3; i++) {
201  trgHistory[i + 1].push_back(trgCount[i]);
202  }
203 
204  trgHistory[4].push_back(trigger.getModules().size());
205 
206  }
207 
208  //-------------------------------------
209  // STEP 3: generating trigger summaries
210  // ------------------------------------
211 
212  for (int i = 0; i < nStages; i++) {
213  for (unsigned j = 0; j < trgHistory[i].size(); j++) {
214  SNTH[i]->Fill(trgHistory[i][j]);
215  }
216  }
217 
218  if (outputFile != "") {
219 
220  TFile out(outputFile.c_str(), "RECREATE");
221 
222  h_vtr->Write();
223  SNTH.Write(out);
224  MUL.Write(out);
225 
226  out.Close();
227 
228  }
229 
230 }
231 
232 
double getLength()
Get length of veto time range.
Definition: JSupernova.hh:102
Utility class to parse command line options.
Definition: JParser.hh:1410
General exception.
Definition: JException.hh:40
Auxiliary methods for geometrical methods.
JModuleSet getModules(JRange< int > A=JRange< int >(6, 10))
Get triggered modules after track veto.
Definition: JSupernova.hh:472
Auxiliary class to define a veto time window on a set of optical modules.
Definition: JSupernova.hh:67
#define STATUS(A)
Definition: JMessage.hh:61
Detector data structure.
Definition: JDetector.hh:77
Auxiliary class to select ROOT class based on class name.
Router for direct addressing of module data in detector data structure.
Structure to store the ToT mean and standard deviation of the hits produced by a nanobeacon in a sour...
Long64_t counter_type
Type definition for counter.
Definition: JCounter.hh:24
Dynamic ROOT object management.
Auxiliary class for a type holder.
Definition: JType.hh:19
string outputFile
Template definition for direct access of elements in ROOT TChain.
Definition: JTreeScanner.hh:89
Data structure for detector geometry and calibration.
Auxiliary interface for direct access of elements in ROOT TChain.
Simple wrapper around JModuleRouter class for direct addressing of PMT data in detector data structur...
int getFrameIndex() const
Get frame index.
Auxiliary class to manage set of compatible ROOT objects (e.g.
Definition: JManager.hh:38
Auxiliary class for defining the range of iterations of objects.
Definition: JLimit.hh:41
I/O formatting auxiliaries.
void setVeto(JVetoSet &vt)
Definition: JSupernova.hh:408
Detector file.
Definition: JHead.hh:126
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1836
Auxiliary class to apply the supernova trigger to SN data.
Definition: JSupernova.hh:385
Data time slice.
void load(const JString &file_name, JDetector &detector)
Load detector from input file.
SN filter based on multiplicity selection optional suppression of multi-module coincidences.
Definition: JSupernova.hh:309
int debug
debug level
Definition: JSirene.cc:59
General purpose messaging.
vector< int > benchmark(JRange< int > A=JRange< int >(6, 10))
Benchmark different trigger steps.
Definition: JSupernova.hh:452
Auxiliary class to select JTreeScanner based on ROOT class name.
#define FATAL(A)
Definition: JMessage.hh:65
Scanning of objects from multiple files according a format that follows from the extension of each fi...
Auxiliary class to build the supernova trigger dataset.
Definition: JSupernova.hh:121
Utility class to parse command line options.
ROOT TTree parameter settings.
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:60
Auxiliary class to manage a set of vetoes.
Definition: JSupernova.hh:237
int main(int argc, char *argv[])
Definition: Main.cpp:15