Jpp - the software that should make you happy
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
JAnglerFish.cc
Go to the documentation of this file.
1 #include <string>
2 #include <iostream>
3 #include <iomanip>
4 
5 #include "TROOT.h"
6 #include "TFile.h"
7 #include "TH1D.h"
8 
13 
14 #include "JDAQ/JDAQTimesliceIO.hh"
16 #include "JDAQ/JDAQEventIO.hh"
18 
21 #include "JSupport/JSupport.hh"
22 #include "JSupport/JMeta.hh"
23 
24 #include "JROOT/JManager.hh"
25 
26 #include "Jeep/JParser.hh"
27 #include "Jeep/JMessage.hh"
28 
30 /**
31  * \file
32  *
33  * Example application to search for sets of consecutive hits on a given PMT,\n
34  * that are separated in time by a constant value. The number of consecutive hits\n
35  * as well as the time difference are selected from the command line.\n
36  * The searches are done within the JDAQSnapshotHits from JDAQEvents\
37  * as well as within the JDAQHits from the JDAQTimeslice data.
38  *
39  * \author rgruiz
40  */
41 int main(int argc, char **argv)
42 {
43  using namespace std;
44  using namespace JPP;
45  using namespace KM3NETDAQ;
46 
48 
49  string inputFile;
50  JLimit_t numberOfEvents;
52  string detectorFile;
53  int deltaT;
54  int multiplicity;
55  int debug;
56 
57  try {
58 
59  JParser<> zap("Program to search for sets of consecutive hits on a given PMT that are equidistant in time.");
60 
61  zap['f'] = make_field(inputFile);
62  zap['a'] = make_field(detectorFile);
63  zap['t'] = make_field(deltaT) = 255;
64  zap['m'] = make_field(multiplicity) = 3;
65  zap['o'] = make_field(outputFile) = "anglerfish.root";
66  zap['n'] = make_field(numberOfEvents) = JLimit::max();
67  zap['d'] = make_field(debug) = 1;
68 
69  zap(argc, argv);
70  }
71  catch(const exception& error) {
72  FATAL(error.what() << endl);
73  }
74 
75  outputFile.open();
76 
77  if (!outputFile.is_open()) FATAL("Error opening file " << outputFile << endl);
78 
79  outputFile.put(JMeta(argc, argv));
80 
82 
83  try {
84  load (detectorFile, detector);
85  }
86  catch(const JException & error) {
87  cerr << "FATAL ERROR. Could not open detector file '" << detectorFile << "'." << endl;
88  exit(1);
89  }
90 
91  const JStringRouter stringRouter(detector);
92 
93  TH1D h1 ("MEvt", "", 50, 0.5, 50.5);
94 
95  TH2D h2 ("FEvt", NULL,
98 
99  TH1D h3 ("ML0" , "", 50, 0.5, 50.5);
100 
101  TH2D h4 ("FL0", NULL,
104 
105  for (Int_t i = 1; i <= h2.GetXaxis()->GetNbins(); ++i) {
106  h2.GetXaxis()->SetBinLabel(i, MAKE_CSTRING(stringRouter.at(i-1)));
107  }
108 
109  for (Int_t i = 1; i <= h2.GetXaxis()->GetNbins(); ++i) {
110  h4.GetXaxis()->SetBinLabel(i, MAKE_CSTRING(stringRouter.at(i-1)));
111  }
112 
113  const JDAQHitRouter hitRouter(detector);
114 
115  {
116  typedef JDAQSnapshotHit JHit_t;
117  typedef vector<JHit_t> JHitBuffer;
118 
119  for (JMultipleFileScanner<JDAQEvent> in(inputFile, numberOfEvents); in.hasNext(); ){
120 
121  JDAQEvent* event = in.next();
122  JHitBuffer& buffer = event->getHits<JHit_t>();
123 
124  std::sort(buffer.begin(),buffer.end(),less<JDAQKeyHit>());
125 
126  for (JHitBuffer::const_iterator p = buffer.begin() ; p != buffer.end() ; ) {
127 
128  JHitBuffer::const_iterator q = p;
129 
130  for (JDAQHit::JTDC_t t1 = p->getT();
131  (++q != buffer.end() &&
132  q->getModuleID() == p->getModuleID() &&
133  q->getPMT() == p->getPMT() &&
134  q->getT() == t1 + deltaT);
135  t1 = q->getT()) {}
136 
137  int d = distance(p,q);
138 
139  const JPMTChannel& channel = hitRouter.getPMTChannel(*p);
140 
141  if (d >= multiplicity)
142  h2.Fill((double) stringRouter.getIndex(channel.getString()), (double) channel.getFloor());
143 
144  if (d >= 2)
145  h1.Fill(d);
146 
147  p=q;
148  }
149  }
150  }
151 
152  {
153  typedef JDAQHit JHit_t;
154  typedef vector<JHit_t> JHitBuffer;
155  typedef JSuperFrame2D<JHit_t> JSuperFrame2D_t;
156 
157  const JModuleRouter moduleRouter(detector);
158 
159  for (JMultipleFileScanner<JDAQTimeslice> in(inputFile, numberOfEvents); in.hasNext(); ) {
160 
161  JDAQTimeslice* slice = in.next();
162 
163  for(JDAQTimeslice::iterator frame = slice->begin(); frame != slice->end() ; ++frame){
164 
165  const JModule& module = moduleRouter.getModule(frame->getModuleID());
166  JSuperFrame2D_t& buffer = JSuperFrame2D_t::demultiplex(*frame, module);
167 
168  for (JSuperFrame2D_t::iterator pmt = buffer.begin() ; pmt != buffer.end() ; ++pmt) {
169 
170  pmt->pop_back();
171 
172  for (JHitBuffer::const_iterator p = pmt->begin() ; p != pmt->end() ;) {
173 
174  JHitBuffer::const_iterator q = p;
175 
176  for (JDAQHit::JTDC_t t1 = p->getT();
177  (++q != pmt->end() &&
178  q->getT() == t1 + deltaT);
179  t1 = q->getT()) {}
180 
181  int d = distance(p,q);
182 
183  const JPMTChannel& channel = hitRouter.getPMTChannel(JDAQKeyHit(frame->getModuleID(),*p));
184 
185  if (d >= multiplicity)
186  h4.Fill((double) stringRouter.getIndex(channel.getString()), (double) channel.getFloor());
187 
188  if (d >= 2)
189  h3.Fill(d);
190 
191  p=q;
192  }
193  }
194  }
195  }
196  }
197 
198  outputFile.put(h1);
199  outputFile.put(h2);
200  outputFile.put(h3);
201  outputFile.put(h4);
202 
203  outputFile.close();
204 }
Auxiliary class to set-up Hit.
Definition: JHit_t.hh:25
Auxiliary class for ROOT I/O of application specific meta data.
Definition: JMeta.hh:70
Object writing to file.
Utility class to parse command line options.
Definition: JParser.hh:1500
General exception.
Definition: JException.hh:23
Direct access to PMT data in detector data structure for DAQ hits.
DAQ key hit.
Definition: JDAQKeyHit.hh:19
Auxiliary class to uniquely identify PMT readout channel.
Definition: JPMTChannel.hh:28
int main(int argc, char *argv[])
Definition: Main.cc:15
ROOT TTree parameter settings of various packages.
int getFloor() const
Get floor number.
Definition: JLocation.hh:145
const JModule & getModule(const JObjectID &id) const
Get module parameters.
Data structure for a composite optical module.
Definition: JModule.hh:57
std::vector< T >::difference_type distance(typename std::vector< T >::const_iterator first, typename PhysicsEvent::const_iterator< T > second)
Specialisation of STL distance.
Detector data structure.
Definition: JDetector.hh:80
Recording of objects on file according a format that follows from the file name extension.
Router for direct addressing of module data in detector data structure.
#define MAKE_CSTRING(A)
Make C-string.
Definition: JPrint.hh:151
exit
Definition: JPizza.sh:36
then for HISTOGRAM in h0 h1
Definition: JMatrixNZ.sh:69
Dynamic ROOT object management.
int getIndex(const T &value) const
Get index of given value.
string outputFile
unsigned int JTDC_t
leading edge [ns]
Definition: JDAQHit.hh:39
Simple wrapper around JModuleRouter class for direct addressing of PMT data in detector data structur...
Type list.
Definition: JTypeList.hh:22
Auxiliary class for defining the range of iterations of objects.
Definition: JLimit.hh:41
Detector file.
Definition: JHead.hh:196
Hit data structure.
Definition: JDAQHit.hh:34
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1961
ROOT I/O of application specific meta data.
int getNumberOfFloors(const JDetector &detector)
Get number of floors.
Data time slice.
int debug
debug level
Definition: JSirene.cc:63
General purpose messaging.
static const JStringCounter getNumberOfStrings
Function object to count unique strings.
Direct access to string in detector data structure.
#define FATAL(A)
Definition: JMessage.hh:67
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.
int getString() const
Get string number.
Definition: JLocation.hh:134
then JMuonMCEvt f $INPUT_FILE o $INTERMEDIATE_FILE d
Definition: JMuonPath.sh:45
void load(const std::string &file_name, JDetector &detector)
Load detector from input file.
Router for mapping of string identifier to index.
General purpose class for object reading from a list of file names.
Utility class to parse command line options.
2-dimensional frame with time calibrated data from one optical module.
do set_variable DETECTOR_TXT $WORKDIR detector
then fatal Wrong number of arguments fi set_variable DETECTOR $argv[1] set_variable INPUT_FILE $argv[2] eval JPrintDetector a $DETECTOR O IDENTIFIER eval JPrintDetector a $DETECTOR O SUMMARY source JAcoustics sh $DETECTOR_ID CHECK_EXIT_CODE typeset A TRIPODS get_tripods $WORKDIR tripod txt TRIPODS for EMITTER in
Definition: JCanberra.sh:38
JPMTChannel getPMTChannel(const JDAQKeyHit &hit) const
Get PMT channel.